/[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.21 - (hide annotations) (download)
Mon May 24 14:15:15 1999 UTC (25 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint23
Changes since 1.20: +2 -2 lines
Moved the final exchange of pressure (cg2d_x or cg3d_x) from the
solve to solve_for_pressure.F so that the adjoint knows whats
going on.

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

  ViewVC Help
Powered by ViewVC 1.1.22