/[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.51 - (hide annotations) (download)
Tue Apr 28 18:01:14 2009 UTC (15 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61n, checkpoint61m
Changes since 1.50: +14 -14 lines
change macros (EXCH & GLOBAL_SUM/MAX) sufix _R4/_R8 to _RS/_RL
 when applied to _RS/_RL variable

1 jmc 1.51 C $Header: /u/gcmpack/MITgcm/model/src/cg2d.F,v 1.50 2008/01/08 23:59:39 jahn Exp $
2 adcroft 1.33 C $Name: $
3 cnh 1.1
4 cnh 1.16 #include "CPP_OPTIONS.h"
5 mlosch 1.49 #ifdef TARGET_NEC_SX
6     C set a sensible default for the outer loop unrolling parameter that can
7     C be overriden in the Makefile with the DEFINES macro or in CPP_OPTIONS.h
8     #ifndef CG2D_OUTERLOOPITERS
9     # define CG2D_OUTERLOOPITERS 10
10     #endif
11     #endif /* TARGET_NEC_SX */
12 cnh 1.1
13 cnh 1.34 CBOP
14     C !ROUTINE: CG2D
15     C !INTERFACE:
16 jmc 1.45 SUBROUTINE CG2D(
17 cnh 1.14 I cg2d_b,
18     U cg2d_x,
19 adcroft 1.33 O firstResidual,
20     O lastResidual,
21 adcroft 1.30 U numIters,
22 cnh 1.1 I myThid )
23 cnh 1.34 C !DESCRIPTION: \bv
24     C *==========================================================*
25 jmc 1.45 C | SUBROUTINE CG2D
26     C | o Two-dimensional grid problem conjugate-gradient
27     C | inverter (with preconditioner).
28 cnh 1.34 C *==========================================================*
29 jmc 1.45 C | Con. grad is an iterative procedure for solving Ax = b.
30     C | It requires the A be symmetric.
31     C | This implementation assumes A is a five-diagonal
32     C | matrix of the form that arises in the discrete
33     C | representation of the del^2 operator in a
34     C | two-dimensional space.
35     C | Notes:
36     C | ======
37     C | This implementation can support shared-memory
38     C | multi-threaded execution. In order to do this COMMON
39     C | blocks are used for many of the arrays - even ones that
40     C | are only used for intermedaite results. This design is
41     C | OK if you want to all the threads to collaborate on
42     C | solving the same problem. On the other hand if you want
43     C | the threads to solve several different problems
44     C | concurrently this implementation will not work.
45 cnh 1.34 C *==========================================================*
46     C \ev
47    
48     C !USES:
49 adcroft 1.18 IMPLICIT NONE
50 cnh 1.1 C === Global data ===
51     #include "SIZE.h"
52     #include "EEPARAMS.h"
53     #include "PARAMS.h"
54 jmc 1.45 #include "CG2D.h"
55 jmc 1.46 c#include "GRID.h"
56     c#include "SURFACE.h"
57 cnh 1.1
58 cnh 1.34 C !INPUT/OUTPUT PARAMETERS:
59 cnh 1.1 C === Routine arguments ===
60 jmc 1.45 C myThid :: Thread on which I am working.
61     C cg2d_b :: The source term or "right hand side"
62     C cg2d_x :: The solution
63     C firstResidual :: the initial residual before any iterations
64     C lastResidual :: the actual residual reached
65     C numIters :: Entry: the maximum number of iterations allowed
66     C Exit: the actual number of iterations used
67 adcroft 1.30 _RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
68     _RL cg2d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
69 adcroft 1.33 _RL firstResidual
70     _RL lastResidual
71 adcroft 1.30 INTEGER numIters
72 cnh 1.1 INTEGER myThid
73    
74 cnh 1.34 C !LOCAL VARIABLES:
75 cnh 1.1 C === Local variables ====
76 jmc 1.45 C actualIts :: Number of iterations taken
77     C actualResidual :: residual
78     C bi, bj :: Block index in X and Y.
79     C eta_qrN :: Used in computing search directions
80 jmc 1.28 C eta_qrNM1 suffix N and NM1 denote current and
81 cnh 1.1 C cgBeta previous iterations respectively.
82 jmc 1.45 C alpha
83     C sumRHS :: Sum of right-hand-side. Sometimes this is a
84 cnh 1.1 C useful debuggin/trouble shooting diagnostic.
85     C For neumann problems sumRHS needs to be ~0.
86     C or they converge at a non-zero residual.
87 jmc 1.45 C err :: Measure of residual of Ax - b, usually the norm.
88     C I, J, it2d :: Loop counters ( it2d counts CG iterations )
89 cnh 1.1 INTEGER actualIts
90 adcroft 1.33 _RL actualResidual
91 jmc 1.45 INTEGER bi, bj
92 cnh 1.1 INTEGER I, J, it2d
93 jmc 1.46 c INTEGER ks
94 jmc 1.48 _RL err, errTile(nSx,nSy)
95     _RL eta_qrN,eta_qrNtile(nSx,nSy)
96 jmc 1.28 _RL eta_qrNM1
97 cnh 1.14 _RL cgBeta
98 jmc 1.48 _RL alpha, alphaTile(nSx,nSy)
99     _RL sumRHS, sumRHStile(nSx,nSy)
100 cnh 1.14 _RL rhsMax
101     _RL rhsNorm
102 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
103     COMMON /CG2D_SINGLECPU_BUF/ localBuf
104     _RL localBuf(1:sNx,1:sNy,nSx,nSy)
105     #endif
106 jmc 1.47 CHARACTER*(MAX_LEN_MBUF) msgBuf
107 cnh 1.34 CEOP
108 cnh 1.13
109    
110 cnh 1.12 CcnhDebugStarts
111 adcroft 1.24 C CHARACTER*(MAX_LEN_FNAM) suff
112 cnh 1.12 CcnhDebugEnds
113    
114    
115 cnh 1.1 C-- Initialise inverter
116 jmc 1.28 eta_qrNM1 = 1. _d 0
117 cnh 1.1
118 cnh 1.10 CcnhDebugStarts
119 jmc 1.51 C _EXCH_XY_RL( cg2d_b, myThid )
120 cnh 1.11 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid )
121 cnh 1.12 C suff = 'unnormalised'
122     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
123 cnh 1.14 C STOP
124 cnh 1.10 CcnhDebugEnds
125    
126 cnh 1.1 C-- Normalise RHS
127     rhsMax = 0. _d 0
128     DO bj=myByLo(myThid),myByHi(myThid)
129     DO bi=myBxLo(myThid),myBxHi(myThid)
130     DO J=1,sNy
131     DO I=1,sNx
132     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm
133     rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax)
134     ENDDO
135     ENDDO
136     ENDDO
137     ENDDO
138 adcroft 1.33
139     IF (cg2dNormaliseRHS) THEN
140     C- Normalise RHS :
141 adcroft 1.23 #ifdef LETS_MAKE_JAM
142 jmc 1.51 C _GLOBAL_MAX_RL( rhsMax, myThid )
143 adcroft 1.25 rhsMax=1.
144 adcroft 1.23 #else
145 jmc 1.51 _GLOBAL_MAX_RL( rhsMax, myThid )
146 adcroft 1.23 #endif
147 cnh 1.1 rhsNorm = 1. _d 0
148     IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
149     DO bj=myByLo(myThid),myByHi(myThid)
150     DO bi=myBxLo(myThid),myBxHi(myThid)
151     DO J=1,sNy
152     DO I=1,sNx
153     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
154     cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
155     ENDDO
156     ENDDO
157     ENDDO
158     ENDDO
159 adcroft 1.33 C- end Normalise RHS
160     ENDIF
161 cnh 1.1
162     C-- Update overlaps
163 jmc 1.51 _EXCH_XY_RL( cg2d_b, myThid )
164     _EXCH_XY_RL( cg2d_x, myThid )
165 cnh 1.1 CcnhDebugStarts
166 cnh 1.11 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
167 cnh 1.12 C suff = 'normalised'
168     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
169 cnh 1.1 CcnhDebugEnds
170    
171     C-- Initial residual calculation
172     err = 0. _d 0
173     sumRHS = 0. _d 0
174     DO bj=myByLo(myThid),myByHi(myThid)
175     DO bi=myBxLo(myThid),myBxHi(myThid)
176 jmc 1.48 sumRHStile(bi,bj) = 0. _d 0
177     errTile(bi,bj) = 0. _d 0
178 mlosch 1.49 #ifdef TARGET_NEC_SX
179     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
180     #endif /* TARGET_NEC_SX */
181 cnh 1.1 DO J=1,sNy
182     DO I=1,sNx
183 jmc 1.46 c ks = ksurfC(I,J,bi,bj)
184 cnh 1.1 cg2d_s(I,J,bi,bj) = 0.
185     cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
186     & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
187     & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
188     & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
189     & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
190 jmc 1.46 & +aC2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
191 cnh 1.4 & )
192 jmc 1.46 c & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
193     c & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
194     c & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
195     c & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
196     c & -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)*recip_Bo(i,j,bi,bj)
197     c & *cg2d_x(I ,J ,bi,bj)
198     c & /deltaTMom/deltaTfreesurf*cg2dNorm
199     c & )
200 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
201     localBuf(I,J,bi,bj) = cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
202     #else
203 jmc 1.48 errTile(bi,bj) = errTile(bi,bj)
204     & + cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
205     sumRHStile(bi,bj) = sumRHStile(bi,bj) + cg2d_b(I,J,bi,bj)
206 jahn 1.50 #endif
207 cnh 1.1 ENDDO
208     ENDDO
209 jmc 1.48 c sumRHS = sumRHS + sumRHStile(bi,bj)
210     c err = err + errTile(bi,bj)
211 cnh 1.1 ENDDO
212     ENDDO
213 adcroft 1.23 #ifdef LETS_MAKE_JAM
214     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
215     #else
216 heimbach 1.35 CALL EXCH_XY_RL( cg2d_r, myThid )
217 adcroft 1.23 #endif
218     #ifdef LETS_MAKE_JAM
219     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
220     #else
221 heimbach 1.35 CALL EXCH_XY_RL( cg2d_s, myThid )
222 adcroft 1.23 #endif
223 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
224     CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, err, myThid)
225     DO bj=myByLo(myThid),myByHi(myThid)
226     DO bi=myBxLo(myThid),myBxHi(myThid)
227     DO J=1,sNy
228     DO I=1,sNx
229     localBuf(I,J,bi,bj) = cg2d_b(I,J,bi,bj)
230     ENDDO
231     ENDDO
232     ENDDO
233     ENDDO
234     CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, sumRHS, myThid)
235     #else
236 jmc 1.51 c _GLOBAL_SUM_RL( sumRHS, myThid )
237     c _GLOBAL_SUM_RL( err , myThid )
238 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( sumRHStile, sumRHS, myThid )
239     CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
240 jahn 1.50 #endif
241 adcroft 1.33 err = SQRT(err)
242     actualIts = 0
243     actualResidual = err
244 dimitri 1.40
245 jmc 1.39 IF ( debugLevel .GE. debLevZero ) THEN
246 heimbach 1.37 _BEGIN_MASTER( myThid )
247 jmc 1.45 WRITE(standardmessageunit,'(A,1P2E22.14)')
248     & ' cg2d: Sum(rhs),rhsMax = ', sumRHS,rhsMax
249 edhill 1.41 _END_MASTER( myThid )
250 heimbach 1.37 ENDIF
251 cnh 1.1 C _BARRIER
252 adcroft 1.30 c _BEGIN_MASTER( myThid )
253 jmc 1.45 c WRITE(standardmessageunit,'(A,I6,1PE30.14)')
254     c & ' CG2D iters, err = ',
255 adcroft 1.30 c & actualIts, actualResidual
256 edhill 1.41 c _END_MASTER( myThid )
257 adcroft 1.33 firstResidual=actualResidual
258    
259     IF ( err .LT. cg2dTolerance ) GOTO 11
260 cnh 1.1
261     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
262 adcroft 1.30 DO 10 it2d=1, numIters
263 cnh 1.1
264     CcnhDebugStarts
265 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' residual = ',
266 cnh 1.14 C & actualResidual
267 cnh 1.1 CcnhDebugEnds
268     C-- Solve preconditioning equation and update
269     C-- conjugate direction vector "s".
270 jmc 1.28 eta_qrN = 0. _d 0
271 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
272     DO bi=myBxLo(myThid),myBxHi(myThid)
273 jmc 1.48 eta_qrNtile(bi,bj) = 0. _d 0
274 mlosch 1.49 #ifdef TARGET_NEC_SX
275     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
276     #endif /* TARGET_NEC_SX */
277 cnh 1.1 DO J=1,sNy
278     DO I=1,sNx
279 jmc 1.45 cg2d_q(I,J,bi,bj) =
280 cnh 1.3 & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
281     & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
282     & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
283     & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
284     & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
285 cnh 1.4 CcnhDebugStarts
286     C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
287     CcnhDebugEnds
288 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
289     localBuf(I,J,bi,bj) =
290     & cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
291     #else
292 jmc 1.48 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
293 cnh 1.1 & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
294 jahn 1.50 #endif
295 cnh 1.1 ENDDO
296     ENDDO
297 jmc 1.48 c eta_qrN = eta_qrN + eta_qrNtile(bi,bj)
298 cnh 1.1 ENDDO
299     ENDDO
300    
301 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
302     CALL GLOBAL_SUM_SINGLECPU_RL( localBuf,eta_qrN,myThid )
303     #else
304 jmc 1.51 c _GLOBAL_SUM_RL(eta_qrN, myThid)
305 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( eta_qrNtile,eta_qrN,myThid )
306 jahn 1.50 #endif
307 cnh 1.1 CcnhDebugStarts
308 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
309 cnh 1.1 CcnhDebugEnds
310 jmc 1.28 cgBeta = eta_qrN/eta_qrNM1
311 cnh 1.1 CcnhDebugStarts
312 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
313 cnh 1.1 CcnhDebugEnds
314 jmc 1.28 eta_qrNM1 = eta_qrN
315 cnh 1.1
316     DO bj=myByLo(myThid),myByHi(myThid)
317     DO bi=myBxLo(myThid),myBxHi(myThid)
318     DO J=1,sNy
319     DO I=1,sNx
320 cnh 1.14 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj)
321     & + cgBeta*cg2d_s(I,J,bi,bj)
322 cnh 1.1 ENDDO
323     ENDDO
324     ENDDO
325     ENDDO
326    
327     C-- Do exchanges that require messages i.e. between
328     C-- processes.
329 jmc 1.51 C _EXCH_XY_RL( cg2d_s, myThid )
330 adcroft 1.23 #ifdef LETS_MAKE_JAM
331     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
332     #else
333 heimbach 1.35 CALL EXCH_XY_RL( cg2d_s, myThid )
334 adcroft 1.23 #endif
335 cnh 1.1
336     C== Evaluate laplace operator on conjugate gradient vector
337     C== q = A.s
338     alpha = 0. _d 0
339     DO bj=myByLo(myThid),myByHi(myThid)
340     DO bi=myBxLo(myThid),myBxHi(myThid)
341 jmc 1.48 alphaTile(bi,bj) = 0. _d 0
342 mlosch 1.49 #ifdef TARGET_NEC_SX
343     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
344     #endif /* TARGET_NEC_SX */
345 cnh 1.1 DO J=1,sNy
346     DO I=1,sNx
347 jmc 1.46 c ks = ksurfC(I,J,bi,bj)
348 jmc 1.45 cg2d_q(I,J,bi,bj) =
349 cnh 1.1 & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
350     & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
351     & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
352     & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
353 jmc 1.46 & +aC2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
354     c & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
355     c & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
356     c & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
357     c & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
358     c & -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)*recip_Bo(i,j,bi,bj)
359     c & *cg2d_s(I ,J ,bi,bj)
360     c & /deltaTMom/deltaTfreesurf*cg2dNorm
361 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
362     localBuf(I,J,bi,bj) = cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
363     #else
364 jmc 1.48 alphaTile(bi,bj) = alphaTile(bi,bj)
365     & + cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
366 jahn 1.50 #endif
367 cnh 1.1 ENDDO
368     ENDDO
369 jmc 1.48 c alpha = alpha + alphaTile(bi,bj)
370 cnh 1.1 ENDDO
371     ENDDO
372 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
373     CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, alpha, myThid)
374     #else
375 jmc 1.51 c _GLOBAL_SUM_RL(alpha,myThid)
376 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( alphaTile, alpha, myThid )
377 jahn 1.50 #endif
378 cnh 1.1 CcnhDebugStarts
379 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
380 cnh 1.1 CcnhDebugEnds
381 jmc 1.28 alpha = eta_qrN/alpha
382 cnh 1.1 CcnhDebugStarts
383 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
384 cnh 1.1 CcnhDebugEnds
385 jmc 1.45
386 cnh 1.1 C== Update solution and residual vectors
387     C Now compute "interior" points.
388     err = 0. _d 0
389     DO bj=myByLo(myThid),myByHi(myThid)
390     DO bi=myBxLo(myThid),myBxHi(myThid)
391 jmc 1.48 errTile(bi,bj) = 0. _d 0
392 cnh 1.1 DO J=1,sNy
393     DO I=1,sNx
394     cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
395     cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
396 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
397     localBuf(I,J,bi,bj) = cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
398     #else
399 jmc 1.48 errTile(bi,bj) = errTile(bi,bj)
400     & + cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
401 jahn 1.50 #endif
402 cnh 1.1 ENDDO
403     ENDDO
404 jmc 1.48 c err = err + errTile(bi,bj)
405 cnh 1.1 ENDDO
406     ENDDO
407    
408 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
409     CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, err, myThid)
410     #else
411 jmc 1.51 c _GLOBAL_SUM_RL( err , myThid )
412 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
413 jahn 1.50 #endif
414 cnh 1.1 err = SQRT(err)
415     actualIts = it2d
416     actualResidual = err
417 jmc 1.47 IF ( debugLevel.GT.debLevB ) THEN
418     c IF ( DIFFERENT_MULTIPLE(monitorFreq,myTime,deltaTClock)
419     c & ) THEN
420     _BEGIN_MASTER( myThid )
421     WRITE(msgBuf,'(A,I6,A,1PE21.14)')
422     & ' cg2d: iter=', actualIts, ' ; resid.= ', actualResidual
423     CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
424     _END_MASTER( myThid )
425     c ENDIF
426     ENDIF
427 adcroft 1.33 IF ( err .LT. cg2dTolerance ) GOTO 11
428 jmc 1.47
429 adcroft 1.23 #ifdef LETS_MAKE_JAM
430     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
431     #else
432 heimbach 1.35 CALL EXCH_XY_RL( cg2d_r, myThid )
433 adcroft 1.23 #endif
434 cnh 1.13
435 cnh 1.1 10 CONTINUE
436     11 CONTINUE
437    
438 adcroft 1.33 IF (cg2dNormaliseRHS) THEN
439 cnh 1.1 C-- Un-normalise the answer
440 adcroft 1.33 DO bj=myByLo(myThid),myByHi(myThid)
441     DO bi=myBxLo(myThid),myBxHi(myThid)
442     DO J=1,sNy
443     DO I=1,sNx
444     cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
445     ENDDO
446     ENDDO
447 cnh 1.1 ENDDO
448     ENDDO
449 adcroft 1.33 ENDIF
450 cnh 1.1
451 adcroft 1.22 C The following exchange was moved up to solve_for_pressure
452     C for compatibility with TAMC.
453 jmc 1.51 C _EXCH_XY_RL(cg2d_x, myThid )
454 adcroft 1.30 c _BEGIN_MASTER( myThid )
455 jmc 1.45 c WRITE(*,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
456 adcroft 1.30 c & actualIts, actualResidual
457 edhill 1.41 c _END_MASTER( myThid )
458 adcroft 1.30
459     C-- Return parameters to caller
460 adcroft 1.33 lastResidual=actualResidual
461 adcroft 1.30 numIters=actualIts
462 cnh 1.1
463     CcnhDebugStarts
464 cnh 1.7 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
465 cnh 1.1 C err = 0. _d 0
466     C DO bj=myByLo(myThid),myByHi(myThid)
467     C DO bi=myBxLo(myThid),myBxHi(myThid)
468     C DO J=1,sNy
469     C DO I=1,sNx
470     C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
471     C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
472     C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
473     C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
474     C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
475     C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
476     C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
477     C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
478     C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
479 jmc 1.45 C err = err +
480 cnh 1.1 C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
481     C ENDDO
482     C ENDDO
483     C ENDDO
484     C ENDDO
485 jmc 1.51 C _GLOBAL_SUM_RL( err , myThid )
486 heimbach 1.31 C write(*,*) 'cg2d: Ax - b = ',SQRT(err)
487 cnh 1.1 CcnhDebugEnds
488    
489 adcroft 1.19 RETURN
490 cnh 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22