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