/[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.13 - (hide annotations) (download)
Tue Sep 29 18:50:56 1998 UTC (25 years, 8 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint15
Changes since 1.12: +63 -5 lines
Changes for new exchange routines which do general tile <-> tile
connectivity, variable width overlap regions and provide
hooks for shared memory  and DMA protocols like Arctic, Memory Channel
etc..

1 cnh 1.13 C $Header: /u/gcmpack/models/MITgcmUV/model/src/cg2d.F,v 1.12 1998/09/08 01:37:49 cnh Exp $
2 cnh 1.1
3     #include "CPP_EEOPTIONS.h"
4    
5     SUBROUTINE CG2D(
6     I myThid )
7     C /==========================================================\
8     C | SUBROUTINE CG2D |
9     C | o Two-dimensional grid problem conjugate-gradient |
10     C | inverter (with preconditioner). |
11     C |==========================================================|
12     C | Con. grad is an iterative procedure for solving Ax = b. |
13     C | It requires the A be symmetric. |
14     C | This implementation assumes A is a five-diagonal |
15     C | matrix of the form that arises in the discrete |
16     C | representation of the del^2 operator in a |
17     C | two-dimensional space. |
18     C | Notes: |
19     C | ====== |
20     C | This implementation can support shared-memory |
21     C | multi-threaded execution. In order to do this COMMON |
22     C | blocks are used for many of the arrays - even ones that |
23     C | are only used for intermedaite results. This design is |
24     C | OK if you want to all the threads to collaborate on |
25     C | solving the same problem. On the other hand if you want |
26     C | the threads to solve several different problems |
27     C | concurrently this implementation will not work. |
28     C \==========================================================/
29    
30     C === Global data ===
31     #include "SIZE.h"
32     #include "EEPARAMS.h"
33     #include "PARAMS.h"
34 cnh 1.4 #include "GRID.h"
35 cnh 1.1 #include "CG2D.h"
36    
37     C === Routine arguments ===
38     C myThid - Thread on which I am working.
39     INTEGER myThid
40    
41     C === Local variables ====
42     C actualIts - Number of iterations taken
43     C actualResidual - residual
44     C bi - Block index in X and Y.
45     C bj
46     C etaN - Used in computing search directions
47     C etaNM1 suffix N and NM1 denote current and
48     C cgBeta previous iterations respectively.
49     C alpha
50     C sumRHS - Sum of right-hand-side. Sometimes this is a
51     C useful debuggin/trouble shooting diagnostic.
52     C For neumann problems sumRHS needs to be ~0.
53     C or they converge at a non-zero residual.
54     C err - Measure of residual of Ax - b, usually the norm.
55     C I, J, N - Loop counters ( N counts CG iterations )
56     INTEGER actualIts
57     REAL actualResidual
58     INTEGER bi, bj
59     INTEGER I, J, it2d
60     REAL err
61     REAL etaN
62     REAL etaNM1
63     REAL cgBeta
64     REAL alpha
65     REAL sumRHS
66     REAL rhsMax
67     REAL rhsNorm
68    
69 cnh 1.13 INTEGER OLw
70     INTEGER OLe
71     INTEGER OLn
72     INTEGER OLs
73     INTEGER exchWidthX
74     INTEGER exchWidthY
75     INTEGER myNz
76    
77    
78 cnh 1.12 CcnhDebugStarts
79     CHARACTER*(MAX_LEN_FNAM) suff
80     CcnhDebugEnds
81    
82    
83 cnh 1.1 C-- Initialise inverter
84     etaNBuf(1,myThid) = 0. D0
85     errBuf(1,myThid) = 0. D0
86     sumRHSBuf(1,myThid) = 0. D0
87     etaNM1 = 1. D0
88    
89 cnh 1.10 CcnhDebugStarts
90 cnh 1.11 C _EXCH_XY_R8( cg2d_b, myThid )
91     C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid )
92 cnh 1.12 C suff = 'unnormalised'
93     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
94 cnh 1.10 CcnhDebugEnds
95    
96 cnh 1.1 C-- Normalise RHS
97     rhsMax = 0. _d 0
98     DO bj=myByLo(myThid),myByHi(myThid)
99     DO bi=myBxLo(myThid),myBxHi(myThid)
100     DO J=1,sNy
101     DO I=1,sNx
102     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm
103     rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax)
104     ENDDO
105     ENDDO
106     ENDDO
107     ENDDO
108     rhsMaxBuf(1,myThid) = rhsMax
109     _GLOBAL_MAX_R8( rhsMaxbuf, rhsMax, myThid )
110     rhsMax = rhsMaxBuf(1,1)
111     rhsNorm = 1. _d 0
112     IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
113     DO bj=myByLo(myThid),myByHi(myThid)
114     DO bi=myBxLo(myThid),myBxHi(myThid)
115     DO J=1,sNy
116     DO I=1,sNx
117     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
118     cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
119     ENDDO
120     ENDDO
121     ENDDO
122     ENDDO
123    
124     C-- Update overlaps
125     _EXCH_XY_R8( cg2d_b, myThid )
126     _EXCH_XY_R8( cg2d_x, myThid )
127     CcnhDebugStarts
128 cnh 1.11 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
129 cnh 1.12 C suff = 'normalised'
130     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
131 cnh 1.1 CcnhDebugEnds
132    
133     C-- Initial residual calculation
134     err = 0. _d 0
135     sumRHS = 0. _d 0
136     DO bj=myByLo(myThid),myByHi(myThid)
137     DO bi=myBxLo(myThid),myBxHi(myThid)
138     DO J=1,sNy
139     DO I=1,sNx
140     cg2d_s(I,J,bi,bj) = 0.
141     cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
142     & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
143     & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
144     & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
145     & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
146     & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
147     & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
148     & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
149 cnh 1.4 & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
150 cnh 1.10 & -freeSurfFac*_rA(i,j,bi,bj)* horiVertRatio*
151 cnh 1.4 & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
152     & )
153 cnh 1.1 err = err +
154     & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
155     sumRHS = sumRHS +
156     & cg2d_b(I,J,bi,bj)
157     ENDDO
158     ENDDO
159     ENDDO
160     ENDDO
161 cnh 1.13 C _EXCH_XY_R8( cg2d_r, myThid )
162     OLw = 1
163     OLe = 1
164     OLn = 1
165     OLs = 1
166     exchWidthX = 1
167     exchWidthY = 1
168     myNz = 1
169     CALL EXCH_RL( cg2d_r,
170     I OLw, OLe, OLs, OLn, myNz,
171     I exchWidthX, exchWidthY,
172     I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
173     C _EXCH_XY_R8( cg2d_s, myThid )
174     OLw = 1
175     OLe = 1
176     OLn = 1
177     OLs = 1
178     exchWidthX = 1
179     exchWidthY = 1
180     myNz = 1
181     CALL EXCH_RL( cg2d_s,
182     I OLw, OLe, OLs, OLn, myNz,
183     I exchWidthX, exchWidthY,
184     I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
185 cnh 1.1 sumRHSBuf(1,myThid) = sumRHS
186     _GLOBAL_SUM_R8( sumRHSBuf , sumRHS, myThid )
187     sumRHS = sumRHSBuf(1,1)
188     errBuf(1,myThid) = err
189     C WRITE(6,*) ' mythid, err = ', mythid, SQRT(err)
190     _GLOBAL_SUM_R8( errBuf , err , myThid )
191     err = errBuf(1,1)
192 cnh 1.13
193     _BEGIN_MASTER( myThid )
194 cnh 1.9 write(0,*) 'cg2d: Sum(rhs) = ',sumRHS
195 cnh 1.13 _END_MASTER( )
196 cnh 1.1
197     actualIts = 0
198     actualResidual = SQRT(err)
199     C _BARRIER
200     _BEGIN_MASTER( myThid )
201 cnh 1.6 WRITE(0,*) ' CG2D iters, err = ', actualIts, actualResidual
202 cnh 1.1 _END_MASTER( )
203    
204     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
205     DO 10 it2d=1, cg2dMaxIters
206    
207     CcnhDebugStarts
208 cnh 1.10 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' residual = ',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     etanBuf(1,myThid) = etaN
235     _GLOBAL_SUM_R8(etaNbuf,etaN, myThid)
236     etaN = etaNBuf(1,1)
237     CcnhDebugStarts
238     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' etaN = ',etaN
239     CcnhDebugEnds
240     cgBeta = etaN/etaNM1
241     CcnhDebugStarts
242     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
243     CcnhDebugEnds
244     etaNM1 = etaN
245    
246     DO bj=myByLo(myThid),myByHi(myThid)
247     DO bi=myBxLo(myThid),myBxHi(myThid)
248     DO J=1,sNy
249     DO I=1,sNx
250     cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj) + cgBeta*cg2d_s(I,J,bi,bj)
251     ENDDO
252     ENDDO
253     ENDDO
254     ENDDO
255    
256     C-- Do exchanges that require messages i.e. between
257     C-- processes.
258 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
259     OLw = 1
260     OLe = 1
261     OLn = 1
262     OLs = 1
263     exchWidthX = 1
264     exchWidthY = 1
265     myNz = 1
266     CALL EXCH_RL( cg2d_s,
267     I OLw, OLe, OLs, OLn, myNz,
268     I exchWidthX, exchWidthY,
269     I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
270    
271 cnh 1.1
272     C== Evaluate laplace operator on conjugate gradient vector
273     C== q = A.s
274     alpha = 0. _d 0
275     DO bj=myByLo(myThid),myByHi(myThid)
276     DO bi=myBxLo(myThid),myBxHi(myThid)
277     DO J=1,sNy
278     DO I=1,sNx
279     cg2d_q(I,J,bi,bj) =
280     & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
281     & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
282     & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
283     & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
284     & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
285     & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
286     & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
287     & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
288 cnh 1.10 & -freeSurfFac*_rA(i,j,bi,bj)* horiVertRatio*
289 cnh 1.4 & cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
290 cnh 1.1 alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
291     ENDDO
292     ENDDO
293     ENDDO
294     ENDDO
295     alphaBuf(1,myThid) = alpha
296     _GLOBAL_SUM_R8(alphaBuf,alpha,myThid)
297     alpha = alphaBuf(1,1)
298     CcnhDebugStarts
299     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
300     CcnhDebugEnds
301     alpha = etaN/alpha
302     CcnhDebugStarts
303     C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
304     CcnhDebugEnds
305    
306     C== Update solution and residual vectors
307     C Now compute "interior" points.
308     err = 0. _d 0
309     DO bj=myByLo(myThid),myByHi(myThid)
310     DO bi=myBxLo(myThid),myBxHi(myThid)
311     DO J=1,sNy
312     DO I=1,sNx
313     cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
314     cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
315     err = err+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
316     ENDDO
317     ENDDO
318     ENDDO
319     ENDDO
320    
321     errBuf(1,myThid) = err
322     _GLOBAL_SUM_R8( errBuf , err , myThid )
323     err = errBuf(1,1)
324     err = SQRT(err)
325     actualIts = it2d
326     actualResidual = err
327     IF ( err .LT. cg2dTargetResidual ) GOTO 11
328 cnh 1.13 C _EXCH_XY_R8(cg2d_r, myThid )
329     OLw = 1
330     OLe = 1
331     OLn = 1
332     OLs = 1
333     exchWidthX = 1
334     exchWidthY = 1
335     myNz = 1
336     CALL EXCH_RL( cg2d_r,
337     I OLw, OLe, OLs, OLn, myNz,
338     I exchWidthX, exchWidthY,
339     I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
340    
341 cnh 1.1 10 CONTINUE
342     11 CONTINUE
343    
344     C-- Un-normalise the answer
345     DO bj=myByLo(myThid),myByHi(myThid)
346     DO bi=myBxLo(myThid),myBxHi(myThid)
347     DO J=1,sNy
348     DO I=1,sNx
349     cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
350     ENDDO
351     ENDDO
352     ENDDO
353     ENDDO
354    
355     _EXCH_XY_R8(cg2d_x, myThid )
356 cnh 1.6 _BEGIN_MASTER( myThid )
357     WRITE(0,*) ' CG2D iters, err = ', actualIts, actualResidual
358     _END_MASTER( )
359 cnh 1.1
360     CcnhDebugStarts
361 cnh 1.7 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
362 cnh 1.1 C err = 0. _d 0
363     C DO bj=myByLo(myThid),myByHi(myThid)
364     C DO bi=myBxLo(myThid),myBxHi(myThid)
365     C DO J=1,sNy
366     C DO I=1,sNx
367     C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
368     C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
369     C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
370     C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
371     C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
372     C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
373     C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
374     C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
375     C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
376     C err = err +
377     C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
378     C ENDDO
379     C ENDDO
380     C ENDDO
381     C ENDDO
382     C errBuf(1,myThid) = err
383     C _GLOBAL_SUM_R8( errBuf , err , myThid )
384     C err = errBuf(1,1)
385     C write(0,*) 'cg2d: Ax - b = ',SQRT(err)
386     CcnhDebugEnds
387    
388    
389     END

  ViewVC Help
Powered by ViewVC 1.1.22