/[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.15 - (show annotations) (download)
Fri Feb 4 19:30:33 2005 UTC (19 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_pre, checkpoint57s_post, checkpoint57g_post, checkpoint57r_post, checkpoint57d_post, checkpoint57i_post, checkpoint57n_post, checkpoint57l_post, checkpoint57t_post, checkpoint57v_post, checkpoint57f_post, checkpoint57h_pre, checkpoint57h_post, checkpoint57e_post, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, eckpoint57e_pre, checkpoint57h_done, checkpoint57j_post, checkpoint57f_pre, checkpoint57o_post, checkpoint57k_post, checkpoint57w_post
Changes since 1.14: +34 -17 lines
compute sum over each tile first, then sum different tiles

1 C $Header: /u/gcmpack/MITgcm/model/src/cg3d.F,v 1.14 2004/12/14 16:54:08 edhill Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #define VERBOSE
8
9 CBOP
10 C !ROUTINE: CG3D
11 C !INTERFACE:
12 SUBROUTINE CG3D(
13 I cg3d_b,
14 U cg3d_x,
15 O firstResidual,
16 O lastResidual,
17 U numIters,
18 I myThid )
19 C !DESCRIPTION: \bv
20 C *==========================================================*
21 C | SUBROUTINE CG3D
22 C | o Three-dimensional grid problem conjugate-gradient
23 C | inverter (with preconditioner).
24 C *==========================================================*
25 C | Con. grad is an iterative procedure for solving Ax = b.
26 C | It requires the A be symmetric.
27 C | This implementation assumes A is a seven-diagonal
28 C | matrix of the form that arises in the discrete
29 C | representation of the del^2 operator in a
30 C | three-dimensional space.
31 C | Notes:
32 C | ======
33 C | This implementation can support shared-memory
34 C | multi-threaded execution. In order to do this COMMON
35 C | blocks are used for many of the arrays - even ones that
36 C | are only used for intermedaite results. This design is
37 C | OK if you want to all the threads to collaborate on
38 C | solving the same problem. On the other hand if you want
39 C | the threads to solve several different problems
40 C | concurrently this implementation will not work.
41 C *==========================================================*
42 C \ev
43
44 C !USES:
45 IMPLICIT NONE
46 C === Global data ===
47 #include "SIZE.h"
48 #include "EEPARAMS.h"
49 #include "PARAMS.h"
50 #include "GRID.h"
51 #include "CG3D.h"
52
53 C !INPUT/OUTPUT PARAMETERS:
54 C === Routine arguments ===
55 C myThid - Thread on which I am working.
56 C cg2d_b - The source term or "right hand side"
57 C cg2d_x - The solution
58 C firstResidual - the initial residual before any iterations
59 C lastResidual - the actual residual reached
60 C numIters - Entry: the maximum number of iterations allowed
61 C Exit: the actual number of iterations used
62 _RL cg3d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
63 _RL cg3d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
64 _RL firstResidual
65 _RL lastResidual
66 INTEGER numIters
67 INTEGER myThid
68
69
70 #ifdef ALLOW_NONHYDROSTATIC
71
72 C !LOCAL VARIABLES:
73 C === Local variables ====
74 C actualIts - Number of iterations taken
75 C actualResidual - residual
76 C bi - Block index in X and Y.
77 C bj
78 C eta_qrN - Used in computing search directions
79 C eta_qrNM1 suffix N and NM1 denote current and
80 C cgBeta previous iterations respectively.
81 C alpha
82 C sumRHS - Sum of right-hand-side. Sometimes this is a
83 C useful debuggin/trouble shooting diagnostic.
84 C For neumann problems sumRHS needs to be ~0.
85 C or they converge at a non-zero residual.
86 C err - Measure of residual of Ax - b, usually the norm.
87 C I, J, N - Loop counters ( N counts CG iterations )
88 INTEGER actualIts
89 _RL actualResidual
90 INTEGER bi, bj
91 INTEGER I, J, K, it3d
92 INTEGER KM1, KP1
93 _RL err, errTile
94 _RL eta_qrN, eta_qrNtile
95 _RL eta_qrNM1
96 _RL cgBeta
97 _RL alpha , alphaTile
98 _RL sumRHS, sumRHStile
99 _RL rhsMax
100 _RL rhsNorm
101
102 INTEGER OLw
103 INTEGER OLe
104 INTEGER OLn
105 INTEGER OLs
106 INTEGER exchWidthX
107 INTEGER exchWidthY
108 INTEGER myNz
109 _RL topLevTerm
110 CEOP
111
112 ceh3 needs an IF ( useNONHYDROSTATIC ) THEN
113
114
115 C-- Initialise inverter
116 eta_qrNM1 = 1. D0
117
118 C-- Normalise RHS
119 rhsMax = 0. _d 0
120 DO bj=myByLo(myThid),myByHi(myThid)
121 DO bi=myBxLo(myThid),myBxHi(myThid)
122 DO K=1,Nr
123 DO J=1,sNy
124 DO I=1,sNx
125 cg3d_b(I,J,K,bi,bj) = cg3d_b(I,J,K,bi,bj)*cg3dNorm
126 rhsMax = MAX(ABS(cg3d_b(I,J,K,bi,bj)),rhsMax)
127 ENDDO
128 ENDDO
129 ENDDO
130 ENDDO
131 ENDDO
132 _GLOBAL_MAX_R8( rhsMax, myThid )
133 rhsNorm = 1. _d 0
134 IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
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)*rhsNorm
141 cg3d_x(I,J,K,bi,bj) = cg3d_x(I,J,K,bi,bj)*rhsNorm
142 ENDDO
143 ENDDO
144 ENDDO
145 ENDDO
146 ENDDO
147
148 C-- Update overlaps
149 _EXCH_XYZ_R8( cg3d_b, myThid )
150 _EXCH_XYZ_R8( cg3d_x, myThid )
151
152 C-- Initial residual calculation (with free-Surface term)
153 err = 0. _d 0
154 sumRHS = 0. _d 0
155 DO bj=myByLo(myThid),myByHi(myThid)
156 DO bi=myBxLo(myThid),myBxHi(myThid)
157 errTile = 0. _d 0
158 sumRHStile = 0. _d 0
159 DO K=1,Nr
160 KM1 = K-1
161 IF ( K .EQ. 1 ) KM1 = 1
162 KP1 = K+1
163 IF ( K .EQ. Nr ) KP1 = 1
164 topLevTerm = 0.
165 IF ( K .EQ. 1) topLevTerm = freeSurfFac*cg3dNorm*
166 & (horiVertRatio/gravity)/deltaTMom/deltaTMom
167 DO J=1,sNy
168 DO I=1,sNx
169 cg3d_s(I,J,K,bi,bj) = 0.
170 cg3d_r(I,J,K,bi,bj) = cg3d_b(I,J,K,bi,bj) -( 0.
171 & +aW3d(I ,J ,K ,bi,bj)*cg3d_x(I-1,J ,K ,bi,bj)
172 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_x(I+1,J ,K ,bi,bj)
173 & +aS3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J-1,K ,bi,bj)
174 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_x(I ,J+1,K ,bi,bj)
175 & +aV3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,KM1,bi,bj)
176 & +aV3d(I ,J ,KP1,bi,bj)*cg3d_x(I ,J ,KP1,bi,bj)
177 & -aW3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
178 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
179 & -aS3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
180 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
181 & -aV3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
182 & -aV3d(I ,J ,KP1,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
183 & -topLevTerm*_rA(I,J,bi,bj)*cg3d_x(I,J,K,bi,bj)
184 & )
185 errTile = errTile
186 & +cg3d_r(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
187 sumRHStile = sumRHStile
188 & +cg3d_b(I,J,K,bi,bj)
189 ENDDO
190 ENDDO
191 ENDDO
192 err = err + errTile
193 sumRHS = sumRHS + sumRHStile
194 ENDDO
195 ENDDO
196 C _EXCH_XYZ_R8( cg3d_r, myThid )
197 OLw = 1
198 OLe = 1
199 OLn = 1
200 OLs = 1
201 exchWidthX = 1
202 exchWidthY = 1
203 myNz = Nr
204 CALL EXCH_RL( cg3d_r,
205 I OLw, OLe, OLs, OLn, myNz,
206 I exchWidthX, exchWidthY,
207 I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
208 C _EXCH_XYZ_R8( cg3d_s, myThid )
209 OLw = 1
210 OLe = 1
211 OLn = 1
212 OLs = 1
213 exchWidthX = 1
214 exchWidthY = 1
215 myNz = Nr
216 CALL EXCH_RL( cg3d_s,
217 I OLw, OLe, OLs, OLn, myNz,
218 I exchWidthX, exchWidthY,
219 I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
220 _GLOBAL_SUM_R8( sumRHS, myThid )
221 _GLOBAL_SUM_R8( err , myThid )
222
223 IF ( debugLevel .GE. debLevZero ) THEN
224 _BEGIN_MASTER( myThid )
225 write(*,'(A,1P2E22.14)')
226 & ' cg3d: Sum(rhs),rhsMax = ',sumRHS,rhsMax
227 _END_MASTER( myThid )
228 ENDIF
229
230 actualIts = 0
231 actualResidual = SQRT(err)
232 C _BARRIER
233 c _BEGIN_MASTER( myThid )
234 c WRITE(*,'(A,I6,1PE30.14)') ' CG3D iters, err = ',
235 c & actualIts, actualResidual
236 c _END_MASTER( myThid )
237 firstResidual=actualResidual
238
239 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
240 DO 10 it3d=1, cg3dMaxIters
241
242 CcnhDebugStarts
243 #ifdef VERBOSE
244 c IF ( mod(it3d-1,10).EQ.0)
245 c & WRITE(*,*) ' CG3D: Iteration ',it3d-1,
246 c & ' residual = ',actualResidual
247 #endif
248 CcnhDebugEnds
249 IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11
250 C-- Solve preconditioning equation and update
251 C-- conjugate direction vector "s".
252 C Note. On the next to loops over all tiles the inner loop ranges
253 C in sNx and sNy are expanded by 1 to avoid a communication
254 C step. However this entails a bit of gynamastics because we only
255 C want eta_qrN for the interior points.
256 eta_qrN = 0. _d 0
257 DO bj=myByLo(myThid),myByHi(myThid)
258 DO bi=myBxLo(myThid),myBxHi(myThid)
259 eta_qrNtile = 0. _d 0
260 DO K=1,1
261 DO J=1-1,sNy+1
262 DO I=1-1,sNx+1
263 cg3d_q(I,J,K,bi,bj) =
264 & zMC(I ,J ,K,bi,bj)*cg3d_r(I ,J ,K,bi,bj)
265 ENDDO
266 ENDDO
267 ENDDO
268 DO K=2,Nr
269 DO J=1-1,sNy+1
270 DO I=1-1,sNx+1
271 cg3d_q(I,J,K,bi,bj) =
272 & zMC(I,J,K,bi,bj)*(cg3d_r(I,J,K ,bi,bj)
273 & -zML(I,J,K,bi,bj)*cg3d_q(I,J,K-1,bi,bj))
274 ENDDO
275 ENDDO
276 ENDDO
277 DO K=Nr,Nr
278 caja IF (Nr .GT. 1) THEN
279 caja DO J=1-1,sNy+1
280 caja DO I=1-1,sNx+1
281 caja cg3d_q(I,J,K,bi,bj) =
282 caja & zMC(i,j,k,bi,bj)*(cg3d_r(i,j,k ,bi,bj)
283 caja & -zML(i,j,k,bi,bj)*cg3d_q(i,j,k-1,bi,bj))
284 caja ENDDO
285 caja ENDDO
286 caja ENDIF
287 DO J=1,sNy
288 DO I=1,sNx
289 eta_qrNtile = eta_qrNtile
290 & +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
291 ENDDO
292 ENDDO
293 ENDDO
294 DO K=Nr-1,1,-1
295 DO J=1-1,sNy+1
296 DO I=1-1,sNx+1
297 cg3d_q(I,J,K,bi,bj) =
298 & cg3d_q(I,J,K,bi,bj)
299 & -zMU(I,J,K,bi,bj)*cg3d_q(I,J,K+1,bi,bj)
300 ENDDO
301 ENDDO
302 DO J=1,sNy
303 DO I=1,sNx
304 eta_qrNtile = eta_qrNtile
305 & +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
306 ENDDO
307 ENDDO
308 ENDDO
309 eta_qrN = eta_qrN + eta_qrNtile
310 ENDDO
311 ENDDO
312 caja
313 caja eta_qrN=0.
314 caja DO bj=myByLo(myThid),myByHi(myThid)
315 caja DO bi=myBxLo(myThid),myBxHi(myThid)
316 caja DO K=1,Nr
317 caja DO J=1,sNy
318 caja DO I=1,sNx
319 caja eta_qrN = eta_qrN
320 caja & +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
321 caja ENDDO
322 caja ENDDO
323 caja ENDDO
324 caja ENDDO
325 caja ENDDO
326 caja
327
328 _GLOBAL_SUM_R8(eta_qrN, myThid)
329 CcnhDebugStarts
330 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' eta_qrN = ',eta_qrN
331 CcnhDebugEnds
332 cgBeta = eta_qrN/eta_qrNM1
333 CcnhDebugStarts
334 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' beta = ',cgBeta
335 CcnhDebugEnds
336 eta_qrNM1 = eta_qrN
337
338 DO bj=myByLo(myThid),myByHi(myThid)
339 DO bi=myBxLo(myThid),myBxHi(myThid)
340 DO K=1,Nr
341 DO J=1-1,sNy+1
342 DO I=1-1,sNx+1
343 cg3d_s(I,J,K,bi,bj) = cg3d_q(I,J,K,bi,bj)
344 & + cgBeta*cg3d_s(I,J,K,bi,bj)
345 ENDDO
346 ENDDO
347 ENDDO
348 ENDDO
349 ENDDO
350
351 C== Evaluate laplace operator on conjugate gradient vector
352 C== q = A.s
353 alpha = 0. _d 0
354 topLevTerm = freeSurfFac*cg3dNorm*
355 & (horiVertRatio/gravity)/deltaTMom/deltaTMom
356 DO bj=myByLo(myThid),myByHi(myThid)
357 DO bi=myBxLo(myThid),myBxHi(myThid)
358 alphaTile = 0. _d 0
359 IF ( Nr .GT. 1 ) THEN
360 DO K=1,1
361 DO J=1,sNy
362 DO I=1,sNx
363 cg3d_q(I,J,K,bi,bj) =
364 & aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj)
365 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj)
366 & +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj)
367 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj)
368 & +aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K+1,bi,bj)
369 & -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
370 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
371 & -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
372 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
373 & -aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
374 & -topLevTerm*_rA(I,J,bi,bj)*cg3d_s(I,J,K,bi,bj)
375 alphaTile = alphaTile
376 & +cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj)
377 ENDDO
378 ENDDO
379 ENDDO
380 ELSE
381 DO K=1,1
382 DO J=1,sNy
383 DO I=1,sNx
384 cg3d_q(I,J,K,bi,bj) =
385 & aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj)
386 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj)
387 & +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj)
388 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj)
389 & -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
390 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
391 & -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
392 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
393 & -topLevTerm*_rA(I,J,bi,bj)*cg3d_s(I,J,K,bi,bj)
394 alphaTile = alphaTile
395 & +cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj)
396 ENDDO
397 ENDDO
398 ENDDO
399 ENDIF
400 DO K=2,Nr-1
401 DO J=1,sNy
402 DO I=1,sNx
403 cg3d_q(I,J,K,bi,bj) =
404 & aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj)
405 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj)
406 & +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj)
407 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj)
408 & +aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K-1,bi,bj)
409 & +aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K+1,bi,bj)
410 & -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
411 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
412 & -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
413 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
414 & -aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
415 & -aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
416 alphaTile = alphaTile
417 & +cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj)
418 ENDDO
419 ENDDO
420 ENDDO
421 IF ( Nr .GT. 1 ) THEN
422 DO K=Nr,Nr
423 DO J=1,sNy
424 DO I=1,sNx
425 cg3d_q(I,J,K,bi,bj) =
426 & aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj)
427 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj)
428 & +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj)
429 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj)
430 & +aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K-1,bi,bj)
431 & -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
432 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
433 & -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
434 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
435 & -aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
436 alphaTile = alphaTile
437 & +cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj)
438 ENDDO
439 ENDDO
440 ENDDO
441 ENDIF
442 alpha = alpha + alphaTile
443 ENDDO
444 ENDDO
445 _GLOBAL_SUM_R8(alpha,myThid)
446 CcnhDebugStarts
447 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' SUM(s*q)= ',alpha
448 CcnhDebugEnds
449 alpha = eta_qrN/alpha
450 CcnhDebugStarts
451 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' alpha= ',alpha
452 CcnhDebugEnds
453
454 C== Update solution and residual vectors
455 C Now compute "interior" points.
456 err = 0. _d 0
457 DO bj=myByLo(myThid),myByHi(myThid)
458 DO bi=myBxLo(myThid),myBxHi(myThid)
459 errTile = 0. _d 0
460 DO K=1,Nr
461 DO J=1,sNy
462 DO I=1,sNx
463 cg3d_x(I,J,K,bi,bj)=cg3d_x(I,J,K,bi,bj)
464 & +alpha*cg3d_s(I,J,K,bi,bj)
465 cg3d_r(I,J,K,bi,bj)=cg3d_r(I,J,K,bi,bj)
466 & -alpha*cg3d_q(I,J,K,bi,bj)
467 errTile = errTile
468 & +cg3d_r(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
469 ENDDO
470 ENDDO
471 ENDDO
472 err = err + errTile
473 ENDDO
474 ENDDO
475
476 _GLOBAL_SUM_R8( err , myThid )
477 err = SQRT(err)
478 actualIts = it3d
479 actualResidual = err
480 IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11
481 C _EXCH_XYZ_R8(cg3d_r, myThid )
482 OLw = 1
483 OLe = 1
484 OLn = 1
485 OLs = 1
486 exchWidthX = 1
487 exchWidthY = 1
488 myNz = Nr
489 CALL EXCH_RL( cg3d_r,
490 I OLw, OLe, OLs, OLn, myNz,
491 I exchWidthX, exchWidthY,
492 I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid )
493
494 10 CONTINUE
495 11 CONTINUE
496
497 C-- Un-normalise the answer
498 DO bj=myByLo(myThid),myByHi(myThid)
499 DO bi=myBxLo(myThid),myBxHi(myThid)
500 DO K=1,Nr
501 DO J=1,sNy
502 DO I=1,sNx
503 cg3d_x(I,J,K,bi,bj) = cg3d_x(I,J,K,bi,bj)/rhsNorm
504 ENDDO
505 ENDDO
506 ENDDO
507 ENDDO
508 ENDDO
509
510 Cadj _EXCH_XYZ_R8(cg3d_x, myThid )
511 c _BEGIN_MASTER( myThid )
512 c WRITE(*,'(A,I6,1PE30.14)') ' CG3D iters, err = ',
513 c & actualIts, actualResidual
514 c _END_MASTER( myThid )
515 lastResidual=actualResidual
516 numIters=actualIts
517
518 #endif /* ALLOW_NONHYDROSTATIC */
519
520 RETURN
521 END

  ViewVC Help
Powered by ViewVC 1.1.22