/[MITgcm]/MITgcm/model/src/cg3d.F
ViewVC logotype

Annotation of /MITgcm/model/src/cg3d.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2 - (hide annotations) (download)
Tue May 18 17:36:55 1999 UTC (25 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint22
Changes since 1.1: +10 -31 lines
Changed number of arguments to GLOBAL_SUM and GLOBAL_MAX to two.
Instigated by Ralf.

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

  ViewVC Help
Powered by ViewVC 1.1.22