/[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.26 - (hide annotations) (download)
Fri Feb 2 21:04:47 2001 UTC (23 years, 3 months ago) by adcroft
Branch: MAIN
Changes since 1.25: +2 -1 lines
Merged changes from branch "branch-atmos-merge" into MAIN (checkpoint34)
 - substantial modifications to algorithm sequence (dynamics.F)
 - packaged OBCS, Shapiro filter, Zonal filter, Atmospheric Physics

1 adcroft 1.26 C $Header: /u/gcmpack/models/MITgcmUV/model/src/cg2d.F,v 1.25.2.1 2001/01/23 15:54:27 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 adcroft 1.24 C CHARACTER*(MAX_LEN_FNAM) suff
86 cnh 1.12 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.25 C _GLOBAL_MAX_R8( rhsMax, myThid )
114     rhsMax=1.
115 adcroft 1.23 #else
116     _GLOBAL_MAX_R8( rhsMax, myThid )
117 adcroft 1.26 Catm rhsMax=1.
118 adcroft 1.23 #endif
119 cnh 1.1 rhsNorm = 1. _d 0
120     IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
121     DO bj=myByLo(myThid),myByHi(myThid)
122     DO bi=myBxLo(myThid),myBxHi(myThid)
123     DO J=1,sNy
124     DO I=1,sNx
125     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
126     cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
127     ENDDO
128     ENDDO
129     ENDDO
130     ENDDO
131    
132     C-- Update overlaps
133     _EXCH_XY_R8( cg2d_b, myThid )
134     _EXCH_XY_R8( cg2d_x, myThid )
135     CcnhDebugStarts
136 cnh 1.11 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
137 cnh 1.12 C suff = 'normalised'
138     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
139 cnh 1.1 CcnhDebugEnds
140    
141     C-- Initial residual calculation
142     err = 0. _d 0
143     sumRHS = 0. _d 0
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_s(I,J,bi,bj) = 0.
149     cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
150     & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
151     & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
152     & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
153     & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
154     & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
155     & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
156     & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
157 cnh 1.4 & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
158 cnh 1.10 & -freeSurfFac*_rA(i,j,bi,bj)* horiVertRatio*
159 cnh 1.4 & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
160     & )
161 cnh 1.1 err = err +
162     & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
163     sumRHS = sumRHS +
164     & cg2d_b(I,J,bi,bj)
165     ENDDO
166     ENDDO
167     ENDDO
168     ENDDO
169 cnh 1.13 C _EXCH_XY_R8( cg2d_r, myThid )
170 adcroft 1.23 #ifdef LETS_MAKE_JAM
171     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
172     #else
173 cnh 1.13 OLw = 1
174     OLe = 1
175     OLn = 1
176     OLs = 1
177     exchWidthX = 1
178     exchWidthY = 1
179     myNz = 1
180     CALL EXCH_RL( cg2d_r,
181     I OLw, OLe, OLs, OLn, myNz,
182     I exchWidthX, exchWidthY,
183     I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
184 adcroft 1.23 #endif
185 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
186 adcroft 1.23 #ifdef LETS_MAKE_JAM
187     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
188     #else
189 cnh 1.13 OLw = 1
190     OLe = 1
191     OLn = 1
192     OLs = 1
193     exchWidthX = 1
194     exchWidthY = 1
195     myNz = 1
196     CALL EXCH_RL( cg2d_s,
197     I OLw, OLe, OLs, OLn, myNz,
198     I exchWidthX, exchWidthY,
199     I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
200 adcroft 1.23 #endif
201 adcroft 1.20 _GLOBAL_SUM_R8( sumRHS, myThid )
202 cnh 1.1 C WRITE(6,*) ' mythid, err = ', mythid, SQRT(err)
203 adcroft 1.20 _GLOBAL_SUM_R8( err , myThid )
204 cnh 1.13
205     _BEGIN_MASTER( myThid )
206 adcroft 1.19 write(0,'(A,1PE30.14)') ' cg2d: Sum(rhs) = ',sumRHS
207 cnh 1.13 _END_MASTER( )
208 cnh 1.1
209     actualIts = 0
210     actualResidual = SQRT(err)
211     C _BARRIER
212     _BEGIN_MASTER( myThid )
213 adcroft 1.17 WRITE(0,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
214 cnh 1.14 & actualIts, actualResidual
215 cnh 1.1 _END_MASTER( )
216    
217     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
218     DO 10 it2d=1, cg2dMaxIters
219    
220     CcnhDebugStarts
221 cnh 1.14 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' residual = ',
222     C & actualResidual
223 cnh 1.1 CcnhDebugEnds
224     IF ( err .LT. cg2dTargetResidual ) GOTO 11
225     C-- Solve preconditioning equation and update
226     C-- conjugate direction vector "s".
227     etaN = 0. _d 0
228     DO bj=myByLo(myThid),myByHi(myThid)
229     DO bi=myBxLo(myThid),myBxHi(myThid)
230     DO J=1,sNy
231     DO I=1,sNx
232     cg2d_q(I,J,bi,bj) =
233 cnh 1.3 & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
234     & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
235     & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
236     & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
237     & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
238 cnh 1.4 CcnhDebugStarts
239     C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
240     CcnhDebugEnds
241 cnh 1.1 etaN = etaN
242     & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
243     ENDDO
244     ENDDO
245     ENDDO
246     ENDDO
247    
248 adcroft 1.20 _GLOBAL_SUM_R8(etaN, myThid)
249 cnh 1.1 CcnhDebugStarts
250     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' etaN = ',etaN
251     CcnhDebugEnds
252     cgBeta = etaN/etaNM1
253     CcnhDebugStarts
254     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
255     CcnhDebugEnds
256     etaNM1 = etaN
257    
258     DO bj=myByLo(myThid),myByHi(myThid)
259     DO bi=myBxLo(myThid),myBxHi(myThid)
260     DO J=1,sNy
261     DO I=1,sNx
262 cnh 1.14 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj)
263     & + cgBeta*cg2d_s(I,J,bi,bj)
264 cnh 1.1 ENDDO
265     ENDDO
266     ENDDO
267     ENDDO
268    
269     C-- Do exchanges that require messages i.e. between
270     C-- processes.
271 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
272 adcroft 1.23 #ifdef LETS_MAKE_JAM
273     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
274     #else
275 cnh 1.13 OLw = 1
276     OLe = 1
277     OLn = 1
278     OLs = 1
279     exchWidthX = 1
280     exchWidthY = 1
281     myNz = 1
282     CALL EXCH_RL( cg2d_s,
283     I OLw, OLe, OLs, OLn, myNz,
284     I exchWidthX, exchWidthY,
285     I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
286 adcroft 1.23 #endif
287 cnh 1.1
288     C== Evaluate laplace operator on conjugate gradient vector
289     C== q = A.s
290     alpha = 0. _d 0
291     DO bj=myByLo(myThid),myByHi(myThid)
292     DO bi=myBxLo(myThid),myBxHi(myThid)
293     DO J=1,sNy
294     DO I=1,sNx
295     cg2d_q(I,J,bi,bj) =
296     & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
297     & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
298     & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
299     & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
300     & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
301     & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
302     & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
303     & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
304 cnh 1.10 & -freeSurfFac*_rA(i,j,bi,bj)* horiVertRatio*
305 cnh 1.4 & cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
306 cnh 1.1 alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
307     ENDDO
308     ENDDO
309     ENDDO
310     ENDDO
311 adcroft 1.20 _GLOBAL_SUM_R8(alpha,myThid)
312 cnh 1.1 CcnhDebugStarts
313     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
314     CcnhDebugEnds
315     alpha = etaN/alpha
316     CcnhDebugStarts
317     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
318     CcnhDebugEnds
319    
320     C== Update solution and residual vectors
321     C Now compute "interior" points.
322     err = 0. _d 0
323     DO bj=myByLo(myThid),myByHi(myThid)
324     DO bi=myBxLo(myThid),myBxHi(myThid)
325     DO J=1,sNy
326     DO I=1,sNx
327     cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
328     cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
329     err = err+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
330     ENDDO
331     ENDDO
332     ENDDO
333     ENDDO
334    
335 adcroft 1.20 _GLOBAL_SUM_R8( err , myThid )
336 cnh 1.1 err = SQRT(err)
337     actualIts = it2d
338     actualResidual = err
339     IF ( err .LT. cg2dTargetResidual ) GOTO 11
340 cnh 1.13 C _EXCH_XY_R8(cg2d_r, myThid )
341 adcroft 1.23 #ifdef LETS_MAKE_JAM
342     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
343     #else
344 cnh 1.13 OLw = 1
345     OLe = 1
346     OLn = 1
347     OLs = 1
348     exchWidthX = 1
349     exchWidthY = 1
350     myNz = 1
351     CALL EXCH_RL( cg2d_r,
352     I OLw, OLe, OLs, OLn, myNz,
353     I exchWidthX, exchWidthY,
354     I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
355 adcroft 1.23 #endif
356 cnh 1.13
357 cnh 1.1 10 CONTINUE
358     11 CONTINUE
359    
360     C-- Un-normalise the answer
361     DO bj=myByLo(myThid),myByHi(myThid)
362     DO bi=myBxLo(myThid),myBxHi(myThid)
363     DO J=1,sNy
364     DO I=1,sNx
365     cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
366     ENDDO
367     ENDDO
368     ENDDO
369     ENDDO
370    
371 adcroft 1.22 C The following exchange was moved up to solve_for_pressure
372     C for compatibility with TAMC.
373     C _EXCH_XY_R8(cg2d_x, myThid )
374 cnh 1.6 _BEGIN_MASTER( myThid )
375 adcroft 1.17 WRITE(0,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
376 cnh 1.14 & actualIts, actualResidual
377 cnh 1.6 _END_MASTER( )
378 cnh 1.1
379     CcnhDebugStarts
380 cnh 1.7 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
381 cnh 1.1 C err = 0. _d 0
382     C DO bj=myByLo(myThid),myByHi(myThid)
383     C DO bi=myBxLo(myThid),myBxHi(myThid)
384     C DO J=1,sNy
385     C DO I=1,sNx
386     C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
387     C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
388     C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
389     C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
390     C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
391     C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
392     C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
393     C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
394     C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
395     C err = err +
396     C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
397     C ENDDO
398     C ENDDO
399     C ENDDO
400     C ENDDO
401 adcroft 1.20 C _GLOBAL_SUM_R8( err , myThid )
402 cnh 1.1 C write(0,*) 'cg2d: Ax - b = ',SQRT(err)
403     CcnhDebugEnds
404    
405 adcroft 1.19 RETURN
406 cnh 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22