/[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.23 - (show annotations) (download)
Sun Jan 17 21:55:48 2010 UTC (14 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62y, checkpoint62x, checkpoint62b
Changes since 1.22: +3 -3 lines
use the right type (RL) of S/R to write residual field (debug)

1 C $Header: /u/gcmpack/MITgcm/model/src/cg3d.F,v 1.22 2009/11/30 15:18:21 mlosch 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 _RL surfFac
111 #ifdef NONLIN_FRSURF
112 INTEGER ks
113 _RL surfTerm(sNx,sNy)
114 #endif /* NONLIN_FRSURF */
115 CEOP
116
117 IF ( select_rStar .NE. 0 ) THEN
118 surfFac = freeSurfFac
119 ELSE
120 surfFac = 0.
121 ENDIF
122 #ifdef NONLIN_FRSURF
123 DO j=1,sNy
124 DO i=1,sNx
125 surfTerm(i,j) = 0.
126 ENDDO
127 ENDDO
128 #endif /* NONLIN_FRSURF */
129
130 C-- Initialise inverter
131 eta_qrNM1 = 1. _d 0
132
133 C-- Normalise RHS
134 rhsMax = 0. _d 0
135 DO bj=myByLo(myThid),myByHi(myThid)
136 DO bi=myBxLo(myThid),myBxHi(myThid)
137 DO k=1,Nr
138 DO j=1,sNy
139 DO i=1,sNx
140 cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj)*cg3dNorm
141 & * maskC(i,j,k,bi,bj)
142 rhsMax = MAX(ABS(cg3d_b(i,j,k,bi,bj)),rhsMax)
143 ENDDO
144 ENDDO
145 ENDDO
146 ENDDO
147 ENDDO
148 _GLOBAL_MAX_RL( rhsMax, myThid )
149 rhsNorm = 1. _d 0
150 IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
151 DO bj=myByLo(myThid),myByHi(myThid)
152 DO bi=myBxLo(myThid),myBxHi(myThid)
153 DO k=1,Nr
154 DO j=1,sNy
155 DO i=1,sNx
156 cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj)*rhsNorm
157 cg3d_x(i,j,k,bi,bj) = cg3d_x(i,j,k,bi,bj)*rhsNorm
158 ENDDO
159 ENDDO
160 ENDDO
161 ENDDO
162 ENDDO
163
164 C-- Update overlaps
165 c _EXCH_XYZ_RL( cg3d_b, myThid )
166 _EXCH_XYZ_RL( cg3d_x, myThid )
167
168 C-- Initial residual calculation (with free-Surface term)
169 err = 0. _d 0
170 sumRHS = 0. _d 0
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=1-1,sNy+1
232 DO I=1-1,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 c CALL EXCH_S3D_RL( cg3d_s, Nr, myThid )
241 CALL GLOBAL_SUM_TILE_RL( sumRHStile, sumRHS, myThid )
242 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
243 IF ( debugLevel.GT.debLevB .AND. diagFreq.GT.0. ) THEN
244 CALL WRITE_FLD_S3D_RL(
245 I 'cg3d_r_I', 'I10', 1, Nr, cg3d_r, myIter, myThid )
246 ENDIF
247
248 IF ( debugLevel .GE. debLevZero ) THEN
249 _BEGIN_MASTER( myThid )
250 WRITE(standardmessageunit,'(A,1P2E22.14)')
251 & ' cg3d: Sum(rhs),rhsMax = ',sumRHS,rhsMax
252 _END_MASTER( myThid )
253 ENDIF
254
255 actualIts = 0
256 actualResidual = SQRT(err)
257 firstResidual=actualResidual
258
259 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
260 DO 10 it3d=1, numIters
261
262 CcnhDebugStarts
263 c IF ( mod(it3d-1,10).EQ.0)
264 c & WRITE(*,*) ' CG3D: Iteration ',it3d-1,
265 c & ' residual = ',actualResidual
266 CcnhDebugEnds
267 IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11
268 C-- Solve preconditioning equation and update
269 C-- conjugate direction vector "s".
270 C Note. On the next two loops over all tiles the inner loop ranges
271 C in sNx and sNy are expanded by 1 to avoid a communication
272 C step. However this entails a bit of gynamastics because we only
273 C want eta_qrN for the interior points.
274 eta_qrN = 0. _d 0
275 DO bj=myByLo(myThid),myByHi(myThid)
276 DO bi=myBxLo(myThid),myBxHi(myThid)
277 eta_qrNtile(bi,bj) = 0. _d 0
278 DO k=1,1
279 #ifdef TARGET_NEC_SX
280 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
281 #endif /* TARGET_NEC_SX */
282 DO j=1-1,sNy+1
283 DO i=1-1,sNx+1
284 cg3d_q(i,j,k,bi,bj) = zMC(i,j,k,bi,bj)
285 & *cg3d_r(i,j,k,bi,bj)
286 ENDDO
287 ENDDO
288 ENDDO
289 DO k=2,Nr
290 #ifdef TARGET_NEC_SX
291 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
292 #endif /* TARGET_NEC_SX */
293 DO j=1-1,sNy+1
294 DO i=1-1,sNx+1
295 cg3d_q(i,j,k,bi,bj) = zMC(i,j,k,bi,bj)
296 & *( cg3d_r(i,j,k,bi,bj)
297 & -zML(i,j,k,bi,bj)*cg3d_q(i,j,k-1,bi,bj)
298 & )
299 ENDDO
300 ENDDO
301 ENDDO
302 DO k=Nr,Nr
303 #ifdef TARGET_NEC_SX
304 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
305 #endif /* TARGET_NEC_SX */
306 DO j=1,sNy
307 DO i=1,sNx
308 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
309 & +cg3d_q(i,j,k,bi,bj)*cg3d_r(i,j,k,bi,bj)
310 ENDDO
311 ENDDO
312 ENDDO
313 DO k=Nr-1,1,-1
314 #ifdef TARGET_NEC_SX
315 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
316 #endif /* TARGET_NEC_SX */
317 DO j=1-1,sNy+1
318 DO i=1-1,sNx+1
319 cg3d_q(i,j,k,bi,bj) = cg3d_q(i,j,k,bi,bj)
320 & -zMU(i,j,k,bi,bj)*cg3d_q(i,j,k+1,bi,bj)
321 ENDDO
322 ENDDO
323 #ifdef TARGET_NEC_SX
324 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
325 #endif /* TARGET_NEC_SX */
326 DO j=1,sNy
327 DO i=1,sNx
328 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
329 & +cg3d_q(i,j,k,bi,bj)*cg3d_r(i,j,k,bi,bj)
330 ENDDO
331 ENDDO
332 ENDDO
333 ENDDO
334 ENDDO
335
336 CALL GLOBAL_SUM_TILE_RL( eta_qrNtile,eta_qrN,myThid )
337 CcnhDebugStarts
338 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' eta_qrN = ',eta_qrN
339 CcnhDebugEnds
340 cgBeta = eta_qrN/eta_qrNM1
341 CcnhDebugStarts
342 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' beta = ',cgBeta
343 CcnhDebugEnds
344 eta_qrNM1 = eta_qrN
345
346 DO bj=myByLo(myThid),myByHi(myThid)
347 DO bi=myBxLo(myThid),myBxHi(myThid)
348 DO k=1,Nr
349 DO j=1-1,sNy+1
350 DO i=1-1,sNx+1
351 cg3d_s(i,j,k,bi,bj) = cg3d_q(i,j,k,bi,bj)
352 & + cgBeta*cg3d_s(i,j,k,bi,bj)
353 ENDDO
354 ENDDO
355 ENDDO
356 ENDDO
357 ENDDO
358
359 C== Evaluate laplace operator on conjugate gradient vector
360 C== q = A.s
361 alpha = 0. _d 0
362 DO bj=myByLo(myThid),myByHi(myThid)
363 DO bi=myBxLo(myThid),myBxHi(myThid)
364 alphaTile(bi,bj) = 0. _d 0
365 #ifdef NONLIN_FRSURF
366 IF ( select_rStar .NE. 0 ) THEN
367 DO j=1,sNy
368 DO i=1,sNx
369 surfTerm(i,j) = 0.
370 ENDDO
371 ENDDO
372 DO k=1,Nr
373 DO j=1,sNy
374 DO i=1,sNx
375 surfTerm(i,j) = surfTerm(i,j)
376 & +cg3d_s(i,j,k,bi,bj)*drF(k)*h0FacC(i,j,k,bi,bj)
377 ENDDO
378 ENDDO
379 ENDDO
380 DO j=1,sNy
381 DO i=1,sNx
382 ks = ksurfC(i,j,bi,bj)
383 surfTerm(i,j) = surfTerm(i,j)*cg3dNorm
384 & *recip_Rcol(i,j,bi,bj)*recip_Rcol(i,j,bi,bj)
385 & *rA(i,j,bi,bj)*deepFac2F(ks)
386 & *recip_Bo(i,j,bi,bj)/deltaTMom/deltaTfreesurf
387 ENDDO
388 ENDDO
389 ENDIF
390 #endif /* NONLIN_FRSURF */
391 IF ( Nr .GT. 1 ) THEN
392 k=1
393 #ifdef TARGET_NEC_SX
394 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
395 #endif /* TARGET_NEC_SX */
396 DO j=1,sNy
397 DO i=1,sNx
398 cg3d_q(i,j,k,bi,bj) =
399 & aW3d( i, j, k, bi,bj)*cg3d_s(i-1,j, k, bi,bj)
400 & +aW3d(i+1,j, k, bi,bj)*cg3d_s(i+1,j, k, bi,bj)
401 & +aS3d( i, j, k, bi,bj)*cg3d_s( i,j-1,k, bi,bj)
402 & +aS3d( i,j+1,k, bi,bj)*cg3d_s( i,j+1,k, bi,bj)
403 & +aV3d( i, j,k+1,bi,bj)*cg3d_s( i, j,k+1,bi,bj)
404 & +aC3d( i, j, k, bi,bj)*cg3d_s( i, j, k, bi,bj)
405 #ifdef NONLIN_FRSURF
406 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
407 #endif /* NONLIN_FRSURF */
408 alphaTile(bi,bj) = alphaTile(bi,bj)
409 & +cg3d_s(i,j,k,bi,bj)*cg3d_q(i,j,k,bi,bj)
410 ENDDO
411 ENDDO
412 ELSE
413 k=1
414 #ifdef TARGET_NEC_SX
415 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
416 #endif /* TARGET_NEC_SX */
417 DO j=1,sNy
418 DO i=1,sNx
419 cg3d_q(i,j,k,bi,bj) =
420 & aW3d( i, j, k, bi,bj)*cg3d_s(i-1,j, k, bi,bj)
421 & +aW3d(i+1,j, k, bi,bj)*cg3d_s(i+1,j, k, bi,bj)
422 & +aS3d( i, j, k, bi,bj)*cg3d_s( i,j-1,k, bi,bj)
423 & +aS3d( i,j+1,k, bi,bj)*cg3d_s( i,j+1,k, bi,bj)
424 & +aC3d( i, j, k, bi,bj)*cg3d_s( i, j, k, bi,bj)
425 #ifdef NONLIN_FRSURF
426 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
427 #endif /* NONLIN_FRSURF */
428 alphaTile(bi,bj) = alphaTile(bi,bj)
429 & +cg3d_s(i,j,k,bi,bj)*cg3d_q(i,j,k,bi,bj)
430 ENDDO
431 ENDDO
432 ENDIF
433 DO k=2,Nr-1
434 #ifdef TARGET_NEC_SX
435 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
436 #endif /* TARGET_NEC_SX */
437 DO j=1,sNy
438 DO i=1,sNx
439 cg3d_q(i,j,k,bi,bj) =
440 & aW3d( i, j, k, bi,bj)*cg3d_s(i-1,j, k, bi,bj)
441 & +aW3d(i+1,j, k, bi,bj)*cg3d_s(i+1,j, k, bi,bj)
442 & +aS3d( i, j, k, bi,bj)*cg3d_s( i,j-1,k, bi,bj)
443 & +aS3d( i,j+1,k, bi,bj)*cg3d_s( i,j+1,k, bi,bj)
444 & +aV3d( i, j, k, bi,bj)*cg3d_s( i, j,k-1,bi,bj)
445 & +aV3d( i, j,k+1,bi,bj)*cg3d_s( i, j,k+1,bi,bj)
446 & +aC3d( i, j, k, bi,bj)*cg3d_s( i, j, k, bi,bj)
447 #ifdef NONLIN_FRSURF
448 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
449 #endif /* NONLIN_FRSURF */
450 alphaTile(bi,bj) = alphaTile(bi,bj)
451 & +cg3d_s(i,j,k,bi,bj)*cg3d_q(i,j,k,bi,bj)
452 ENDDO
453 ENDDO
454 ENDDO
455 IF ( Nr .GT. 1 ) THEN
456 k=Nr
457 #ifdef TARGET_NEC_SX
458 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
459 #endif /* TARGET_NEC_SX */
460 DO j=1,sNy
461 DO i=1,sNx
462 cg3d_q(i,j,k,bi,bj) =
463 & aW3d( i, j, k, bi,bj)*cg3d_s(i-1,j, k, bi,bj)
464 & +aW3d(i+1,j, k, bi,bj)*cg3d_s(i+1,j, k, bi,bj)
465 & +aS3d( i, j, k, bi,bj)*cg3d_s( i,j-1,k, bi,bj)
466 & +aS3d( i,j+1,k, bi,bj)*cg3d_s( i,j+1,k, bi,bj)
467 & +aV3d( i, j, k, bi,bj)*cg3d_s( i, j,k-1,bi,bj)
468 & +aC3d( i, j, k, bi,bj)*cg3d_s( i, j, k, bi,bj)
469 #ifdef NONLIN_FRSURF
470 & -surfFac*surfTerm(i,j)*drF(k)*h0FacC(i,j,k,bi,bj)
471 #endif /* NONLIN_FRSURF */
472 alphaTile(bi,bj) = alphaTile(bi,bj)
473 & +cg3d_s(i,j,k,bi,bj)*cg3d_q(i,j,k,bi,bj)
474 ENDDO
475 ENDDO
476 ENDIF
477 ENDDO
478 ENDDO
479 CALL GLOBAL_SUM_TILE_RL( alphaTile, alpha, myThid )
480 CcnhDebugStarts
481 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' SUM(s*q)= ',alpha
482 CcnhDebugEnds
483 alpha = eta_qrN/alpha
484 CcnhDebugStarts
485 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' alpha= ',alpha
486 CcnhDebugEnds
487
488 C== Update solution and residual vectors
489 C Now compute "interior" points.
490 err = 0. _d 0
491 DO bj=myByLo(myThid),myByHi(myThid)
492 DO bi=myBxLo(myThid),myBxHi(myThid)
493 errTile(bi,bj) = 0. _d 0
494 DO k=1,Nr
495 #ifdef TARGET_NEC_SX
496 !CDIR OUTERUNROLL=CG3D_OUTERLOOPITERS
497 #endif /* TARGET_NEC_SX */
498 DO j=1,sNy
499 DO i=1,sNx
500 cg3d_x(i,j,k,bi,bj)=cg3d_x(i,j,k,bi,bj)
501 & +alpha*cg3d_s(i,j,k,bi,bj)
502 cg3d_r(i,j,k,bi,bj)=cg3d_r(i,j,k,bi,bj)
503 & -alpha*cg3d_q(i,j,k,bi,bj)
504 errTile(bi,bj) = errTile(bi,bj)
505 & +cg3d_r(i,j,k,bi,bj)*cg3d_r(i,j,k,bi,bj)
506 ENDDO
507 ENDDO
508 ENDDO
509 ENDDO
510 ENDDO
511
512 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
513 err = SQRT(err)
514 actualIts = it3d
515 actualResidual = err
516 IF ( debugLevel.GT.debLevB ) THEN
517 c IF ( DIFFERENT_MULTIPLE(monitorFreq,myTime,deltaTClock)
518 c & ) THEN
519 _BEGIN_MASTER( myThid )
520 WRITE(msgBuf,'(A,I6,A,1PE21.14)')
521 & ' cg3d: iter=', actualIts, ' ; resid.= ', actualResidual
522 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
523 & SQUEEZE_RIGHT, myThid )
524 _END_MASTER( myThid )
525 c ENDIF
526 ENDIF
527 IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11
528 CALL EXCH_S3D_RL( cg3d_r, Nr, myThid )
529
530 10 CONTINUE
531 11 CONTINUE
532
533 IF ( debugLevel.GT.debLevB .AND. diagFreq.GT.0. ) THEN
534 CALL WRITE_FLD_S3D_RL(
535 I 'cg3d_r_F', 'I10', 1, Nr, cg3d_r, myIter, myThid )
536 ENDIF
537
538 C-- Un-normalise the answer
539 DO bj=myByLo(myThid),myByHi(myThid)
540 DO bi=myBxLo(myThid),myBxHi(myThid)
541 DO k=1,Nr
542 DO j=1,sNy
543 DO i=1,sNx
544 cg3d_x(i,j,k,bi,bj) = cg3d_x(i,j,k,bi,bj)/rhsNorm
545 ENDDO
546 ENDDO
547 ENDDO
548 ENDDO
549 ENDDO
550
551 lastResidual = actualResidual
552 numIters = actualIts
553
554 #endif /* ALLOW_NONHYDROSTATIC */
555
556 RETURN
557 END

  ViewVC Help
Powered by ViewVC 1.1.22