/[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.5 - (show annotations) (download)
Sun Feb 4 14:38:46 2001 UTC (23 years, 4 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint36, checkpoint35
Changes since 1.4: +2 -1 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

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

  ViewVC Help
Powered by ViewVC 1.1.22