/[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.1 - (hide annotations) (download)
Mon Mar 22 15:54:03 1999 UTC (25 years, 2 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint20, checkpoint21
Modifications for non-hydrostatic ability + updates for open-boundaries.

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

  ViewVC Help
Powered by ViewVC 1.1.22