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