/[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.25 - (show annotations) (download)
Fri May 11 23:34:06 2012 UTC (11 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint64, checkpoint65, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63n, checkpoint63o, checkpoint65o, checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, HEAD
Changes since 1.24: +56 -66 lines
import some minor changes from cg2d.F

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

  ViewVC Help
Powered by ViewVC 1.1.22