/[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.54 - (hide annotations) (download)
Wed Jun 8 01:46:34 2011 UTC (13 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62z, checkpoint63g, checkpoint63, checkpoint63l, checkpoint63m, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.53: +15 -29 lines
use new parameter "printResidualFreq" to print residual in CG iterations

1 jmc 1.54 C $Header: /u/gcmpack/MITgcm/model/src/cg2d.F,v 1.53 2009/07/11 22:00:40 jmc 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 cg2d_b :: The source term or "right hand side"
61     C cg2d_x :: The solution
62     C firstResidual :: the initial residual before any iterations
63     C lastResidual :: the actual residual reached
64     C numIters :: Entry: the maximum number of iterations allowed
65     C Exit: the actual number of iterations used
66 jmc 1.53 C myThid :: Thread on which I am working.
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     _RL localBuf(1:sNx,1:sNy,nSx,nSy)
104     #endif
105 jmc 1.47 CHARACTER*(MAX_LEN_MBUF) msgBuf
106 jmc 1.54 LOGICAL printResidual
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.53 c CALL EXCH_XY_RL( cg2d_b, myThid )
164     CALL 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     DO bj=myByLo(myThid),myByHi(myThid)
173     DO bi=myBxLo(myThid),myBxHi(myThid)
174 jmc 1.53 DO J=1-1,sNy+1
175     DO I=1-1,sNx+1
176     cg2d_s(I,J,bi,bj) = 0.
177     ENDDO
178     ENDDO
179 jmc 1.48 sumRHStile(bi,bj) = 0. _d 0
180     errTile(bi,bj) = 0. _d 0
181 mlosch 1.49 #ifdef TARGET_NEC_SX
182     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
183     #endif /* TARGET_NEC_SX */
184 cnh 1.1 DO J=1,sNy
185     DO I=1,sNx
186 jmc 1.54 c ks = kSurfC(I,J,bi,bj)
187 cnh 1.1 cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
188     & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
189     & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
190     & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
191     & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
192 jmc 1.46 & +aC2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
193 cnh 1.4 & )
194 jmc 1.46 c & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
195     c & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
196     c & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
197     c & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
198     c & -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)*recip_Bo(i,j,bi,bj)
199     c & *cg2d_x(I ,J ,bi,bj)
200     c & /deltaTMom/deltaTfreesurf*cg2dNorm
201     c & )
202 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
203     localBuf(I,J,bi,bj) = cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
204     #else
205 jmc 1.48 errTile(bi,bj) = errTile(bi,bj)
206     & + cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
207     sumRHStile(bi,bj) = sumRHStile(bi,bj) + cg2d_b(I,J,bi,bj)
208 jahn 1.50 #endif
209 cnh 1.1 ENDDO
210     ENDDO
211     ENDDO
212     ENDDO
213 jmc 1.53 CALL EXCH_S3D_RL( cg2d_r, 1, myThid )
214 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
215 jmc 1.52 CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, err, 0, 0, myThid)
216     CALL GLOBAL_SUM_SINGLECPU_RL(cg2d_b, sumRHS, OLx, OLy, myThid)
217 jahn 1.50 #else
218 jmc 1.53 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
219     CALL GLOBAL_SUM_TILE_RL( sumRHStile, sumRHS, myThid )
220 jahn 1.50 #endif
221 adcroft 1.33 err = SQRT(err)
222     actualIts = 0
223     actualResidual = err
224 dimitri 1.40
225 jmc 1.54 printResidual = .FALSE.
226 jmc 1.53 IF ( debugLevel .GE. debLevZero ) THEN
227 heimbach 1.37 _BEGIN_MASTER( myThid )
228 jmc 1.54 printResidual = printResidualFreq.GE.1
229 jmc 1.45 WRITE(standardmessageunit,'(A,1P2E22.14)')
230     & ' cg2d: Sum(rhs),rhsMax = ', sumRHS,rhsMax
231 edhill 1.41 _END_MASTER( myThid )
232 jmc 1.53 ENDIF
233 adcroft 1.33 firstResidual=actualResidual
234    
235     IF ( err .LT. cg2dTolerance ) GOTO 11
236 cnh 1.1
237     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
238 adcroft 1.30 DO 10 it2d=1, numIters
239 cnh 1.1
240     C-- Solve preconditioning equation and update
241     C-- conjugate direction vector "s".
242     DO bj=myByLo(myThid),myByHi(myThid)
243     DO bi=myBxLo(myThid),myBxHi(myThid)
244 jmc 1.48 eta_qrNtile(bi,bj) = 0. _d 0
245 mlosch 1.49 #ifdef TARGET_NEC_SX
246     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
247     #endif /* TARGET_NEC_SX */
248 cnh 1.1 DO J=1,sNy
249     DO I=1,sNx
250 jmc 1.45 cg2d_q(I,J,bi,bj) =
251 cnh 1.3 & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
252     & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
253     & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
254     & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
255     & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
256 cnh 1.4 CcnhDebugStarts
257     C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
258     CcnhDebugEnds
259 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
260     localBuf(I,J,bi,bj) =
261     & cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
262     #else
263 jmc 1.48 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
264 cnh 1.1 & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
265 jahn 1.50 #endif
266 cnh 1.1 ENDDO
267     ENDDO
268     ENDDO
269     ENDDO
270    
271 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
272 jmc 1.52 CALL GLOBAL_SUM_SINGLECPU_RL( localBuf,eta_qrN,0,0,myThid )
273 jahn 1.50 #else
274 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( eta_qrNtile,eta_qrN,myThid )
275 jahn 1.50 #endif
276 cnh 1.1 CcnhDebugStarts
277 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
278 cnh 1.1 CcnhDebugEnds
279 jmc 1.28 cgBeta = eta_qrN/eta_qrNM1
280 cnh 1.1 CcnhDebugStarts
281 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
282 cnh 1.1 CcnhDebugEnds
283 jmc 1.28 eta_qrNM1 = eta_qrN
284 cnh 1.1
285     DO bj=myByLo(myThid),myByHi(myThid)
286     DO bi=myBxLo(myThid),myBxHi(myThid)
287     DO J=1,sNy
288     DO I=1,sNx
289 cnh 1.14 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj)
290     & + cgBeta*cg2d_s(I,J,bi,bj)
291 cnh 1.1 ENDDO
292     ENDDO
293     ENDDO
294     ENDDO
295    
296 jmc 1.54 C-- Do exchanges that require messages i.e. between processes.
297 jmc 1.53 CALL EXCH_S3D_RL( cg2d_s, 1, myThid )
298 cnh 1.1
299     C== Evaluate laplace operator on conjugate gradient vector
300     C== q = A.s
301     DO bj=myByLo(myThid),myByHi(myThid)
302     DO bi=myBxLo(myThid),myBxHi(myThid)
303 jmc 1.48 alphaTile(bi,bj) = 0. _d 0
304 mlosch 1.49 #ifdef TARGET_NEC_SX
305     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
306     #endif /* TARGET_NEC_SX */
307 cnh 1.1 DO J=1,sNy
308     DO I=1,sNx
309 jmc 1.54 c ks = kSurfC(I,J,bi,bj)
310 jmc 1.45 cg2d_q(I,J,bi,bj) =
311 cnh 1.1 & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
312     & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
313     & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
314     & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
315 jmc 1.46 & +aC2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
316     c & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
317     c & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
318     c & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
319     c & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
320     c & -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)*recip_Bo(i,j,bi,bj)
321     c & *cg2d_s(I ,J ,bi,bj)
322     c & /deltaTMom/deltaTfreesurf*cg2dNorm
323 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
324     localBuf(I,J,bi,bj) = cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
325     #else
326 jmc 1.48 alphaTile(bi,bj) = alphaTile(bi,bj)
327     & + cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
328 jahn 1.50 #endif
329 cnh 1.1 ENDDO
330     ENDDO
331     ENDDO
332     ENDDO
333 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
334 jmc 1.52 CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, alpha, 0, 0, myThid)
335 jahn 1.50 #else
336 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( alphaTile, alpha, myThid )
337 jahn 1.50 #endif
338 cnh 1.1 CcnhDebugStarts
339 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
340 cnh 1.1 CcnhDebugEnds
341 jmc 1.28 alpha = eta_qrN/alpha
342 cnh 1.1 CcnhDebugStarts
343 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
344 cnh 1.1 CcnhDebugEnds
345 jmc 1.45
346 cnh 1.1 C== Update solution and residual vectors
347     C Now compute "interior" points.
348     DO bj=myByLo(myThid),myByHi(myThid)
349     DO bi=myBxLo(myThid),myBxHi(myThid)
350 jmc 1.48 errTile(bi,bj) = 0. _d 0
351 cnh 1.1 DO J=1,sNy
352     DO I=1,sNx
353     cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
354     cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
355 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
356     localBuf(I,J,bi,bj) = cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
357     #else
358 jmc 1.48 errTile(bi,bj) = errTile(bi,bj)
359     & + cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
360 jahn 1.50 #endif
361 cnh 1.1 ENDDO
362     ENDDO
363     ENDDO
364     ENDDO
365    
366 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
367 jmc 1.52 CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, err, 0, 0, myThid)
368 jahn 1.50 #else
369 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
370 jahn 1.50 #endif
371 cnh 1.1 err = SQRT(err)
372     actualIts = it2d
373     actualResidual = err
374 jmc 1.54 IF ( printResidual ) THEN
375     IF ( MOD( it2d-1, printResidualFreq ).EQ.0 ) THEN
376     WRITE(msgBuf,'(A,I6,A,1PE21.14)')
377 jmc 1.47 & ' cg2d: iter=', actualIts, ' ; resid.= ', actualResidual
378 jmc 1.54 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
379     & SQUEEZE_RIGHT, myThid )
380     ENDIF
381 jmc 1.47 ENDIF
382 adcroft 1.33 IF ( err .LT. cg2dTolerance ) GOTO 11
383 jmc 1.47
384 jmc 1.53 CALL EXCH_S3D_RL( cg2d_r, 1, myThid )
385 cnh 1.13
386 cnh 1.1 10 CONTINUE
387     11 CONTINUE
388    
389 adcroft 1.33 IF (cg2dNormaliseRHS) THEN
390 cnh 1.1 C-- Un-normalise the answer
391 adcroft 1.33 DO bj=myByLo(myThid),myByHi(myThid)
392     DO bi=myBxLo(myThid),myBxHi(myThid)
393     DO J=1,sNy
394     DO I=1,sNx
395     cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
396     ENDDO
397     ENDDO
398 cnh 1.1 ENDDO
399     ENDDO
400 adcroft 1.33 ENDIF
401 cnh 1.1
402 adcroft 1.22 C The following exchange was moved up to solve_for_pressure
403     C for compatibility with TAMC.
404 jmc 1.51 C _EXCH_XY_RL(cg2d_x, myThid )
405 adcroft 1.30
406     C-- Return parameters to caller
407 jmc 1.54 lastResidual = actualResidual
408     numIters = actualIts
409 cnh 1.1
410     CcnhDebugStarts
411 cnh 1.7 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
412 cnh 1.1 C err = 0. _d 0
413     C DO bj=myByLo(myThid),myByHi(myThid)
414     C DO bi=myBxLo(myThid),myBxHi(myThid)
415     C DO J=1,sNy
416     C DO I=1,sNx
417     C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
418     C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
419     C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
420     C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
421     C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
422     C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
423     C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
424     C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
425     C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
426 jmc 1.45 C err = err +
427 cnh 1.1 C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
428     C ENDDO
429     C ENDDO
430     C ENDDO
431     C ENDDO
432 jmc 1.51 C _GLOBAL_SUM_RL( err , myThid )
433 heimbach 1.31 C write(*,*) 'cg2d: Ax - b = ',SQRT(err)
434 cnh 1.1 CcnhDebugEnds
435    
436 adcroft 1.19 RETURN
437 cnh 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22