/[MITgcm]/MITgcm/model/src/cg2d.F
ViewVC logotype

Contents of /MITgcm/model/src/cg2d.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.29 - (show annotations) (download)
Thu Mar 8 20:49:08 2001 UTC (23 years, 2 months ago) by jmc
Branch: MAIN
Changes since 1.28: +4 -3 lines
change (*g) units of cg2d_x to have same units as all other potential Phi

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/cg2d.F,v 1.28 2001/03/06 17:02:57 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 SUBROUTINE CG2D(
7 I cg2d_b,
8 U cg2d_x,
9 I myThid )
10 C /==========================================================\
11 C | SUBROUTINE CG2D |
12 C | o Two-dimensional grid problem conjugate-gradient |
13 C | inverter (with preconditioner). |
14 C |==========================================================|
15 C | Con. grad is an iterative procedure for solving Ax = b. |
16 C | It requires the A be symmetric. |
17 C | This implementation assumes A is a five-diagonal |
18 C | matrix of the form that arises in the discrete |
19 C | representation of the del^2 operator in a |
20 C | two-dimensional space. |
21 C | Notes: |
22 C | ====== |
23 C | This implementation can support shared-memory |
24 C | multi-threaded execution. In order to do this COMMON |
25 C | blocks are used for many of the arrays - even ones that |
26 C | are only used for intermedaite results. This design is |
27 C | OK if you want to all the threads to collaborate on |
28 C | solving the same problem. On the other hand if you want |
29 C | the threads to solve several different problems |
30 C | concurrently this implementation will not work. |
31 C \==========================================================/
32 IMPLICIT NONE
33
34 C === Global data ===
35 #include "SIZE.h"
36 #include "EEPARAMS.h"
37 #include "PARAMS.h"
38 #include "GRID.h"
39 #include "CG2D_INTERNAL.h"
40 #include "SURFACE.h"
41
42 C === Routine arguments ===
43 C myThid - Thread on which I am working.
44 INTEGER myThid
45 _RL cg2d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
46 _RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
47
48
49 C === Local variables ====
50 C actualIts - Number of iterations taken
51 C actualResidual - residual
52 C bi - Block index in X and Y.
53 C bj
54 C eta_qrN - Used in computing search directions
55 C eta_qrNM1 suffix N and NM1 denote current and
56 C cgBeta previous iterations respectively.
57 C alpha
58 C sumRHS - Sum of right-hand-side. Sometimes this is a
59 C useful debuggin/trouble shooting diagnostic.
60 C For neumann problems sumRHS needs to be ~0.
61 C or they converge at a non-zero residual.
62 C err - Measure of residual of Ax - b, usually the norm.
63 C I, J, N - Loop counters ( N counts CG iterations )
64 INTEGER actualIts
65 _RL actualResidual
66 INTEGER bi, bj
67 INTEGER I, J, it2d
68 _RL err
69 _RL eta_qrN
70 _RL eta_qrNM1
71 _RL cgBeta
72 _RL alpha
73 _RL sumRHS
74 _RL rhsMax
75 _RL rhsNorm
76
77 INTEGER OLw
78 INTEGER OLe
79 INTEGER OLn
80 INTEGER OLs
81 INTEGER exchWidthX
82 INTEGER exchWidthY
83 INTEGER myNz
84
85
86 CcnhDebugStarts
87 C CHARACTER*(MAX_LEN_FNAM) suff
88 CcnhDebugEnds
89
90
91 C-- Initialise inverter
92 eta_qrNM1 = 1. _d 0
93
94 CcnhDebugStarts
95 C _EXCH_XY_R8( cg2d_b, myThid )
96 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid )
97 C suff = 'unnormalised'
98 C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
99 C STOP
100 CcnhDebugEnds
101
102 C-- Normalise RHS
103 rhsMax = 0. _d 0
104 DO bj=myByLo(myThid),myByHi(myThid)
105 DO bi=myBxLo(myThid),myBxHi(myThid)
106 DO J=1,sNy
107 DO I=1,sNx
108 cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm
109 rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax)
110 ENDDO
111 ENDDO
112 ENDDO
113 ENDDO
114 #ifdef LETS_MAKE_JAM
115 C _GLOBAL_MAX_R8( rhsMax, myThid )
116 rhsMax=1.
117 #else
118 _GLOBAL_MAX_R8( rhsMax, myThid )
119 Catm rhsMax=1.
120 #endif
121 rhsNorm = 1. _d 0
122 IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
123 DO bj=myByLo(myThid),myByHi(myThid)
124 DO bi=myBxLo(myThid),myBxHi(myThid)
125 DO J=1,sNy
126 DO I=1,sNx
127 cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
128 cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
129 ENDDO
130 ENDDO
131 ENDDO
132 ENDDO
133
134 C-- Update overlaps
135 _EXCH_XY_R8( cg2d_b, myThid )
136 _EXCH_XY_R8( cg2d_x, myThid )
137 CcnhDebugStarts
138 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
139 C suff = 'normalised'
140 C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
141 CcnhDebugEnds
142
143 C-- Initial residual calculation
144 err = 0. _d 0
145 sumRHS = 0. _d 0
146 DO bj=myByLo(myThid),myByHi(myThid)
147 DO bi=myBxLo(myThid),myBxHi(myThid)
148 DO J=1,sNy
149 DO I=1,sNx
150 cg2d_s(I,J,bi,bj) = 0.
151 cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
152 & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
153 & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
154 & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
155 & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
156 & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
157 & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
158 & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
159 & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
160 & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
161 & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
162 & )
163 err = err +
164 & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
165 sumRHS = sumRHS +
166 & cg2d_b(I,J,bi,bj)
167 ENDDO
168 ENDDO
169 ENDDO
170 ENDDO
171 C _EXCH_XY_R8( cg2d_r, myThid )
172 #ifdef LETS_MAKE_JAM
173 CALL EXCH_XY_O1_R8_JAM( cg2d_r )
174 #else
175 OLw = 1
176 OLe = 1
177 OLn = 1
178 OLs = 1
179 exchWidthX = 1
180 exchWidthY = 1
181 myNz = 1
182 CALL EXCH_RL( cg2d_r,
183 I OLw, OLe, OLs, OLn, myNz,
184 I exchWidthX, exchWidthY,
185 I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
186 #endif
187 C _EXCH_XY_R8( cg2d_s, myThid )
188 #ifdef LETS_MAKE_JAM
189 CALL EXCH_XY_O1_R8_JAM( cg2d_s )
190 #else
191 OLw = 1
192 OLe = 1
193 OLn = 1
194 OLs = 1
195 exchWidthX = 1
196 exchWidthY = 1
197 myNz = 1
198 CALL EXCH_RL( cg2d_s,
199 I OLw, OLe, OLs, OLn, myNz,
200 I exchWidthX, exchWidthY,
201 I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
202 #endif
203 _GLOBAL_SUM_R8( sumRHS, myThid )
204 C WRITE(6,*) ' mythid, err = ', mythid, SQRT(err)
205 _GLOBAL_SUM_R8( err , myThid )
206
207 _BEGIN_MASTER( myThid )
208 write(0,'(A,1PE30.14)') ' cg2d: Sum(rhs) = ',sumRHS
209 _END_MASTER( )
210
211 actualIts = 0
212 actualResidual = SQRT(err)
213 C _BARRIER
214 _BEGIN_MASTER( myThid )
215 WRITE(0,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
216 & actualIts, actualResidual
217 _END_MASTER( )
218
219 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
220 DO 10 it2d=1, cg2dMaxIters
221
222 CcnhDebugStarts
223 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' residual = ',
224 C & actualResidual
225 CcnhDebugEnds
226 IF ( err .LT. cg2dTargetResidual ) GOTO 11
227 C-- Solve preconditioning equation and update
228 C-- conjugate direction vector "s".
229 eta_qrN = 0. _d 0
230 DO bj=myByLo(myThid),myByHi(myThid)
231 DO bi=myBxLo(myThid),myBxHi(myThid)
232 DO J=1,sNy
233 DO I=1,sNx
234 cg2d_q(I,J,bi,bj) =
235 & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
236 & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
237 & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
238 & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
239 & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
240 CcnhDebugStarts
241 C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
242 CcnhDebugEnds
243 eta_qrN = eta_qrN
244 & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
245 ENDDO
246 ENDDO
247 ENDDO
248 ENDDO
249
250 _GLOBAL_SUM_R8(eta_qrN, myThid)
251 CcnhDebugStarts
252 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
253 CcnhDebugEnds
254 cgBeta = eta_qrN/eta_qrNM1
255 CcnhDebugStarts
256 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
257 CcnhDebugEnds
258 eta_qrNM1 = eta_qrN
259
260 DO bj=myByLo(myThid),myByHi(myThid)
261 DO bi=myBxLo(myThid),myBxHi(myThid)
262 DO J=1,sNy
263 DO I=1,sNx
264 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj)
265 & + cgBeta*cg2d_s(I,J,bi,bj)
266 ENDDO
267 ENDDO
268 ENDDO
269 ENDDO
270
271 C-- Do exchanges that require messages i.e. between
272 C-- processes.
273 C _EXCH_XY_R8( cg2d_s, myThid )
274 #ifdef LETS_MAKE_JAM
275 CALL EXCH_XY_O1_R8_JAM( cg2d_s )
276 #else
277 OLw = 1
278 OLe = 1
279 OLn = 1
280 OLs = 1
281 exchWidthX = 1
282 exchWidthY = 1
283 myNz = 1
284 CALL EXCH_RL( cg2d_s,
285 I OLw, OLe, OLs, OLn, myNz,
286 I exchWidthX, exchWidthY,
287 I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
288 #endif
289
290 C== Evaluate laplace operator on conjugate gradient vector
291 C== q = A.s
292 alpha = 0. _d 0
293 DO bj=myByLo(myThid),myByHi(myThid)
294 DO bi=myBxLo(myThid),myBxHi(myThid)
295 DO J=1,sNy
296 DO I=1,sNx
297 cg2d_q(I,J,bi,bj) =
298 & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
299 & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
300 & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
301 & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
302 & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
303 & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
304 & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
305 & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
306 & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
307 & cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
308 alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
309 ENDDO
310 ENDDO
311 ENDDO
312 ENDDO
313 _GLOBAL_SUM_R8(alpha,myThid)
314 CcnhDebugStarts
315 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
316 CcnhDebugEnds
317 alpha = eta_qrN/alpha
318 CcnhDebugStarts
319 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
320 CcnhDebugEnds
321
322 C== Update solution and residual vectors
323 C Now compute "interior" points.
324 err = 0. _d 0
325 DO bj=myByLo(myThid),myByHi(myThid)
326 DO bi=myBxLo(myThid),myBxHi(myThid)
327 DO J=1,sNy
328 DO I=1,sNx
329 cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
330 cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
331 err = err+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
332 ENDDO
333 ENDDO
334 ENDDO
335 ENDDO
336
337 _GLOBAL_SUM_R8( err , myThid )
338 err = SQRT(err)
339 actualIts = it2d
340 actualResidual = err
341 IF ( err .LT. cg2dTargetResidual ) GOTO 11
342 C _EXCH_XY_R8(cg2d_r, myThid )
343 #ifdef LETS_MAKE_JAM
344 CALL EXCH_XY_O1_R8_JAM( cg2d_r )
345 #else
346 OLw = 1
347 OLe = 1
348 OLn = 1
349 OLs = 1
350 exchWidthX = 1
351 exchWidthY = 1
352 myNz = 1
353 CALL EXCH_RL( cg2d_r,
354 I OLw, OLe, OLs, OLn, myNz,
355 I exchWidthX, exchWidthY,
356 I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
357 #endif
358
359 10 CONTINUE
360 11 CONTINUE
361
362 C-- Un-normalise the answer
363 DO bj=myByLo(myThid),myByHi(myThid)
364 DO bi=myBxLo(myThid),myBxHi(myThid)
365 DO J=1,sNy
366 DO I=1,sNx
367 cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
368 ENDDO
369 ENDDO
370 ENDDO
371 ENDDO
372
373 C The following exchange was moved up to solve_for_pressure
374 C for compatibility with TAMC.
375 C _EXCH_XY_R8(cg2d_x, myThid )
376 _BEGIN_MASTER( myThid )
377 WRITE(0,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
378 & actualIts, actualResidual
379 _END_MASTER( )
380
381 CcnhDebugStarts
382 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
383 C err = 0. _d 0
384 C DO bj=myByLo(myThid),myByHi(myThid)
385 C DO bi=myBxLo(myThid),myBxHi(myThid)
386 C DO J=1,sNy
387 C DO I=1,sNx
388 C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
389 C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
390 C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
391 C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
392 C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
393 C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
394 C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
395 C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
396 C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
397 C err = err +
398 C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
399 C ENDDO
400 C ENDDO
401 C ENDDO
402 C ENDDO
403 C _GLOBAL_SUM_R8( err , myThid )
404 C write(0,*) 'cg2d: Ax - b = ',SQRT(err)
405 CcnhDebugEnds
406
407 RETURN
408 END

  ViewVC Help
Powered by ViewVC 1.1.22