/[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.4 - (show annotations) (download)
Mon Mar 27 22:25:43 2000 UTC (24 years, 2 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint28, checkpoint29, checkpoint27, branch-atmos-merge-freeze, branch-atmos-merge-start, checkpoint26, branch-atmos-merge-shapiro, checkpoint33, checkpoint32, checkpoint31, checkpoint30, checkpoint34, branch-atmos-merge-zonalfilt, branch-atmos-merge-phase5, branch-atmos-merge-phase4, branch-atmos-merge-phase7, branch-atmos-merge-phase6, branch-atmos-merge-phase1, branch-atmos-merge-phase3, branch-atmos-merge-phase2
Branch point for: branch-atmos-merge
Changes since 1.3: +3 -8 lines
Removed unused variables and fixed some unitialized variables.

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

  ViewVC Help
Powered by ViewVC 1.1.22