/[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.34.4.1 - (show annotations) (download)
Thu May 30 02:46:13 2002 UTC (22 years ago) by heimbach
Branch: release1
CVS Tags: release1_p13_pre, release1_p13, release1_p8, release1_p9, release1_p3, release1_p4, release1_p5, release1_p6, release1_p7, release1_p11, release1_p12, release1_p10, release1_p16, release1_p17, release1_p14, release1_p15, release1_p12_pre
Branch point for: release1_50yr
Changes since 1.34: +4 -74 lines
Cleaned exchange calls in cg2d.

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

  ViewVC Help
Powered by ViewVC 1.1.22