1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/cg3d.F,v 1.9 2001/05/14 21:33:30 heimbach 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 eta_qrN - Used in computing search directions |
53 |
C eta_qrNM1 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 eta_qrN |
69 |
_RL eta_qrNM1 |
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 topLevTerm |
84 |
|
85 |
C-- Initialise inverter |
86 |
eta_qrNM1 = 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 |
C-- Initial residual calculation (with free-Surface term) |
123 |
err = 0. _d 0 |
124 |
sumRHS = 0. _d 0 |
125 |
DO bj=myByLo(myThid),myByHi(myThid) |
126 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
127 |
DO K=1,Nr |
128 |
KM1 = K-1 |
129 |
IF ( K .EQ. 1 ) KM1 = 1 |
130 |
KP1 = K+1 |
131 |
IF ( K .EQ. Nr ) KP1 = 1 |
132 |
topLevTerm = 0. |
133 |
IF ( K .EQ. 1) topLevTerm = freeSurfFac*cg3dNorm* |
134 |
& (horiVertRatio/gravity)/deltaTMom/deltaTMom |
135 |
DO J=1,sNy |
136 |
DO I=1,sNx |
137 |
cg3d_s(I,J,K,bi,bj) = 0. |
138 |
cg3d_r(I,J,K,bi,bj) = cg3d_b(I,J,K,bi,bj) -( 0. |
139 |
& +aW3d(I ,J ,K ,bi,bj)*cg3d_x(I-1,J ,K ,bi,bj) |
140 |
& +aW3d(I+1,J ,K ,bi,bj)*cg3d_x(I+1,J ,K ,bi,bj) |
141 |
& +aS3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J-1,K ,bi,bj) |
142 |
& +aS3d(I ,J+1,K ,bi,bj)*cg3d_x(I ,J+1,K ,bi,bj) |
143 |
& +aV3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,KM1,bi,bj) |
144 |
& +aV3d(I ,J ,KP1,bi,bj)*cg3d_x(I ,J ,KP1,bi,bj) |
145 |
& -aW3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj) |
146 |
& -aW3d(I+1,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj) |
147 |
& -aS3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj) |
148 |
& -aS3d(I ,J+1,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj) |
149 |
& -aV3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj) |
150 |
& -aV3d(I ,J ,KP1,bi,bj)*cg3d_x(I ,J ,K ,bi,bj) |
151 |
& -topLevTerm*_rA(I,J,bi,bj)*cg3d_x(I,J,K,bi,bj) |
152 |
& ) |
153 |
err = err |
154 |
& +cg3d_r(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj) |
155 |
sumRHS = sumRHS |
156 |
& +cg3d_b(I,J,K,bi,bj) |
157 |
ENDDO |
158 |
ENDDO |
159 |
ENDDO |
160 |
ENDDO |
161 |
ENDDO |
162 |
C _EXCH_XYZ_R8( cg3d_r, myThid ) |
163 |
OLw = 1 |
164 |
OLe = 1 |
165 |
OLn = 1 |
166 |
OLs = 1 |
167 |
exchWidthX = 1 |
168 |
exchWidthY = 1 |
169 |
myNz = Nr |
170 |
CALL EXCH_RL( cg3d_r, |
171 |
I OLw, OLe, OLs, OLn, myNz, |
172 |
I exchWidthX, exchWidthY, |
173 |
I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid ) |
174 |
C _EXCH_XYZ_R8( cg3d_s, myThid ) |
175 |
OLw = 1 |
176 |
OLe = 1 |
177 |
OLn = 1 |
178 |
OLs = 1 |
179 |
exchWidthX = 1 |
180 |
exchWidthY = 1 |
181 |
myNz = Nr |
182 |
CALL EXCH_RL( cg3d_s, |
183 |
I OLw, OLe, OLs, OLn, myNz, |
184 |
I exchWidthX, exchWidthY, |
185 |
I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid ) |
186 |
_GLOBAL_SUM_R8( sumRHS, myThid ) |
187 |
_GLOBAL_SUM_R8( err , myThid ) |
188 |
|
189 |
_BEGIN_MASTER( myThid ) |
190 |
write(*,'(A,1PE30.14)') ' cg3d: Sum(rhs) = ',sumRHS |
191 |
_END_MASTER( ) |
192 |
|
193 |
actualIts = 0 |
194 |
actualResidual = SQRT(err) |
195 |
C _BARRIER |
196 |
_BEGIN_MASTER( myThid ) |
197 |
WRITE(*,'(A,I6,1PE30.14)') ' CG3D iters, err = ', |
198 |
& actualIts, actualResidual |
199 |
_END_MASTER( ) |
200 |
|
201 |
C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
202 |
DO 10 it3d=1, cg3dMaxIters |
203 |
|
204 |
CcnhDebugStarts |
205 |
#ifdef VERBOSE |
206 |
IF ( mod(it3d-1,10).EQ.0) |
207 |
& WRITE(*,*) ' CG3D: Iteration ',it3d-1, |
208 |
& ' residual = ',actualResidual |
209 |
#endif |
210 |
CcnhDebugEnds |
211 |
IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11 |
212 |
C-- Solve preconditioning equation and update |
213 |
C-- conjugate direction vector "s". |
214 |
C Note. On the next to loops over all tiles the inner loop ranges |
215 |
C in sNx and sNy are expanded by 1 to avoid a communication |
216 |
C step. However this entails a bit of gynamastics because we only |
217 |
C want eta_qrN for the interior points. |
218 |
eta_qrN = 0. _d 0 |
219 |
DO bj=myByLo(myThid),myByHi(myThid) |
220 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
221 |
DO K=1,1 |
222 |
DO J=1-1,sNy+1 |
223 |
DO I=1-1,sNx+1 |
224 |
cg3d_q(I,J,K,bi,bj) = |
225 |
& zMC(I ,J ,K,bi,bj)*cg3d_r(I ,J ,K,bi,bj) |
226 |
ENDDO |
227 |
ENDDO |
228 |
ENDDO |
229 |
DO K=2,Nr |
230 |
DO J=1-1,sNy+1 |
231 |
DO I=1-1,sNx+1 |
232 |
cg3d_q(I,J,K,bi,bj) = |
233 |
& zMC(I,J,K,bi,bj)*(cg3d_r(I,J,K ,bi,bj) |
234 |
& -zML(I,J,K,bi,bj)*cg3d_q(I,J,K-1,bi,bj)) |
235 |
ENDDO |
236 |
ENDDO |
237 |
ENDDO |
238 |
DO K=Nr,Nr |
239 |
caja IF (Nr .GT. 1) THEN |
240 |
caja DO J=1-1,sNy+1 |
241 |
caja DO I=1-1,sNx+1 |
242 |
caja cg3d_q(I,J,K,bi,bj) = |
243 |
caja & zMC(i,j,k,bi,bj)*(cg3d_r(i,j,k ,bi,bj) |
244 |
caja & -zML(i,j,k,bi,bj)*cg3d_q(i,j,k-1,bi,bj)) |
245 |
caja ENDDO |
246 |
caja ENDDO |
247 |
caja ENDIF |
248 |
DO J=1,sNy |
249 |
DO I=1,sNx |
250 |
eta_qrN = eta_qrN |
251 |
& +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj) |
252 |
ENDDO |
253 |
ENDDO |
254 |
ENDDO |
255 |
DO K=Nr-1,1,-1 |
256 |
DO J=1-1,sNy+1 |
257 |
DO I=1-1,sNx+1 |
258 |
cg3d_q(I,J,K,bi,bj) = |
259 |
& cg3d_q(I,J,K,bi,bj) |
260 |
& -zMU(I,J,K,bi,bj)*cg3d_q(I,J,K+1,bi,bj) |
261 |
ENDDO |
262 |
ENDDO |
263 |
DO J=1,sNy |
264 |
DO I=1,sNx |
265 |
eta_qrN = eta_qrN |
266 |
& +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj) |
267 |
ENDDO |
268 |
ENDDO |
269 |
ENDDO |
270 |
ENDDO |
271 |
ENDDO |
272 |
caja |
273 |
caja eta_qrN=0. |
274 |
caja DO bj=myByLo(myThid),myByHi(myThid) |
275 |
caja DO bi=myBxLo(myThid),myBxHi(myThid) |
276 |
caja DO K=1,Nr |
277 |
caja DO J=1,sNy |
278 |
caja DO I=1,sNx |
279 |
caja eta_qrN = eta_qrN |
280 |
caja & +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj) |
281 |
caja ENDDO |
282 |
caja ENDDO |
283 |
caja ENDDO |
284 |
caja ENDDO |
285 |
caja ENDDO |
286 |
caja |
287 |
|
288 |
_GLOBAL_SUM_R8(eta_qrN, myThid) |
289 |
CcnhDebugStarts |
290 |
C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' eta_qrN = ',eta_qrN |
291 |
CcnhDebugEnds |
292 |
cgBeta = eta_qrN/eta_qrNM1 |
293 |
CcnhDebugStarts |
294 |
C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' beta = ',cgBeta |
295 |
CcnhDebugEnds |
296 |
eta_qrNM1 = eta_qrN |
297 |
|
298 |
DO bj=myByLo(myThid),myByHi(myThid) |
299 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
300 |
DO K=1,Nr |
301 |
DO J=1-1,sNy+1 |
302 |
DO I=1-1,sNx+1 |
303 |
cg3d_s(I,J,K,bi,bj) = cg3d_q(I,J,K,bi,bj) |
304 |
& + cgBeta*cg3d_s(I,J,K,bi,bj) |
305 |
ENDDO |
306 |
ENDDO |
307 |
ENDDO |
308 |
ENDDO |
309 |
ENDDO |
310 |
|
311 |
C== Evaluate laplace operator on conjugate gradient vector |
312 |
C== q = A.s |
313 |
alpha = 0. _d 0 |
314 |
topLevTerm = freeSurfFac*cg3dNorm* |
315 |
& (horiVertRatio/gravity)/deltaTMom/deltaTMom |
316 |
DO bj=myByLo(myThid),myByHi(myThid) |
317 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
318 |
IF ( Nr .GT. 1 ) THEN |
319 |
DO K=1,1 |
320 |
DO J=1,sNy |
321 |
DO I=1,sNx |
322 |
cg3d_q(I,J,K,bi,bj) = |
323 |
& aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj) |
324 |
& +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj) |
325 |
& +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj) |
326 |
& +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj) |
327 |
& +aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K+1,bi,bj) |
328 |
& -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
329 |
& -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
330 |
& -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
331 |
& -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
332 |
& -aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
333 |
& -topLevTerm*_rA(I,J,bi,bj)*cg3d_s(I,J,K,bi,bj) |
334 |
alpha = alpha+cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj) |
335 |
ENDDO |
336 |
ENDDO |
337 |
ENDDO |
338 |
ELSE |
339 |
DO K=1,1 |
340 |
DO J=1,sNy |
341 |
DO I=1,sNx |
342 |
cg3d_q(I,J,K,bi,bj) = |
343 |
& aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj) |
344 |
& +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj) |
345 |
& +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj) |
346 |
& +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,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 |
& -topLevTerm*_rA(I,J,bi,bj)*cg3d_s(I,J,K,bi,bj) |
352 |
alpha = alpha+cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj) |
353 |
ENDDO |
354 |
ENDDO |
355 |
ENDDO |
356 |
ENDIF |
357 |
DO K=2,Nr-1 |
358 |
DO J=1,sNy |
359 |
DO I=1,sNx |
360 |
cg3d_q(I,J,K,bi,bj) = |
361 |
& aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj) |
362 |
& +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj) |
363 |
& +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj) |
364 |
& +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj) |
365 |
& +aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K-1,bi,bj) |
366 |
& +aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K+1,bi,bj) |
367 |
& -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
368 |
& -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
369 |
& -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
370 |
& -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
371 |
& -aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
372 |
& -aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
373 |
alpha = alpha+cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj) |
374 |
ENDDO |
375 |
ENDDO |
376 |
ENDDO |
377 |
IF ( Nr .GT. 1 ) THEN |
378 |
DO K=Nr,Nr |
379 |
DO J=1,sNy |
380 |
DO I=1,sNx |
381 |
cg3d_q(I,J,K,bi,bj) = |
382 |
& aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj) |
383 |
& +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj) |
384 |
& +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj) |
385 |
& +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj) |
386 |
& +aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K-1,bi,bj) |
387 |
& -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
388 |
& -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
389 |
& -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
390 |
& -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj) |
391 |
& -aV3d(I ,J ,K ,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 |
ENDIF |
397 |
ENDDO |
398 |
ENDDO |
399 |
_GLOBAL_SUM_R8(alpha,myThid) |
400 |
CcnhDebugStarts |
401 |
C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' SUM(s*q)= ',alpha |
402 |
CcnhDebugEnds |
403 |
alpha = eta_qrN/alpha |
404 |
CcnhDebugStarts |
405 |
C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' alpha= ',alpha |
406 |
CcnhDebugEnds |
407 |
|
408 |
C== Update solution and residual vectors |
409 |
C Now compute "interior" points. |
410 |
err = 0. _d 0 |
411 |
DO bj=myByLo(myThid),myByHi(myThid) |
412 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
413 |
DO K=1,Nr |
414 |
DO J=1,sNy |
415 |
DO I=1,sNx |
416 |
cg3d_x(I,J,K,bi,bj)=cg3d_x(I,J,K,bi,bj) |
417 |
& +alpha*cg3d_s(I,J,K,bi,bj) |
418 |
cg3d_r(I,J,K,bi,bj)=cg3d_r(I,J,K,bi,bj) |
419 |
& -alpha*cg3d_q(I,J,K,bi,bj) |
420 |
err = err+cg3d_r(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj) |
421 |
ENDDO |
422 |
ENDDO |
423 |
ENDDO |
424 |
ENDDO |
425 |
ENDDO |
426 |
|
427 |
_GLOBAL_SUM_R8( err , myThid ) |
428 |
err = SQRT(err) |
429 |
actualIts = it3d |
430 |
actualResidual = err |
431 |
IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11 |
432 |
C _EXCH_XYZ_R8(cg3d_r, myThid ) |
433 |
OLw = 1 |
434 |
OLe = 1 |
435 |
OLn = 1 |
436 |
OLs = 1 |
437 |
exchWidthX = 1 |
438 |
exchWidthY = 1 |
439 |
myNz = Nr |
440 |
CALL EXCH_RL( cg3d_r, |
441 |
I OLw, OLe, OLs, OLn, myNz, |
442 |
I exchWidthX, exchWidthY, |
443 |
I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid ) |
444 |
|
445 |
10 CONTINUE |
446 |
11 CONTINUE |
447 |
|
448 |
C-- Un-normalise the answer |
449 |
DO bj=myByLo(myThid),myByHi(myThid) |
450 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
451 |
DO K=1,Nr |
452 |
DO J=1,sNy |
453 |
DO I=1,sNx |
454 |
cg3d_x(I,J,K,bi,bj) = cg3d_x(I,J,K,bi,bj)/rhsNorm |
455 |
ENDDO |
456 |
ENDDO |
457 |
ENDDO |
458 |
ENDDO |
459 |
ENDDO |
460 |
|
461 |
Cadj _EXCH_XYZ_R8(cg3d_x, myThid ) |
462 |
_BEGIN_MASTER( myThid ) |
463 |
WRITE(*,'(A,I6,1PE30.14)') ' CG3D iters, err = ', |
464 |
& actualIts, actualResidual |
465 |
_END_MASTER( ) |
466 |
|
467 |
#endif /* ALLOW_NONHYDROSTATIC */ |
468 |
|
469 |
RETURN |
470 |
END |