/[MITgcm]/MITgcm/model/src/cg3d.F
ViewVC logotype

Contents of /MITgcm/model/src/cg3d.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.24 - (show 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.23: +11 -16 lines
use new parameter "printResidualFreq" to print residual in CG iterations

1 C $Header: /u/gcmpack/MITgcm/model/src/cg3d.F,v 1.23 2010/01/17 21:55:48 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5 #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 CG3D_OUTERLOOPITERS
9 # define CG3D_OUTERLOOPITERS 10
10 #endif
11 #endif /* TARGET_NEC_SX */
12
13 CBOP
14 C !ROUTINE: CG3D
15 C !INTERFACE:
16 SUBROUTINE CG3D(
17 I cg3d_b,
18 U cg3d_x,
19 O firstResidual,
20 O lastResidual,
21 U numIters,
22 I myIter, myThid )
23 C !DESCRIPTION: \bv
24 C *==========================================================*
25 C | SUBROUTINE CG3D
26 C | o Three-dimensional grid problem conjugate-gradient
27 C | inverter (with preconditioner).
28 C *==========================================================*
29 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 seven-diagonal
32 C | matrix of the form that arises in the discrete
33 C | representation of the del^2 operator in a
34 C | three-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 C *==========================================================*
46 C \ev
47
48 C !USES:
49 IMPLICIT NONE
50 C === Global data ===
51 #include "SIZE.h"
52 #include "EEPARAMS.h"
53 #include "PARAMS.h"
54 #include "GRID.h"
55 #include "SURFACE.h"
56 #include "CG3D.h"
57
58 C !INPUT/OUTPUT PARAMETERS:
59 C === Routine arguments ===
60 C cg3d_b :: The source term or "right hand side"
61 C cg3d_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 C myIter :: Current iteration number in simulation
67 C myThid :: my Thread Id number
68 _RL cg3d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
69 _RL cg3d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
70 _RL firstResidual
71 _RL lastResidual
72 INTEGER numIters
73 INTEGER myIter
74 INTEGER myThid
75
76 #ifdef ALLOW_NONHYDROSTATIC
77
78 C !LOCAL VARIABLES:
79 C === Local variables ====
80 C actualIts :: Number of iterations taken
81 C actualResidual :: residual
82 C bi,bj :: tile indices
83 C eta_qrN :: Used in computing search directions
84 C eta_qrNM1 suffix N and NM1 denote current and
85 C cgBeta previous iterations respectively.
86 C alpha
87 C sumRHS :: Sum of right-hand-side. Sometimes this is a
88 C useful debuggin/trouble shooting diagnostic.
89 C For neumann problems sumRHS needs to be ~0.
90 C or they converge at a non-zero residual.
91 C err :: Measure of residual of Ax - b, usually the norm.
92 C i, j, k :: Loop counters
93 C it3d :: Loop counter for CG iterations
94 C msgBuf :: Informational/error message buffer
95 INTEGER actualIts
96 _RL actualResidual
97 INTEGER bi, bj
98 INTEGER i, j, k, it3d
99 INTEGER km1, kp1
100 _RL maskM1, maskP1
101 _RL err, errTile(nSx,nSy)
102 _RL eta_qrN,eta_qrNtile(nSx,nSy)
103 _RL eta_qrNM1
104 _RL cgBeta
105 _RL alpha , alphaTile(nSx,nSy)
106 _RL sumRHS, sumRHStile(nSx,nSy)
107 _RL rhsMax
108 _RL rhsNorm
109 CHARACTER*(MAX_LEN_MBUF) msgBuf
110 LOGICAL printResidual
111 _RL surfFac
112 #ifdef NONLIN_FRSURF
113 INTEGER ks
114 _RL surfTerm(sNx,sNy)
115 #endif /* NONLIN_FRSURF */
116 CEOP
117
118 IF ( select_rStar .NE. 0 ) THEN
119 surfFac = freeSurfFac
120 ELSE
121 surfFac = 0.
122 ENDIF
123 #ifdef NONLIN_FRSURF
124 DO j=1,sNy
125 DO i=1,sNx
126 surfTerm(i,j) = 0.
127 ENDDO
128 ENDDO
129 #endif /* NONLIN_FRSURF */
130
131 C-- Initialise inverter
132 eta_qrNM1 = 1. _d 0
133
134 C-- Normalise RHS
135 rhsMax = 0. _d 0
136 DO bj=myByLo(myThid),myByHi(myThid)
137 DO bi=myBxLo(myThid),myBxHi(myThid)
138 DO k=1,Nr
139 DO j=1,sNy
140 DO i=1,sNx
141 cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj)*cg3dNorm
142 & * maskC(i,j,k,bi,bj)
143 rhsMax = MAX(ABS(cg3d_b(i,j,k,bi,bj)),rhsMax)
144 ENDDO
145 ENDDO
146 ENDDO
147 ENDDO
148 ENDDO
149 _GLOBAL_MAX_RL( rhsMax, myThid )
150 rhsNorm = 1. _d 0
151 IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
152 DO bj=myByLo(myThid),myByHi(myThid)
153 DO bi=myBxLo(myThid),myBxHi(myThid)
154 DO k=1,Nr
155 DO j=1,sNy
156 DO i=1,sNx
157 cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj)*rhsNorm
158 cg3d_x(i,j,k,bi,bj) = cg3d_x(i,j,k,bi,bj)*rhsNorm
159 ENDDO
160 ENDDO
161 ENDDO
162 ENDDO
163 ENDDO
164
165 C-- Update overlaps
166 c _EXCH_XYZ_RL( cg3d_b, myThid )
167 _EXCH_XYZ_RL( cg3d_x, myThid )
168
169 C-- Initial residual calculation (with free-Surface term)
170 err = 0. _d 0
171 sumRHS = 0. _d 0
172 DO bj=myByLo(myThid),myByHi(myThid)
173 DO bi=myBxLo(myThid),myBxHi(myThid)
174 errTile(bi,bj) = 0. _d 0
175 sumRHStile(bi,bj) = 0. _d 0
176 #ifdef NONLIN_FRSURF
177 IF ( select_rStar .NE. 0 ) THEN
178 DO j=1,sNy
179 DO i=1,sNx
180 surfTerm(i,j) = 0.
181 ENDDO
182 ENDDO
183 DO k=1,Nr
184 DO j=1,sNy
185 DO i=1,sNx
186 surfTerm(i,j) = surfTerm(i,j)
187 & +cg3d_x(i,j,k,bi,bj)*drF(k)*h0FacC(i,j,k,bi,bj)
188 ENDDO
189 ENDDO
190 ENDDO
191 DO j=1,sNy
192 DO i=1,sNx
193 ks = kSurfC(i,j,bi,bj)
194 surfTerm(i,j) = surfTerm(i,j)*cg3dNorm
195 & *recip_Rcol(i,j,bi,bj)*recip_Rcol(i,j,bi,bj)
196 & *rA(i,j,bi,bj)*deepFac2F(ks)
197 & *recip_Bo(i,j,bi,bj)/deltaTMom/deltaTfreesurf
198 ENDDO
199 ENDDO
200 ENDIF
201 #endif /* NONLIN_FRSURF */
202 DO k=1,Nr
203 km1 = MAX(k-1, 1 )
204 kp1 = MIN(k+1, Nr)
205 maskM1 = 1. _d 0
206 maskP1 = 1. _d 0
207 IF ( k .EQ. 1 ) maskM1 = 0. _d 0
208 IF ( k .EQ. Nr) maskP1 = 0. _d 0
209 #ifdef TARGET_NEC_SX
210 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
211 #endif /* TARGET_NEC_SX */
212 DO j=1,sNy
213 DO i=1,sNx
214 cg3d_r(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj)
215 & -( 0.
216 & +aW3d( i, j, k, bi,bj)*cg3d_x(i-1,j, k, bi,bj)
217 & +aW3d(i+1,j, k, bi,bj)*cg3d_x(i+1,j, k, bi,bj)
218 & +aS3d( i, j, k, bi,bj)*cg3d_x( i,j-1,k, bi,bj)
219 & +aS3d( i,j+1,k, bi,bj)*cg3d_x( i,j+1,k, bi,bj)
220 & +aV3d( i, j, k, bi,bj)*cg3d_x( i, j,km1,bi,bj)*maskM1
221 & +aV3d( i, j,kp1,bi,bj)*cg3d_x( i, j,kp1,bi,bj)*maskP1
222 & +aC3d( i, j, k, bi,bj)*cg3d_x( i, j, k, bi,bj)
223 #ifdef NONLIN_FRSURF
224 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
225 #endif /* NONLIN_FRSURF */
226 & )
227 errTile(bi,bj) = errTile(bi,bj)
228 & +cg3d_r(i,j,k,bi,bj)*cg3d_r(i,j,k,bi,bj)
229 sumRHStile(bi,bj) = sumRHStile(bi,bj)+cg3d_b(i,j,k,bi,bj)
230 ENDDO
231 ENDDO
232 DO J=1-1,sNy+1
233 DO I=1-1,sNx+1
234 cg3d_s(i,j,k,bi,bj) = 0.
235 ENDDO
236 ENDDO
237 ENDDO
238 ENDDO
239 ENDDO
240 CALL EXCH_S3D_RL( cg3d_r, Nr, myThid )
241 c CALL EXCH_S3D_RL( cg3d_s, Nr, myThid )
242 CALL GLOBAL_SUM_TILE_RL( sumRHStile, sumRHS, myThid )
243 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
244 IF ( debugLevel.GE.debLevC .AND. diagFreq.GT.0. ) THEN
245 CALL WRITE_FLD_S3D_RL(
246 I 'cg3d_r_I', 'I10', 1, Nr, cg3d_r, myIter, myThid )
247 ENDIF
248
249 printResidual = .FALSE.
250 IF ( debugLevel .GE. debLevZero ) THEN
251 _BEGIN_MASTER( myThid )
252 printResidual = printResidualFreq.GE.1
253 WRITE(standardmessageunit,'(A,1P2E22.14)')
254 & ' cg3d: Sum(rhs),rhsMax = ',sumRHS,rhsMax
255 _END_MASTER( myThid )
256 ENDIF
257
258 actualIts = 0
259 actualResidual = SQRT(err)
260 firstResidual=actualResidual
261
262 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
263 DO 10 it3d=1, numIters
264
265 IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11
266 C-- Solve preconditioning equation and update
267 C-- conjugate direction vector "s".
268 C Note. On the next two loops over all tiles the inner loop ranges
269 C in sNx and sNy are expanded by 1 to avoid a communication
270 C step. However this entails a bit of gynamastics because we only
271 C want eta_qrN for the interior points.
272 eta_qrN = 0. _d 0
273 DO bj=myByLo(myThid),myByHi(myThid)
274 DO bi=myBxLo(myThid),myBxHi(myThid)
275 eta_qrNtile(bi,bj) = 0. _d 0
276 DO k=1,1
277 #ifdef TARGET_NEC_SX
278 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
279 #endif /* TARGET_NEC_SX */
280 DO j=1-1,sNy+1
281 DO i=1-1,sNx+1
282 cg3d_q(i,j,k,bi,bj) = zMC(i,j,k,bi,bj)
283 & *cg3d_r(i,j,k,bi,bj)
284 ENDDO
285 ENDDO
286 ENDDO
287 DO k=2,Nr
288 #ifdef TARGET_NEC_SX
289 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
290 #endif /* TARGET_NEC_SX */
291 DO j=1-1,sNy+1
292 DO i=1-1,sNx+1
293 cg3d_q(i,j,k,bi,bj) = zMC(i,j,k,bi,bj)
294 & *( cg3d_r(i,j,k,bi,bj)
295 & -zML(i,j,k,bi,bj)*cg3d_q(i,j,k-1,bi,bj)
296 & )
297 ENDDO
298 ENDDO
299 ENDDO
300 DO k=Nr,Nr
301 #ifdef TARGET_NEC_SX
302 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
303 #endif /* TARGET_NEC_SX */
304 DO j=1,sNy
305 DO i=1,sNx
306 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
307 & +cg3d_q(i,j,k,bi,bj)*cg3d_r(i,j,k,bi,bj)
308 ENDDO
309 ENDDO
310 ENDDO
311 DO k=Nr-1,1,-1
312 #ifdef TARGET_NEC_SX
313 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
314 #endif /* TARGET_NEC_SX */
315 DO j=1-1,sNy+1
316 DO i=1-1,sNx+1
317 cg3d_q(i,j,k,bi,bj) = cg3d_q(i,j,k,bi,bj)
318 & -zMU(i,j,k,bi,bj)*cg3d_q(i,j,k+1,bi,bj)
319 ENDDO
320 ENDDO
321 #ifdef TARGET_NEC_SX
322 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
323 #endif /* TARGET_NEC_SX */
324 DO j=1,sNy
325 DO i=1,sNx
326 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
327 & +cg3d_q(i,j,k,bi,bj)*cg3d_r(i,j,k,bi,bj)
328 ENDDO
329 ENDDO
330 ENDDO
331 ENDDO
332 ENDDO
333
334 CALL GLOBAL_SUM_TILE_RL( eta_qrNtile,eta_qrN,myThid )
335 CcnhDebugStarts
336 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' eta_qrN = ',eta_qrN
337 CcnhDebugEnds
338 cgBeta = eta_qrN/eta_qrNM1
339 CcnhDebugStarts
340 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' beta = ',cgBeta
341 CcnhDebugEnds
342 eta_qrNM1 = eta_qrN
343
344 DO bj=myByLo(myThid),myByHi(myThid)
345 DO bi=myBxLo(myThid),myBxHi(myThid)
346 DO k=1,Nr
347 DO j=1-1,sNy+1
348 DO i=1-1,sNx+1
349 cg3d_s(i,j,k,bi,bj) = cg3d_q(i,j,k,bi,bj)
350 & + cgBeta*cg3d_s(i,j,k,bi,bj)
351 ENDDO
352 ENDDO
353 ENDDO
354 ENDDO
355 ENDDO
356
357 C== Evaluate laplace operator on conjugate gradient vector
358 C== q = A.s
359 alpha = 0. _d 0
360 DO bj=myByLo(myThid),myByHi(myThid)
361 DO bi=myBxLo(myThid),myBxHi(myThid)
362 alphaTile(bi,bj) = 0. _d 0
363 #ifdef NONLIN_FRSURF
364 IF ( select_rStar .NE. 0 ) THEN
365 DO j=1,sNy
366 DO i=1,sNx
367 surfTerm(i,j) = 0.
368 ENDDO
369 ENDDO
370 DO k=1,Nr
371 DO j=1,sNy
372 DO i=1,sNx
373 surfTerm(i,j) = surfTerm(i,j)
374 & +cg3d_s(i,j,k,bi,bj)*drF(k)*h0FacC(i,j,k,bi,bj)
375 ENDDO
376 ENDDO
377 ENDDO
378 DO j=1,sNy
379 DO i=1,sNx
380 ks = kSurfC(i,j,bi,bj)
381 surfTerm(i,j) = surfTerm(i,j)*cg3dNorm
382 & *recip_Rcol(i,j,bi,bj)*recip_Rcol(i,j,bi,bj)
383 & *rA(i,j,bi,bj)*deepFac2F(ks)
384 & *recip_Bo(i,j,bi,bj)/deltaTMom/deltaTfreesurf
385 ENDDO
386 ENDDO
387 ENDIF
388 #endif /* NONLIN_FRSURF */
389 IF ( Nr .GT. 1 ) THEN
390 k=1
391 #ifdef TARGET_NEC_SX
392 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
393 #endif /* TARGET_NEC_SX */
394 DO j=1,sNy
395 DO i=1,sNx
396 cg3d_q(i,j,k,bi,bj) =
397 & aW3d( i, j, k, bi,bj)*cg3d_s(i-1,j, k, bi,bj)
398 & +aW3d(i+1,j, k, bi,bj)*cg3d_s(i+1,j, k, bi,bj)
399 & +aS3d( i, j, k, bi,bj)*cg3d_s( i,j-1,k, bi,bj)
400 & +aS3d( i,j+1,k, bi,bj)*cg3d_s( i,j+1,k, bi,bj)
401 & +aV3d( i, j,k+1,bi,bj)*cg3d_s( i, j,k+1,bi,bj)
402 & +aC3d( i, j, k, bi,bj)*cg3d_s( i, j, k, bi,bj)
403 #ifdef NONLIN_FRSURF
404 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
405 #endif /* NONLIN_FRSURF */
406 alphaTile(bi,bj) = alphaTile(bi,bj)
407 & +cg3d_s(i,j,k,bi,bj)*cg3d_q(i,j,k,bi,bj)
408 ENDDO
409 ENDDO
410 ELSE
411 k=1
412 #ifdef TARGET_NEC_SX
413 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
414 #endif /* TARGET_NEC_SX */
415 DO j=1,sNy
416 DO i=1,sNx
417 cg3d_q(i,j,k,bi,bj) =
418 & aW3d( i, j, k, bi,bj)*cg3d_s(i-1,j, k, bi,bj)
419 & +aW3d(i+1,j, k, bi,bj)*cg3d_s(i+1,j, k, bi,bj)
420 & +aS3d( i, j, k, bi,bj)*cg3d_s( i,j-1,k, bi,bj)
421 & +aS3d( i,j+1,k, bi,bj)*cg3d_s( i,j+1,k, bi,bj)
422 & +aC3d( i, j, k, bi,bj)*cg3d_s( i, j, k, bi,bj)
423 #ifdef NONLIN_FRSURF
424 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
425 #endif /* NONLIN_FRSURF */
426 alphaTile(bi,bj) = alphaTile(bi,bj)
427 & +cg3d_s(i,j,k,bi,bj)*cg3d_q(i,j,k,bi,bj)
428 ENDDO
429 ENDDO
430 ENDIF
431 DO k=2,Nr-1
432 #ifdef TARGET_NEC_SX
433 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
434 #endif /* TARGET_NEC_SX */
435 DO j=1,sNy
436 DO i=1,sNx
437 cg3d_q(i,j,k,bi,bj) =
438 & aW3d( i, j, k, bi,bj)*cg3d_s(i-1,j, k, bi,bj)
439 & +aW3d(i+1,j, k, bi,bj)*cg3d_s(i+1,j, k, bi,bj)
440 & +aS3d( i, j, k, bi,bj)*cg3d_s( i,j-1,k, bi,bj)
441 & +aS3d( i,j+1,k, bi,bj)*cg3d_s( i,j+1,k, bi,bj)
442 & +aV3d( i, j, k, bi,bj)*cg3d_s( i, j,k-1,bi,bj)
443 & +aV3d( i, j,k+1,bi,bj)*cg3d_s( i, j,k+1,bi,bj)
444 & +aC3d( i, j, k, bi,bj)*cg3d_s( i, j, k, bi,bj)
445 #ifdef NONLIN_FRSURF
446 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
447 #endif /* NONLIN_FRSURF */
448 alphaTile(bi,bj) = alphaTile(bi,bj)
449 & +cg3d_s(i,j,k,bi,bj)*cg3d_q(i,j,k,bi,bj)
450 ENDDO
451 ENDDO
452 ENDDO
453 IF ( Nr .GT. 1 ) THEN
454 k=Nr
455 #ifdef TARGET_NEC_SX
456 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
457 #endif /* TARGET_NEC_SX */
458 DO j=1,sNy
459 DO i=1,sNx
460 cg3d_q(i,j,k,bi,bj) =
461 & aW3d( i, j, k, bi,bj)*cg3d_s(i-1,j, k, bi,bj)
462 & +aW3d(i+1,j, k, bi,bj)*cg3d_s(i+1,j, k, bi,bj)
463 & +aS3d( i, j, k, bi,bj)*cg3d_s( i,j-1,k, bi,bj)
464 & +aS3d( i,j+1,k, bi,bj)*cg3d_s( i,j+1,k, bi,bj)
465 & +aV3d( i, j, k, bi,bj)*cg3d_s( i, j,k-1,bi,bj)
466 & +aC3d( i, j, k, bi,bj)*cg3d_s( i, j, k, bi,bj)
467 #ifdef NONLIN_FRSURF
468 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
469 #endif /* NONLIN_FRSURF */
470 alphaTile(bi,bj) = alphaTile(bi,bj)
471 & +cg3d_s(i,j,k,bi,bj)*cg3d_q(i,j,k,bi,bj)
472 ENDDO
473 ENDDO
474 ENDIF
475 ENDDO
476 ENDDO
477 CALL GLOBAL_SUM_TILE_RL( alphaTile, alpha, myThid )
478 CcnhDebugStarts
479 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' SUM(s*q)= ',alpha
480 CcnhDebugEnds
481 alpha = eta_qrN/alpha
482 CcnhDebugStarts
483 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' alpha= ',alpha
484 CcnhDebugEnds
485
486 C== Update solution and residual vectors
487 C Now compute "interior" points.
488 err = 0. _d 0
489 DO bj=myByLo(myThid),myByHi(myThid)
490 DO bi=myBxLo(myThid),myBxHi(myThid)
491 errTile(bi,bj) = 0. _d 0
492 DO k=1,Nr
493 #ifdef TARGET_NEC_SX
494 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
495 #endif /* TARGET_NEC_SX */
496 DO j=1,sNy
497 DO i=1,sNx
498 cg3d_x(i,j,k,bi,bj)=cg3d_x(i,j,k,bi,bj)
499 & +alpha*cg3d_s(i,j,k,bi,bj)
500 cg3d_r(i,j,k,bi,bj)=cg3d_r(i,j,k,bi,bj)
501 & -alpha*cg3d_q(i,j,k,bi,bj)
502 errTile(bi,bj) = errTile(bi,bj)
503 & +cg3d_r(i,j,k,bi,bj)*cg3d_r(i,j,k,bi,bj)
504 ENDDO
505 ENDDO
506 ENDDO
507 ENDDO
508 ENDDO
509
510 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
511 err = SQRT(err)
512 actualIts = it3d
513 actualResidual = err
514 IF ( printResidual ) THEN
515 IF ( MOD( it3d-1, printResidualFreq ).EQ.0 ) THEN
516 WRITE(msgBuf,'(A,I6,A,1PE21.14)')
517 & ' cg3d: iter=', actualIts, ' ; resid.= ', actualResidual
518 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
519 & SQUEEZE_RIGHT, myThid )
520 ENDIF
521 ENDIF
522 IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11
523 CALL EXCH_S3D_RL( cg3d_r, Nr, myThid )
524
525 10 CONTINUE
526 11 CONTINUE
527
528 IF ( debugLevel.GE.debLevC .AND. diagFreq.GT.0. ) THEN
529 CALL WRITE_FLD_S3D_RL(
530 I 'cg3d_r_F', 'I10', 1, Nr, cg3d_r, myIter, myThid )
531 ENDIF
532
533 C-- Un-normalise the answer
534 DO bj=myByLo(myThid),myByHi(myThid)
535 DO bi=myBxLo(myThid),myBxHi(myThid)
536 DO k=1,Nr
537 DO j=1,sNy
538 DO i=1,sNx
539 cg3d_x(i,j,k,bi,bj) = cg3d_x(i,j,k,bi,bj)/rhsNorm
540 ENDDO
541 ENDDO
542 ENDDO
543 ENDDO
544 ENDDO
545
546 lastResidual = actualResidual
547 numIters = actualIts
548
549 #endif /* ALLOW_NONHYDROSTATIC */
550
551 RETURN
552 END

  ViewVC Help
Powered by ViewVC 1.1.22