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

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

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


Revision 1.1 - (show annotations) (download)
Mon May 14 13:21:59 2012 UTC (12 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64o, checkpoint64a, checkpoint63r, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint64n, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint64b, checkpoint64e, checkpoint63q, checkpoint64d, checkpoint64c, checkpoint64g, checkpoint64f, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint63n, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint64i, checkpoint63o, checkpoint63p, checkpoint64h, checkpoint63s, checkpoint64k, checkpoint64, checkpoint65, checkpoint64j, checkpoint64m, checkpoint64l, HEAD
new CG-solver version (_EX0) for disconnected-tiles special case.

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

  ViewVC Help
Powered by ViewVC 1.1.22