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