/[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.23 - (hide annotations) (download)
Tue Mar 14 16:09:41 2000 UTC (24 years, 2 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint25
Changes since 1.22: +22 -2 lines
Added "JAM" routines for use with Artic network (Hyades cluster).

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

  ViewVC Help
Powered by ViewVC 1.1.22