/[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.53 - (hide annotations) (download)
Sat Jul 11 22:00:40 2009 UTC (14 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62y, checkpoint62x, checkpoint62, checkpoint62b, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.52: +16 -30 lines
use simple EXCH (overlap 1 and no Corner Exch): this reduces number of
 EXCH calls by 2 if using exch2).

1 jmc 1.53 C $Header: /u/gcmpack/MITgcm/model/src/cg2d.F,v 1.52 2009/05/16 13:42:15 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 cnh 1.34 CEOP
107 cnh 1.13
108    
109 cnh 1.12 CcnhDebugStarts
110 adcroft 1.24 C CHARACTER*(MAX_LEN_FNAM) suff
111 cnh 1.12 CcnhDebugEnds
112    
113    
114 cnh 1.1 C-- Initialise inverter
115 jmc 1.28 eta_qrNM1 = 1. _d 0
116 cnh 1.1
117 cnh 1.10 CcnhDebugStarts
118 jmc 1.51 C _EXCH_XY_RL( cg2d_b, myThid )
119 cnh 1.11 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid )
120 cnh 1.12 C suff = 'unnormalised'
121     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
122 cnh 1.14 C STOP
123 cnh 1.10 CcnhDebugEnds
124    
125 cnh 1.1 C-- Normalise RHS
126     rhsMax = 0. _d 0
127     DO bj=myByLo(myThid),myByHi(myThid)
128     DO bi=myBxLo(myThid),myBxHi(myThid)
129     DO J=1,sNy
130     DO I=1,sNx
131     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm
132     rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax)
133     ENDDO
134     ENDDO
135     ENDDO
136     ENDDO
137 adcroft 1.33
138     IF (cg2dNormaliseRHS) THEN
139     C- Normalise RHS :
140 adcroft 1.23 #ifdef LETS_MAKE_JAM
141 jmc 1.51 C _GLOBAL_MAX_RL( rhsMax, myThid )
142 adcroft 1.25 rhsMax=1.
143 adcroft 1.23 #else
144 jmc 1.51 _GLOBAL_MAX_RL( rhsMax, myThid )
145 adcroft 1.23 #endif
146 cnh 1.1 rhsNorm = 1. _d 0
147     IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
148     DO bj=myByLo(myThid),myByHi(myThid)
149     DO bi=myBxLo(myThid),myBxHi(myThid)
150     DO J=1,sNy
151     DO I=1,sNx
152     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
153     cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
154     ENDDO
155     ENDDO
156     ENDDO
157     ENDDO
158 adcroft 1.33 C- end Normalise RHS
159     ENDIF
160 cnh 1.1
161     C-- Update overlaps
162 jmc 1.53 c CALL EXCH_XY_RL( cg2d_b, myThid )
163     CALL EXCH_XY_RL( cg2d_x, myThid )
164 cnh 1.1 CcnhDebugStarts
165 cnh 1.11 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
166 cnh 1.12 C suff = 'normalised'
167     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
168 cnh 1.1 CcnhDebugEnds
169    
170     C-- Initial residual calculation
171     DO bj=myByLo(myThid),myByHi(myThid)
172     DO bi=myBxLo(myThid),myBxHi(myThid)
173 jmc 1.53 DO J=1-1,sNy+1
174     DO I=1-1,sNx+1
175     cg2d_s(I,J,bi,bj) = 0.
176     ENDDO
177     ENDDO
178 jmc 1.48 sumRHStile(bi,bj) = 0. _d 0
179     errTile(bi,bj) = 0. _d 0
180 mlosch 1.49 #ifdef TARGET_NEC_SX
181     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
182     #endif /* TARGET_NEC_SX */
183 cnh 1.1 DO J=1,sNy
184     DO I=1,sNx
185 jmc 1.46 c ks = ksurfC(I,J,bi,bj)
186 cnh 1.1 cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
187     & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
188     & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
189     & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
190     & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
191 jmc 1.46 & +aC2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
192 cnh 1.4 & )
193 jmc 1.46 c & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
194     c & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
195     c & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
196     c & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
197     c & -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)*recip_Bo(i,j,bi,bj)
198     c & *cg2d_x(I ,J ,bi,bj)
199     c & /deltaTMom/deltaTfreesurf*cg2dNorm
200     c & )
201 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
202     localBuf(I,J,bi,bj) = cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
203     #else
204 jmc 1.48 errTile(bi,bj) = errTile(bi,bj)
205     & + cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
206     sumRHStile(bi,bj) = sumRHStile(bi,bj) + cg2d_b(I,J,bi,bj)
207 jahn 1.50 #endif
208 cnh 1.1 ENDDO
209     ENDDO
210     ENDDO
211     ENDDO
212 jmc 1.53 CALL EXCH_S3D_RL( cg2d_r, 1, myThid )
213 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
214 jmc 1.52 CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, err, 0, 0, myThid)
215     CALL GLOBAL_SUM_SINGLECPU_RL(cg2d_b, sumRHS, OLx, OLy, myThid)
216 jahn 1.50 #else
217 jmc 1.53 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
218     CALL GLOBAL_SUM_TILE_RL( sumRHStile, sumRHS, myThid )
219 jahn 1.50 #endif
220 adcroft 1.33 err = SQRT(err)
221     actualIts = 0
222     actualResidual = err
223 dimitri 1.40
224 jmc 1.53 IF ( debugLevel .GE. debLevZero ) THEN
225 heimbach 1.37 _BEGIN_MASTER( myThid )
226 jmc 1.45 WRITE(standardmessageunit,'(A,1P2E22.14)')
227     & ' cg2d: Sum(rhs),rhsMax = ', sumRHS,rhsMax
228 edhill 1.41 _END_MASTER( myThid )
229 jmc 1.53 ENDIF
230 cnh 1.1 C _BARRIER
231 adcroft 1.30 c _BEGIN_MASTER( myThid )
232 jmc 1.45 c WRITE(standardmessageunit,'(A,I6,1PE30.14)')
233     c & ' CG2D iters, err = ',
234 adcroft 1.30 c & actualIts, actualResidual
235 edhill 1.41 c _END_MASTER( myThid )
236 adcroft 1.33 firstResidual=actualResidual
237    
238     IF ( err .LT. cg2dTolerance ) GOTO 11
239 cnh 1.1
240     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
241 adcroft 1.30 DO 10 it2d=1, numIters
242 cnh 1.1
243     CcnhDebugStarts
244 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' residual = ',
245 cnh 1.14 C & actualResidual
246 cnh 1.1 CcnhDebugEnds
247     C-- Solve preconditioning equation and update
248     C-- conjugate direction vector "s".
249     DO bj=myByLo(myThid),myByHi(myThid)
250     DO bi=myBxLo(myThid),myBxHi(myThid)
251 jmc 1.48 eta_qrNtile(bi,bj) = 0. _d 0
252 mlosch 1.49 #ifdef TARGET_NEC_SX
253     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
254     #endif /* TARGET_NEC_SX */
255 cnh 1.1 DO J=1,sNy
256     DO I=1,sNx
257 jmc 1.45 cg2d_q(I,J,bi,bj) =
258 cnh 1.3 & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
259     & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
260     & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
261     & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
262     & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
263 cnh 1.4 CcnhDebugStarts
264     C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
265     CcnhDebugEnds
266 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
267     localBuf(I,J,bi,bj) =
268     & cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
269     #else
270 jmc 1.48 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
271 cnh 1.1 & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
272 jahn 1.50 #endif
273 cnh 1.1 ENDDO
274     ENDDO
275     ENDDO
276     ENDDO
277    
278 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
279 jmc 1.52 CALL GLOBAL_SUM_SINGLECPU_RL( localBuf,eta_qrN,0,0,myThid )
280 jahn 1.50 #else
281 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( eta_qrNtile,eta_qrN,myThid )
282 jahn 1.50 #endif
283 cnh 1.1 CcnhDebugStarts
284 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
285 cnh 1.1 CcnhDebugEnds
286 jmc 1.28 cgBeta = eta_qrN/eta_qrNM1
287 cnh 1.1 CcnhDebugStarts
288 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
289 cnh 1.1 CcnhDebugEnds
290 jmc 1.28 eta_qrNM1 = eta_qrN
291 cnh 1.1
292     DO bj=myByLo(myThid),myByHi(myThid)
293     DO bi=myBxLo(myThid),myBxHi(myThid)
294     DO J=1,sNy
295     DO I=1,sNx
296 cnh 1.14 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj)
297     & + cgBeta*cg2d_s(I,J,bi,bj)
298 cnh 1.1 ENDDO
299     ENDDO
300     ENDDO
301     ENDDO
302    
303     C-- Do exchanges that require messages i.e. between
304     C-- processes.
305 jmc 1.53 CALL EXCH_S3D_RL( cg2d_s, 1, myThid )
306 cnh 1.1
307     C== Evaluate laplace operator on conjugate gradient vector
308     C== q = A.s
309     DO bj=myByLo(myThid),myByHi(myThid)
310     DO bi=myBxLo(myThid),myBxHi(myThid)
311 jmc 1.48 alphaTile(bi,bj) = 0. _d 0
312 mlosch 1.49 #ifdef TARGET_NEC_SX
313     !CDIR OUTERUNROLL=CG2D_OUTERLOOPITERS
314     #endif /* TARGET_NEC_SX */
315 cnh 1.1 DO J=1,sNy
316     DO I=1,sNx
317 jmc 1.46 c ks = ksurfC(I,J,bi,bj)
318 jmc 1.45 cg2d_q(I,J,bi,bj) =
319 cnh 1.1 & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
320     & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
321     & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
322     & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
323 jmc 1.46 & +aC2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
324     c & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
325     c & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
326     c & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
327     c & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
328     c & -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)*recip_Bo(i,j,bi,bj)
329     c & *cg2d_s(I ,J ,bi,bj)
330     c & /deltaTMom/deltaTfreesurf*cg2dNorm
331 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
332     localBuf(I,J,bi,bj) = cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
333     #else
334 jmc 1.48 alphaTile(bi,bj) = alphaTile(bi,bj)
335     & + cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
336 jahn 1.50 #endif
337 cnh 1.1 ENDDO
338     ENDDO
339     ENDDO
340     ENDDO
341 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
342 jmc 1.52 CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, alpha, 0, 0, myThid)
343 jahn 1.50 #else
344 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( alphaTile, alpha, myThid )
345 jahn 1.50 #endif
346 cnh 1.1 CcnhDebugStarts
347 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
348 cnh 1.1 CcnhDebugEnds
349 jmc 1.28 alpha = eta_qrN/alpha
350 cnh 1.1 CcnhDebugStarts
351 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
352 cnh 1.1 CcnhDebugEnds
353 jmc 1.45
354 cnh 1.1 C== Update solution and residual vectors
355     C Now compute "interior" points.
356     DO bj=myByLo(myThid),myByHi(myThid)
357     DO bi=myBxLo(myThid),myBxHi(myThid)
358 jmc 1.48 errTile(bi,bj) = 0. _d 0
359 cnh 1.1 DO J=1,sNy
360     DO I=1,sNx
361     cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
362     cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
363 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
364     localBuf(I,J,bi,bj) = cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
365     #else
366 jmc 1.48 errTile(bi,bj) = errTile(bi,bj)
367     & + cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
368 jahn 1.50 #endif
369 cnh 1.1 ENDDO
370     ENDDO
371     ENDDO
372     ENDDO
373    
374 jahn 1.50 #ifdef CG2D_SINGLECPU_SUM
375 jmc 1.52 CALL GLOBAL_SUM_SINGLECPU_RL(localBuf, err, 0, 0, myThid)
376 jahn 1.50 #else
377 jmc 1.48 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
378 jahn 1.50 #endif
379 cnh 1.1 err = SQRT(err)
380     actualIts = it2d
381     actualResidual = err
382 jmc 1.47 IF ( debugLevel.GT.debLevB ) THEN
383     c IF ( DIFFERENT_MULTIPLE(monitorFreq,myTime,deltaTClock)
384     c & ) THEN
385     _BEGIN_MASTER( myThid )
386     WRITE(msgBuf,'(A,I6,A,1PE21.14)')
387     & ' cg2d: iter=', actualIts, ' ; resid.= ', actualResidual
388     CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
389     _END_MASTER( myThid )
390     c ENDIF
391     ENDIF
392 adcroft 1.33 IF ( err .LT. cg2dTolerance ) GOTO 11
393 jmc 1.47
394 jmc 1.53 CALL EXCH_S3D_RL( cg2d_r, 1, myThid )
395 cnh 1.13
396 cnh 1.1 10 CONTINUE
397     11 CONTINUE
398    
399 adcroft 1.33 IF (cg2dNormaliseRHS) THEN
400 cnh 1.1 C-- Un-normalise the answer
401 adcroft 1.33 DO bj=myByLo(myThid),myByHi(myThid)
402     DO bi=myBxLo(myThid),myBxHi(myThid)
403     DO J=1,sNy
404     DO I=1,sNx
405     cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
406     ENDDO
407     ENDDO
408 cnh 1.1 ENDDO
409     ENDDO
410 adcroft 1.33 ENDIF
411 cnh 1.1
412 adcroft 1.22 C The following exchange was moved up to solve_for_pressure
413     C for compatibility with TAMC.
414 jmc 1.51 C _EXCH_XY_RL(cg2d_x, myThid )
415 adcroft 1.30 c _BEGIN_MASTER( myThid )
416 jmc 1.45 c WRITE(*,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
417 adcroft 1.30 c & actualIts, actualResidual
418 edhill 1.41 c _END_MASTER( myThid )
419 adcroft 1.30
420     C-- Return parameters to caller
421 adcroft 1.33 lastResidual=actualResidual
422 adcroft 1.30 numIters=actualIts
423 cnh 1.1
424     CcnhDebugStarts
425 cnh 1.7 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
426 cnh 1.1 C err = 0. _d 0
427     C DO bj=myByLo(myThid),myByHi(myThid)
428     C DO bi=myBxLo(myThid),myBxHi(myThid)
429     C DO J=1,sNy
430     C DO I=1,sNx
431     C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
432     C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
433     C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
434     C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
435     C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
436     C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
437     C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
438     C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
439     C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
440 jmc 1.45 C err = err +
441 cnh 1.1 C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
442     C ENDDO
443     C ENDDO
444     C ENDDO
445     C ENDDO
446 jmc 1.51 C _GLOBAL_SUM_RL( err , myThid )
447 heimbach 1.31 C write(*,*) 'cg2d: Ax - b = ',SQRT(err)
448 cnh 1.1 CcnhDebugEnds
449    
450 adcroft 1.19 RETURN
451 cnh 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22