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

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

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


Revision 1.39 - (hide annotations) (download)
Thu Apr 22 22:43:37 2004 UTC (20 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint56b_post, checkpoint54d_post, checkpoint54e_post, checkpoint55, checkpoint54, checkpoint57, checkpoint56, checkpoint53, checkpoint54f_post, checkpoint55i_post, checkpoint55c_post, checkpoint53d_post, checkpoint57a_post, checkpoint54b_post, checkpoint55g_post, checkpoint55d_post, checkpoint54a_pre, checkpoint53c_post, checkpoint55d_pre, checkpoint55j_post, checkpoint54a_post, checkpoint55h_post, checkpoint55b_post, checkpoint53a_post, checkpoint55f_post, checkpoint53g_post, checkpoint56a_post, checkpoint53f_post, checkpoint52n_post, checkpoint53b_pre, checkpoint56c_post, checkpoint57a_pre, checkpoint55a_post, checkpoint53b_post, checkpoint53d_pre, checkpoint55e_post, checkpoint54c_post
Changes since 1.38: +2 -2 lines
skip writing "Sum(rhs),rhsMax = .."  only if debugLevel < 0

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

  ViewVC Help
Powered by ViewVC 1.1.22