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

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

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


Revision 1.33 - (hide annotations) (download)
Tue May 29 14:01:36 2001 UTC (23 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre3, checkpoint40pre1, checkpoint40pre7, checkpoint40pre6, checkpoint40pre9, checkpoint40pre8, checkpoint40pre2, checkpoint40pre4, checkpoint40pre5, checkpoint40
Changes since 1.32: +73 -39 lines
Merge from branch pre38:
 o essential mods for cubed sphere
 o debugged atmosphere, dynamcis + physics (aim)
 o new packages (mom_vecinv, mom_fluxform, ...)

1 adcroft 1.33 C $Header: /u/gcmpack/models/MITgcmUV/model/src/cg2d.F,v 1.30.2.3 2001/04/06 13:06:10 jmc Exp $
2     C $Name: $
3 cnh 1.1
4 cnh 1.16 #include "CPP_OPTIONS.h"
5 cnh 1.1
6     SUBROUTINE CG2D(
7 cnh 1.14 I cg2d_b,
8     U cg2d_x,
9 adcroft 1.33 O firstResidual,
10     O lastResidual,
11 adcroft 1.30 U numIters,
12 cnh 1.1 I myThid )
13     C /==========================================================\
14     C | SUBROUTINE CG2D |
15     C | o Two-dimensional grid problem conjugate-gradient |
16     C | inverter (with preconditioner). |
17     C |==========================================================|
18     C | Con. grad is an iterative procedure for solving Ax = b. |
19     C | It requires the A be symmetric. |
20     C | This implementation assumes A is a five-diagonal |
21     C | matrix of the form that arises in the discrete |
22     C | representation of the del^2 operator in a |
23     C | two-dimensional space. |
24     C | Notes: |
25     C | ====== |
26     C | This implementation can support shared-memory |
27     C | multi-threaded execution. In order to do this COMMON |
28     C | blocks are used for many of the arrays - even ones that |
29     C | are only used for intermedaite results. This design is |
30     C | OK if you want to all the threads to collaborate on |
31     C | solving the same problem. On the other hand if you want |
32     C | the threads to solve several different problems |
33     C | concurrently this implementation will not work. |
34     C \==========================================================/
35 adcroft 1.18 IMPLICIT NONE
36 cnh 1.1
37     C === Global data ===
38     #include "SIZE.h"
39     #include "EEPARAMS.h"
40     #include "PARAMS.h"
41 cnh 1.4 #include "GRID.h"
42 adcroft 1.33 #include "CG2D.h"
43 jmc 1.29 #include "SURFACE.h"
44 cnh 1.1
45     C === Routine arguments ===
46 adcroft 1.30 C myThid - Thread on which I am working.
47     C cg2d_b - The source term or "right hand side"
48     C cg2d_x - The solution
49 adcroft 1.33 C firstResidual - the initial residual before any iterations
50     C lastResidual - the actual residual reached
51 adcroft 1.30 C numIters - Entry: the maximum number of iterations allowed
52     C Exit: the actual number of iterations used
53     _RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
54     _RL cg2d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
55 adcroft 1.33 _RL firstResidual
56     _RL lastResidual
57 adcroft 1.30 INTEGER numIters
58 cnh 1.1 INTEGER myThid
59    
60     C === Local variables ====
61     C actualIts - Number of iterations taken
62     C actualResidual - residual
63     C bi - Block index in X and Y.
64     C bj
65 jmc 1.28 C eta_qrN - Used in computing search directions
66     C eta_qrNM1 suffix N and NM1 denote current and
67 cnh 1.1 C cgBeta previous iterations respectively.
68     C alpha
69     C sumRHS - Sum of right-hand-side. Sometimes this is a
70     C useful debuggin/trouble shooting diagnostic.
71     C For neumann problems sumRHS needs to be ~0.
72     C or they converge at a non-zero residual.
73     C err - Measure of residual of Ax - b, usually the norm.
74     C I, J, N - Loop counters ( N counts CG iterations )
75     INTEGER actualIts
76 adcroft 1.33 _RL actualResidual
77 cnh 1.1 INTEGER bi, bj
78     INTEGER I, J, it2d
79 cnh 1.14 _RL err
80 jmc 1.28 _RL eta_qrN
81     _RL eta_qrNM1
82 cnh 1.14 _RL cgBeta
83     _RL alpha
84     _RL sumRHS
85     _RL rhsMax
86     _RL rhsNorm
87 cnh 1.1
88 cnh 1.13 INTEGER OLw
89     INTEGER OLe
90     INTEGER OLn
91     INTEGER OLs
92     INTEGER exchWidthX
93     INTEGER exchWidthY
94     INTEGER myNz
95    
96    
97 cnh 1.12 CcnhDebugStarts
98 adcroft 1.24 C CHARACTER*(MAX_LEN_FNAM) suff
99 cnh 1.12 CcnhDebugEnds
100    
101    
102 cnh 1.1 C-- Initialise inverter
103 jmc 1.28 eta_qrNM1 = 1. _d 0
104 cnh 1.1
105 cnh 1.10 CcnhDebugStarts
106 cnh 1.11 C _EXCH_XY_R8( cg2d_b, myThid )
107     C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid )
108 cnh 1.12 C suff = 'unnormalised'
109     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
110 cnh 1.14 C STOP
111 cnh 1.10 CcnhDebugEnds
112    
113 cnh 1.1 C-- Normalise RHS
114     rhsMax = 0. _d 0
115     DO bj=myByLo(myThid),myByHi(myThid)
116     DO bi=myBxLo(myThid),myBxHi(myThid)
117     DO J=1,sNy
118     DO I=1,sNx
119     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm
120     rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax)
121     ENDDO
122     ENDDO
123     ENDDO
124     ENDDO
125 adcroft 1.33
126     IF (cg2dNormaliseRHS) THEN
127     C- Normalise RHS :
128 adcroft 1.23 #ifdef LETS_MAKE_JAM
129 adcroft 1.25 C _GLOBAL_MAX_R8( rhsMax, myThid )
130     rhsMax=1.
131 adcroft 1.23 #else
132     _GLOBAL_MAX_R8( rhsMax, myThid )
133 adcroft 1.26 Catm rhsMax=1.
134 adcroft 1.23 #endif
135 cnh 1.1 rhsNorm = 1. _d 0
136     IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
137     DO bj=myByLo(myThid),myByHi(myThid)
138     DO bi=myBxLo(myThid),myBxHi(myThid)
139     DO J=1,sNy
140     DO I=1,sNx
141     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
142     cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
143     ENDDO
144     ENDDO
145     ENDDO
146     ENDDO
147 adcroft 1.33 C- end Normalise RHS
148     ENDIF
149 cnh 1.1
150     C-- Update overlaps
151     _EXCH_XY_R8( cg2d_b, myThid )
152     _EXCH_XY_R8( cg2d_x, myThid )
153     CcnhDebugStarts
154 cnh 1.11 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
155 cnh 1.12 C suff = 'normalised'
156     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
157 cnh 1.1 CcnhDebugEnds
158    
159     C-- Initial residual calculation
160     err = 0. _d 0
161     sumRHS = 0. _d 0
162     DO bj=myByLo(myThid),myByHi(myThid)
163     DO bi=myBxLo(myThid),myBxHi(myThid)
164     DO J=1,sNy
165     DO I=1,sNx
166     cg2d_s(I,J,bi,bj) = 0.
167     cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
168     & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
169     & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
170     & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
171     & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
172     & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
173     & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
174     & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
175 cnh 1.4 & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
176 jmc 1.29 & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
177 cnh 1.4 & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
178     & )
179 cnh 1.1 err = err +
180     & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
181     sumRHS = sumRHS +
182     & cg2d_b(I,J,bi,bj)
183     ENDDO
184     ENDDO
185     ENDDO
186     ENDDO
187 cnh 1.13 C _EXCH_XY_R8( cg2d_r, myThid )
188 adcroft 1.23 #ifdef LETS_MAKE_JAM
189     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
190     #else
191 cnh 1.13 OLw = 1
192     OLe = 1
193     OLn = 1
194     OLs = 1
195     exchWidthX = 1
196     exchWidthY = 1
197     myNz = 1
198 adcroft 1.33 IF (useCubedSphereExchange) THEN
199     CALL EXCH_RL_CUBE( cg2d_r,
200     I OLw, OLe, OLs, OLn, myNz,
201     I exchWidthX, exchWidthY,
202     I FORWARD_SIMULATION, EXCH_UPDATE_CORNERS, myThid )
203     ELSE
204 cnh 1.13 CALL EXCH_RL( cg2d_r,
205 adcroft 1.33 I OLw, OLe, OLs, OLn, myNz,
206 cnh 1.13 I exchWidthX, exchWidthY,
207 adcroft 1.33 I FORWARD_SIMULATION, EXCH_UPDATE_CORNERS, myThid )
208     ENDIF
209 adcroft 1.23 #endif
210 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
211 adcroft 1.23 #ifdef LETS_MAKE_JAM
212     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
213     #else
214 cnh 1.13 OLw = 1
215     OLe = 1
216     OLn = 1
217     OLs = 1
218     exchWidthX = 1
219     exchWidthY = 1
220     myNz = 1
221 adcroft 1.33 IF (useCubedSphereExchange) THEN
222     CALL EXCH_RL_CUBE( cg2d_s,
223     I OLw, OLe, OLs, OLn, myNz,
224     I exchWidthX, exchWidthY,
225     I FORWARD_SIMULATION, EXCH_UPDATE_CORNERS, myThid )
226     ELSE
227 cnh 1.13 CALL EXCH_RL( cg2d_s,
228 adcroft 1.33 I OLw, OLe, OLs, OLn, myNz,
229 cnh 1.13 I exchWidthX, exchWidthY,
230 adcroft 1.33 I FORWARD_SIMULATION, EXCH_UPDATE_CORNERS, myThid )
231     ENDIF
232 adcroft 1.23 #endif
233 adcroft 1.33 _GLOBAL_SUM_R8( sumRHS, myThid )
234     _GLOBAL_SUM_R8( err , myThid )
235     err = SQRT(err)
236     actualIts = 0
237     actualResidual = err
238 cnh 1.13
239     _BEGIN_MASTER( myThid )
240 adcroft 1.33 write(*,'(A,1P2E22.14)')' cg2d: Sum(rhs),rhsMax = ',
241     & sumRHS,rhsMax
242 cnh 1.13 _END_MASTER( )
243 cnh 1.1 C _BARRIER
244 adcroft 1.30 c _BEGIN_MASTER( myThid )
245 heimbach 1.31 c WRITE(*,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
246 adcroft 1.30 c & actualIts, actualResidual
247     c _END_MASTER( )
248 adcroft 1.33 firstResidual=actualResidual
249    
250     IF ( err .LT. cg2dTolerance ) GOTO 11
251 cnh 1.1
252     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
253 adcroft 1.30 DO 10 it2d=1, numIters
254 cnh 1.1
255     CcnhDebugStarts
256 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' residual = ',
257 cnh 1.14 C & actualResidual
258 cnh 1.1 CcnhDebugEnds
259     C-- Solve preconditioning equation and update
260     C-- conjugate direction vector "s".
261 jmc 1.28 eta_qrN = 0. _d 0
262 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
263     DO bi=myBxLo(myThid),myBxHi(myThid)
264     DO J=1,sNy
265     DO I=1,sNx
266     cg2d_q(I,J,bi,bj) =
267 cnh 1.3 & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
268     & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
269     & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
270     & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
271     & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
272 cnh 1.4 CcnhDebugStarts
273     C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
274     CcnhDebugEnds
275 jmc 1.28 eta_qrN = eta_qrN
276 cnh 1.1 & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
277     ENDDO
278     ENDDO
279     ENDDO
280     ENDDO
281    
282 jmc 1.28 _GLOBAL_SUM_R8(eta_qrN, myThid)
283 cnh 1.1 CcnhDebugStarts
284 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
285 cnh 1.1 CcnhDebugEnds
286 jmc 1.28 cgBeta = eta_qrN/eta_qrNM1
287 cnh 1.1 CcnhDebugStarts
288 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
289 cnh 1.1 CcnhDebugEnds
290 jmc 1.28 eta_qrNM1 = eta_qrN
291 cnh 1.1
292     DO bj=myByLo(myThid),myByHi(myThid)
293     DO bi=myBxLo(myThid),myBxHi(myThid)
294     DO J=1,sNy
295     DO I=1,sNx
296 cnh 1.14 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj)
297     & + cgBeta*cg2d_s(I,J,bi,bj)
298 cnh 1.1 ENDDO
299     ENDDO
300     ENDDO
301     ENDDO
302    
303     C-- Do exchanges that require messages i.e. between
304     C-- processes.
305 cnh 1.13 C _EXCH_XY_R8( cg2d_s, myThid )
306 adcroft 1.23 #ifdef LETS_MAKE_JAM
307     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
308     #else
309 cnh 1.13 OLw = 1
310     OLe = 1
311     OLn = 1
312     OLs = 1
313     exchWidthX = 1
314     exchWidthY = 1
315     myNz = 1
316 adcroft 1.33 IF (useCubedSphereExchange) THEN
317     CALL EXCH_RL_CUBE( cg2d_s,
318     I OLw, OLe, OLs, OLn, myNz,
319     I exchWidthX, exchWidthY,
320     I FORWARD_SIMULATION, EXCH_UPDATE_CORNERS, myThid )
321     ELSE
322 cnh 1.13 CALL EXCH_RL( cg2d_s,
323 adcroft 1.33 I OLw, OLe, OLs, OLn, myNz,
324 cnh 1.13 I exchWidthX, exchWidthY,
325 adcroft 1.33 I FORWARD_SIMULATION, EXCH_UPDATE_CORNERS, myThid )
326     ENDIF
327 adcroft 1.23 #endif
328 cnh 1.1
329     C== Evaluate laplace operator on conjugate gradient vector
330     C== q = A.s
331     alpha = 0. _d 0
332     DO bj=myByLo(myThid),myByHi(myThid)
333     DO bi=myBxLo(myThid),myBxHi(myThid)
334     DO J=1,sNy
335     DO I=1,sNx
336     cg2d_q(I,J,bi,bj) =
337     & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
338     & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
339     & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
340     & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
341     & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
342     & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
343     & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
344     & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
345 jmc 1.29 & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
346 cnh 1.4 & cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
347 cnh 1.1 alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
348     ENDDO
349     ENDDO
350     ENDDO
351     ENDDO
352 adcroft 1.20 _GLOBAL_SUM_R8(alpha,myThid)
353 cnh 1.1 CcnhDebugStarts
354 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
355 cnh 1.1 CcnhDebugEnds
356 jmc 1.28 alpha = eta_qrN/alpha
357 cnh 1.1 CcnhDebugStarts
358 heimbach 1.31 C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
359 cnh 1.1 CcnhDebugEnds
360    
361     C== Update solution and residual vectors
362     C Now compute "interior" points.
363     err = 0. _d 0
364     DO bj=myByLo(myThid),myByHi(myThid)
365     DO bi=myBxLo(myThid),myBxHi(myThid)
366     DO J=1,sNy
367     DO I=1,sNx
368     cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
369     cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
370     err = err+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
371     ENDDO
372     ENDDO
373     ENDDO
374     ENDDO
375    
376 adcroft 1.20 _GLOBAL_SUM_R8( err , myThid )
377 cnh 1.1 err = SQRT(err)
378     actualIts = it2d
379     actualResidual = err
380 adcroft 1.33 IF ( err .LT. cg2dTolerance ) GOTO 11
381 cnh 1.13 C _EXCH_XY_R8(cg2d_r, myThid )
382 adcroft 1.23 #ifdef LETS_MAKE_JAM
383     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
384     #else
385 cnh 1.13 OLw = 1
386     OLe = 1
387     OLn = 1
388     OLs = 1
389     exchWidthX = 1
390     exchWidthY = 1
391     myNz = 1
392 adcroft 1.33 IF (useCubedSphereExchange) THEN
393     CALL EXCH_RL_CUBE( cg2d_r,
394     I OLw, OLe, OLs, OLn, myNz,
395     I exchWidthX, exchWidthY,
396     I FORWARD_SIMULATION, EXCH_UPDATE_CORNERS, myThid )
397     ELSE
398 cnh 1.13 CALL EXCH_RL( cg2d_r,
399 adcroft 1.33 I OLw, OLe, OLs, OLn, myNz,
400     I exchWidthX, exchWidthY,
401     I FORWARD_SIMULATION, EXCH_UPDATE_CORNERS, myThid )
402     ENDIF
403 adcroft 1.23 #endif
404 cnh 1.13
405 cnh 1.1 10 CONTINUE
406     11 CONTINUE
407    
408 adcroft 1.33 IF (cg2dNormaliseRHS) THEN
409 cnh 1.1 C-- Un-normalise the answer
410 adcroft 1.33 DO bj=myByLo(myThid),myByHi(myThid)
411     DO bi=myBxLo(myThid),myBxHi(myThid)
412     DO J=1,sNy
413     DO I=1,sNx
414     cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
415     ENDDO
416     ENDDO
417 cnh 1.1 ENDDO
418     ENDDO
419 adcroft 1.33 ENDIF
420 cnh 1.1
421 adcroft 1.22 C The following exchange was moved up to solve_for_pressure
422     C for compatibility with TAMC.
423     C _EXCH_XY_R8(cg2d_x, myThid )
424 adcroft 1.30 c _BEGIN_MASTER( myThid )
425 heimbach 1.31 c WRITE(*,'(A,I6,1PE30.14)') ' CG2D iters, err = ',
426 adcroft 1.30 c & actualIts, actualResidual
427     c _END_MASTER( )
428    
429     C-- Return parameters to caller
430 adcroft 1.33 lastResidual=actualResidual
431 adcroft 1.30 numIters=actualIts
432 cnh 1.1
433     CcnhDebugStarts
434 cnh 1.7 C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
435 cnh 1.1 C err = 0. _d 0
436     C DO bj=myByLo(myThid),myByHi(myThid)
437     C DO bi=myBxLo(myThid),myBxHi(myThid)
438     C DO J=1,sNy
439     C DO I=1,sNx
440     C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
441     C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
442     C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
443     C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
444     C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
445     C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
446     C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
447     C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
448     C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
449     C err = err +
450     C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
451     C ENDDO
452     C ENDDO
453     C ENDDO
454     C ENDDO
455 adcroft 1.20 C _GLOBAL_SUM_R8( err , myThid )
456 heimbach 1.31 C write(*,*) 'cg2d: Ax - b = ',SQRT(err)
457 cnh 1.1 CcnhDebugEnds
458    
459 adcroft 1.19 RETURN
460 cnh 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22