/[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.35 - (hide annotations) (download)
Thu May 30 02:26:01 2002 UTC (22 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint45b_post, checkpoint45c_post
Changes since 1.34: +5 -73 lines
Cleaned exchange calls.

1 heimbach 1.35 C $Header: /u/gcmpack/MITgcm/model/src/cg2d.F,v 1.34 2001/09/26 18:09:14 cnh 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 cnh 1.14 _RL err
88 jmc 1.28 _RL eta_qrN
89     _RL eta_qrNM1
90 cnh 1.14 _RL cgBeta
91     _RL alpha
92     _RL sumRHS
93     _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     DO J=1,sNy
174     DO I=1,sNx
175     cg2d_s(I,J,bi,bj) = 0.
176     cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
177     & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
178     & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
179     & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
180     & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
181     & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
182     & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
183     & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
184 cnh 1.4 & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
185 jmc 1.29 & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
186 cnh 1.4 & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
187     & )
188 cnh 1.1 err = err +
189     & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
190     sumRHS = sumRHS +
191     & cg2d_b(I,J,bi,bj)
192     ENDDO
193     ENDDO
194     ENDDO
195     ENDDO
196 cnh 1.13 C _EXCH_XY_R8( cg2d_r, myThid )
197 adcroft 1.23 #ifdef LETS_MAKE_JAM
198     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
199     #else
200 heimbach 1.35 CALL EXCH_XY_RL( cg2d_r, myThid )
201 adcroft 1.23 #endif
202 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
203 adcroft 1.23 #ifdef LETS_MAKE_JAM
204     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
205     #else
206 heimbach 1.35 CALL EXCH_XY_RL( cg2d_s, myThid )
207 adcroft 1.23 #endif
208 adcroft 1.33 _GLOBAL_SUM_R8( sumRHS, myThid )
209     _GLOBAL_SUM_R8( err , myThid )
210     err = SQRT(err)
211     actualIts = 0
212     actualResidual = err
213 cnh 1.13
214     _BEGIN_MASTER( myThid )
215 adcroft 1.33 write(*,'(A,1P2E22.14)')' cg2d: Sum(rhs),rhsMax = ',
216     & sumRHS,rhsMax
217 cnh 1.13 _END_MASTER( )
218 cnh 1.1 C _BARRIER
219 adcroft 1.30 c _BEGIN_MASTER( myThid )
220 heimbach 1.31 c WRITE(*,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
221 adcroft 1.30 c & actualIts, actualResidual
222     c _END_MASTER( )
223 adcroft 1.33 firstResidual=actualResidual
224    
225     IF ( err .LT. cg2dTolerance ) GOTO 11
226 cnh 1.1
227     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
228 adcroft 1.30 DO 10 it2d=1, numIters
229 cnh 1.1
230     CcnhDebugStarts
231 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' residual = ',
232 cnh 1.14 C & actualResidual
233 cnh 1.1 CcnhDebugEnds
234     C-- Solve preconditioning equation and update
235     C-- conjugate direction vector "s".
236 jmc 1.28 eta_qrN = 0. _d 0
237 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
238     DO bi=myBxLo(myThid),myBxHi(myThid)
239     DO J=1,sNy
240     DO I=1,sNx
241     cg2d_q(I,J,bi,bj) =
242 cnh 1.3 & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
243     & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
244     & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
245     & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
246     & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
247 cnh 1.4 CcnhDebugStarts
248     C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
249     CcnhDebugEnds
250 jmc 1.28 eta_qrN = eta_qrN
251 cnh 1.1 & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
252     ENDDO
253     ENDDO
254     ENDDO
255     ENDDO
256    
257 jmc 1.28 _GLOBAL_SUM_R8(eta_qrN, myThid)
258 cnh 1.1 CcnhDebugStarts
259 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
260 cnh 1.1 CcnhDebugEnds
261 jmc 1.28 cgBeta = eta_qrN/eta_qrNM1
262 cnh 1.1 CcnhDebugStarts
263 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
264 cnh 1.1 CcnhDebugEnds
265 jmc 1.28 eta_qrNM1 = eta_qrN
266 cnh 1.1
267     DO bj=myByLo(myThid),myByHi(myThid)
268     DO bi=myBxLo(myThid),myBxHi(myThid)
269     DO J=1,sNy
270     DO I=1,sNx
271 cnh 1.14 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj)
272     & + cgBeta*cg2d_s(I,J,bi,bj)
273 cnh 1.1 ENDDO
274     ENDDO
275     ENDDO
276     ENDDO
277    
278     C-- Do exchanges that require messages i.e. between
279     C-- processes.
280 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
281 adcroft 1.23 #ifdef LETS_MAKE_JAM
282     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
283     #else
284 heimbach 1.35 CALL EXCH_XY_RL( cg2d_s, myThid )
285 adcroft 1.23 #endif
286 cnh 1.1
287     C== Evaluate laplace operator on conjugate gradient vector
288     C== q = A.s
289     alpha = 0. _d 0
290     DO bj=myByLo(myThid),myByHi(myThid)
291     DO bi=myBxLo(myThid),myBxHi(myThid)
292     DO J=1,sNy
293     DO I=1,sNx
294     cg2d_q(I,J,bi,bj) =
295     & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
296     & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
297     & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
298     & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
299     & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
300     & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
301     & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
302     & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
303 jmc 1.29 & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
304 cnh 1.4 & cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
305 cnh 1.1 alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
306     ENDDO
307     ENDDO
308     ENDDO
309     ENDDO
310 adcroft 1.20 _GLOBAL_SUM_R8(alpha,myThid)
311 cnh 1.1 CcnhDebugStarts
312 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
313 cnh 1.1 CcnhDebugEnds
314 jmc 1.28 alpha = eta_qrN/alpha
315 cnh 1.1 CcnhDebugStarts
316 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
317 cnh 1.1 CcnhDebugEnds
318    
319     C== Update solution and residual vectors
320     C Now compute "interior" points.
321     err = 0. _d 0
322     DO bj=myByLo(myThid),myByHi(myThid)
323     DO bi=myBxLo(myThid),myBxHi(myThid)
324     DO J=1,sNy
325     DO I=1,sNx
326     cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
327     cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
328     err = err+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
329     ENDDO
330     ENDDO
331     ENDDO
332     ENDDO
333    
334 adcroft 1.20 _GLOBAL_SUM_R8( err , myThid )
335 cnh 1.1 err = SQRT(err)
336     actualIts = it2d
337     actualResidual = err
338 adcroft 1.33 IF ( err .LT. cg2dTolerance ) GOTO 11
339 cnh 1.13 C _EXCH_XY_R8(cg2d_r, myThid )
340 adcroft 1.23 #ifdef LETS_MAKE_JAM
341     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
342     #else
343 heimbach 1.35 CALL EXCH_XY_RL( cg2d_r, myThid )
344 adcroft 1.23 #endif
345 cnh 1.13
346 cnh 1.1 10 CONTINUE
347     11 CONTINUE
348    
349 adcroft 1.33 IF (cg2dNormaliseRHS) THEN
350 cnh 1.1 C-- Un-normalise the answer
351 adcroft 1.33 DO bj=myByLo(myThid),myByHi(myThid)
352     DO bi=myBxLo(myThid),myBxHi(myThid)
353     DO J=1,sNy
354     DO I=1,sNx
355     cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
356     ENDDO
357     ENDDO
358 cnh 1.1 ENDDO
359     ENDDO
360 adcroft 1.33 ENDIF
361 cnh 1.1
362 adcroft 1.22 C The following exchange was moved up to solve_for_pressure
363     C for compatibility with TAMC.
364     C _EXCH_XY_R8(cg2d_x, myThid )
365 adcroft 1.30 c _BEGIN_MASTER( myThid )
366 heimbach 1.31 c WRITE(*,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
367 adcroft 1.30 c & actualIts, actualResidual
368     c _END_MASTER( )
369    
370     C-- Return parameters to caller
371 adcroft 1.33 lastResidual=actualResidual
372 adcroft 1.30 numIters=actualIts
373 cnh 1.1
374     CcnhDebugStarts
375 cnh 1.7 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
376 cnh 1.1 C err = 0. _d 0
377     C DO bj=myByLo(myThid),myByHi(myThid)
378     C DO bi=myBxLo(myThid),myBxHi(myThid)
379     C DO J=1,sNy
380     C DO I=1,sNx
381     C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
382     C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
383     C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
384     C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
385     C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
386     C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
387     C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
388     C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
389     C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
390     C err = err +
391     C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
392     C ENDDO
393     C ENDDO
394     C ENDDO
395     C ENDDO
396 adcroft 1.20 C _GLOBAL_SUM_R8( err , myThid )
397 heimbach 1.31 C write(*,*) 'cg2d: Ax - b = ',SQRT(err)
398 cnh 1.1 CcnhDebugEnds
399    
400 adcroft 1.19 RETURN
401 cnh 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22