/[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.11 - (show annotations) (download)
Fri Jun 29 17:14:49 2001 UTC (22 years, 11 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre3, checkpoint40pre1, checkpoint40pre7, checkpoint40pre6, checkpoint40pre9, checkpoint40pre8, checkpoint40pre2, checkpoint40pre4, checkpoint40pre5, checkpoint40
Changes since 1.10: +35 -14 lines
Moved cg3d_x into DYNVARS.h and renamed it to phi_nh.
 - cg3d and cg2d now look more similar
 - output formatted to fit Chris's tastes (I think)

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

  ViewVC Help
Powered by ViewVC 1.1.22