/[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.29 - (hide 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 jmc 1.29 C $Header: /u/gcmpack/models/MITgcmUV/model/src/cg2d.F,v 1.28 2001/03/06 17:02:57 jmc Exp $
2 jmc 1.28 C $Name: $
3 cnh 1.1
4 cnh 1.16 #include "CPP_OPTIONS.h"
5 cnh 1.1
6     SUBROUTINE CG2D(
7 cnh 1.14 I cg2d_b,
8     U cg2d_x,
9 cnh 1.1 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 adcroft 1.18 IMPLICIT NONE
33 cnh 1.1
34     C === Global data ===
35     #include "SIZE.h"
36     #include "EEPARAMS.h"
37     #include "PARAMS.h"
38 cnh 1.4 #include "GRID.h"
39 cnh 1.14 #include "CG2D_INTERNAL.h"
40 jmc 1.29 #include "SURFACE.h"
41 cnh 1.1
42     C === Routine arguments ===
43     C myThid - Thread on which I am working.
44     INTEGER myThid
45 cnh 1.14 _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 cnh 1.1
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 jmc 1.28 C eta_qrN - Used in computing search directions
55     C eta_qrNM1 suffix N and NM1 denote current and
56 cnh 1.1 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 cnh 1.14 _RL actualResidual
66 cnh 1.1 INTEGER bi, bj
67     INTEGER I, J, it2d
68 cnh 1.14 _RL err
69 jmc 1.28 _RL eta_qrN
70     _RL eta_qrNM1
71 cnh 1.14 _RL cgBeta
72     _RL alpha
73     _RL sumRHS
74     _RL rhsMax
75     _RL rhsNorm
76 cnh 1.1
77 cnh 1.13 INTEGER OLw
78     INTEGER OLe
79     INTEGER OLn
80     INTEGER OLs
81     INTEGER exchWidthX
82     INTEGER exchWidthY
83     INTEGER myNz
84    
85    
86 cnh 1.12 CcnhDebugStarts
87 adcroft 1.24 C CHARACTER*(MAX_LEN_FNAM) suff
88 cnh 1.12 CcnhDebugEnds
89    
90    
91 cnh 1.1 C-- Initialise inverter
92 jmc 1.28 eta_qrNM1 = 1. _d 0
93 cnh 1.1
94 cnh 1.10 CcnhDebugStarts
95 cnh 1.11 C _EXCH_XY_R8( cg2d_b, myThid )
96     C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid )
97 cnh 1.12 C suff = 'unnormalised'
98     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
99 cnh 1.14 C STOP
100 cnh 1.10 CcnhDebugEnds
101    
102 cnh 1.1 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 adcroft 1.23 #ifdef LETS_MAKE_JAM
115 adcroft 1.25 C _GLOBAL_MAX_R8( rhsMax, myThid )
116     rhsMax=1.
117 adcroft 1.23 #else
118     _GLOBAL_MAX_R8( rhsMax, myThid )
119 adcroft 1.26 Catm rhsMax=1.
120 adcroft 1.23 #endif
121 cnh 1.1 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 cnh 1.11 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
139 cnh 1.12 C suff = 'normalised'
140     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
141 cnh 1.1 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 cnh 1.4 & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
160 jmc 1.29 & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
161 cnh 1.4 & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
162     & )
163 cnh 1.1 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 cnh 1.13 C _EXCH_XY_R8( cg2d_r, myThid )
172 adcroft 1.23 #ifdef LETS_MAKE_JAM
173     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
174     #else
175 cnh 1.13 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 adcroft 1.23 #endif
187 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
188 adcroft 1.23 #ifdef LETS_MAKE_JAM
189     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
190     #else
191 cnh 1.13 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 adcroft 1.23 #endif
203 adcroft 1.20 _GLOBAL_SUM_R8( sumRHS, myThid )
204 cnh 1.1 C WRITE(6,*) ' mythid, err = ', mythid, SQRT(err)
205 adcroft 1.20 _GLOBAL_SUM_R8( err , myThid )
206 cnh 1.13
207     _BEGIN_MASTER( myThid )
208 adcroft 1.19 write(0,'(A,1PE30.14)') ' cg2d: Sum(rhs) = ',sumRHS
209 cnh 1.13 _END_MASTER( )
210 cnh 1.1
211     actualIts = 0
212     actualResidual = SQRT(err)
213     C _BARRIER
214     _BEGIN_MASTER( myThid )
215 adcroft 1.17 WRITE(0,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
216 cnh 1.14 & actualIts, actualResidual
217 cnh 1.1 _END_MASTER( )
218    
219     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
220     DO 10 it2d=1, cg2dMaxIters
221    
222     CcnhDebugStarts
223 cnh 1.14 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' residual = ',
224     C & actualResidual
225 cnh 1.1 CcnhDebugEnds
226     IF ( err .LT. cg2dTargetResidual ) GOTO 11
227     C-- Solve preconditioning equation and update
228     C-- conjugate direction vector "s".
229 jmc 1.28 eta_qrN = 0. _d 0
230 cnh 1.1 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 cnh 1.3 & 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 cnh 1.4 CcnhDebugStarts
241     C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
242     CcnhDebugEnds
243 jmc 1.28 eta_qrN = eta_qrN
244 cnh 1.1 & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
245     ENDDO
246     ENDDO
247     ENDDO
248     ENDDO
249    
250 jmc 1.28 _GLOBAL_SUM_R8(eta_qrN, myThid)
251 cnh 1.1 CcnhDebugStarts
252 jmc 1.28 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
253 cnh 1.1 CcnhDebugEnds
254 jmc 1.28 cgBeta = eta_qrN/eta_qrNM1
255 cnh 1.1 CcnhDebugStarts
256     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
257     CcnhDebugEnds
258 jmc 1.28 eta_qrNM1 = eta_qrN
259 cnh 1.1
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 cnh 1.14 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj)
265     & + cgBeta*cg2d_s(I,J,bi,bj)
266 cnh 1.1 ENDDO
267     ENDDO
268     ENDDO
269     ENDDO
270    
271     C-- Do exchanges that require messages i.e. between
272     C-- processes.
273 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
274 adcroft 1.23 #ifdef LETS_MAKE_JAM
275     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
276     #else
277 cnh 1.13 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 adcroft 1.23 #endif
289 cnh 1.1
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 jmc 1.29 & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
307 cnh 1.4 & cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
308 cnh 1.1 alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
309     ENDDO
310     ENDDO
311     ENDDO
312     ENDDO
313 adcroft 1.20 _GLOBAL_SUM_R8(alpha,myThid)
314 cnh 1.1 CcnhDebugStarts
315     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
316     CcnhDebugEnds
317 jmc 1.28 alpha = eta_qrN/alpha
318 cnh 1.1 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 adcroft 1.20 _GLOBAL_SUM_R8( err , myThid )
338 cnh 1.1 err = SQRT(err)
339     actualIts = it2d
340     actualResidual = err
341     IF ( err .LT. cg2dTargetResidual ) GOTO 11
342 cnh 1.13 C _EXCH_XY_R8(cg2d_r, myThid )
343 adcroft 1.23 #ifdef LETS_MAKE_JAM
344     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
345     #else
346 cnh 1.13 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 adcroft 1.23 #endif
358 cnh 1.13
359 cnh 1.1 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 adcroft 1.22 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 cnh 1.6 _BEGIN_MASTER( myThid )
377 adcroft 1.17 WRITE(0,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
378 cnh 1.14 & actualIts, actualResidual
379 cnh 1.6 _END_MASTER( )
380 cnh 1.1
381     CcnhDebugStarts
382 cnh 1.7 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
383 cnh 1.1 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 adcroft 1.20 C _GLOBAL_SUM_R8( err , myThid )
404 cnh 1.1 C write(0,*) 'cg2d: Ax - b = ',SQRT(err)
405     CcnhDebugEnds
406    
407 adcroft 1.19 RETURN
408 cnh 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22