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

Contents of /MITgcm/model/src/cg2d_nsa.F

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


Revision 1.3 - (show annotations) (download)
Sat Jul 11 21:45:44 2009 UTC (14 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63l, checkpoint63m, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint62, checkpoint63, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.2: +117 -122 lines
- import modifs from standard version (cg2d.F, v1.16):
 use pre-computed solver main-diagonal term (stored in common block)
 which affects truncation error.

1 C $Header: /u/gcmpack/MITgcm/model/src/cg2d_nsa.F,v 1.2 2009/04/28 18:01:14 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5 #ifdef ALLOW_USE_MPI
6 C HACK to avoid global_max
7 # define ALLOW_CONST_RHSMAX
8 #endif
9
10 CML THIS DOES NOT WORK +++++
11 #undef ALLOW_LOOP_DIRECTIVE
12 CBOP
13 C !ROUTINE: CG2D_NSA
14 C !INTERFACE:
15 SUBROUTINE CG2D_NSA(
16 I cg2d_b,
17 U cg2d_x,
18 O firstResidual,
19 O lastResidual,
20 U numIters,
21 I myThid )
22 C !DESCRIPTION: \bv
23 C *==========================================================*
24 C | SUBROUTINE CG2D_NSA
25 C | o Two-dimensional grid problem conjugate-gradient
26 C | inverter (with preconditioner).
27 C | o This version is used only in the case when the matrix
28 C | operator is not "self-adjoint" (NSA). Any remaining
29 C | residuals will immediately reported to the department
30 C | of homeland security.
31 C *==========================================================*
32 C | Con. grad is an iterative procedure for solving Ax = b.
33 C | It requires the A be symmetric.
34 C | This implementation assumes A is a five-diagonal
35 C | matrix of the form that arises in the discrete
36 C | representation of the del^2 operator in a
37 C | two-dimensional space.
38 C | Notes:
39 C | ======
40 C | This implementation can support shared-memory
41 C | multi-threaded execution. In order to do this COMMON
42 C | blocks are used for many of the arrays - even ones that
43 C | are only used for intermedaite results. This design is
44 C | OK if you want to all the threads to collaborate on
45 C | solving the same problem. On the other hand if you want
46 C | the threads to solve several different problems
47 C | concurrently this implementation will not work.
48 C *==========================================================*
49 C \ev
50
51 C !USES:
52 IMPLICIT NONE
53 C === Global data ===
54 #include "SIZE.h"
55 #include "EEPARAMS.h"
56 #include "PARAMS.h"
57 #include "CG2D.h"
58 c#include "GRID.h"
59 c#include "SURFACE.h"
60 #ifdef ALLOW_AUTODIFF_TAMC
61 # include "tamc.h"
62 # include "tamc_keys.h"
63 #endif
64
65 C !INPUT/OUTPUT PARAMETERS:
66 C === Routine arguments ===
67 C cg2d_b :: The source term or "right hand side"
68 C cg2d_x :: The solution
69 C firstResidual :: the initial residual before any iterations
70 C lastResidual :: the actual residual reached
71 C numIters :: Entry: the maximum number of iterations allowed
72 C Exit: the actual number of iterations used
73 C myThid :: Thread on which I am working.
74 _RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
75 _RL cg2d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
76 _RL firstResidual
77 _RL lastResidual
78 INTEGER numIters
79 INTEGER myThid
80
81 #ifdef ALLOW_CG2D_NSA
82 C !LOCAL VARIABLES:
83 C === Local variables ====
84 C actualIts :: Number of iterations taken
85 C actualResidual :: residual
86 C bi, bj :: Block index in X and Y.
87 C eta_qrN :: Used in computing search directions
88 C eta_qrNM1 suffix N and NM1 denote current and
89 C cgBeta previous iterations respectively.
90 C recip_eta_qrNM1 :: reciprocal of eta_qrNM1
91 C alpha
92 C alpha_aux :: to avoid the statement: alpha = 1./alpha (for TAMC/TAF)
93 C sumRHS :: Sum of right-hand-side. Sometimes this is a
94 C useful debuggin/trouble shooting diagnostic.
95 C For neumann problems sumRHS needs to be ~0.
96 C or they converge at a non-zero residual.
97 C err :: Measure of residual of Ax - b, usually the norm.
98 C err_sq :: square of err (for TAMC/TAF)
99 C I, J, it2d :: Loop counters ( it2d counts CG iterations )
100 INTEGER actualIts
101 _RL actualResidual
102 INTEGER bi, bj
103 INTEGER I, J, it2d
104
105 _RL err
106 _RL err_sq
107 _RL eta_qrN
108 _RL eta_qrNM1
109 _RL recip_eta_qrNM1
110 _RL cgBeta
111 _RL alpha
112 _RL alpha_aux
113 _RL sumRHS
114 _RL rhsMax, rhsMaxGlobal
115 _RL rhsNorm
116 _RL cg2dTolerance_sq
117 CEOP
118
119 #ifdef ALLOW_AUTODIFF_TAMC
120 IF ( numIters .GT. numItersMax ) THEN
121 WRITE(standardMessageUnit,'(A,I10)')
122 & 'CG2D_NSA: numIters > numItersMax = ', numItersMax
123 STOP 'NON-NORMAL in CG2D_NSA'
124 ENDIF
125 #endif /* ALLOW_AUTODIFF_TAMC */
126
127 CcnhDebugStarts
128 C CHARACTER*(MAX_LEN_FNAM) suff
129 CcnhDebugEnds
130
131 #ifdef ALLOW_AUTODIFF_TAMC
132 act1 = myThid - 1
133 max1 = nTx*nTy
134 act2 = ikey_dynamics - 1
135 ikey = (act1 + 1) + act2*max1
136 #endif /* ALLOW_AUTODIFF_TAMC */
137
138 C-- Initialise inverter
139 eta_qrNM1 = 1. _d 0
140 recip_eta_qrNM1 = 1./eta_qrNM1
141
142 CcnhDebugStarts
143 C _EXCH_XY_RL( cg2d_b, myThid )
144 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid )
145 C suff = 'unnormalised'
146 C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
147 C STOP
148 CcnhDebugEnds
149
150 C-- Normalise RHS
151 #ifdef ALLOW_AUTODIFF_TAMC
152 CADJ STORE cg2d_b = comlev1_cg2d, key = ikey, byte = isbyte
153 #endif /* ALLOW_AUTODIFF_TAMC */
154 rhsMax = 0. _d 0
155 DO bj=myByLo(myThid),myByHi(myThid)
156 DO bi=myBxLo(myThid),myBxHi(myThid)
157 DO J=1,sNy
158 DO I=1,sNx
159 cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm
160 rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax)
161 ENDDO
162 ENDDO
163 ENDDO
164 ENDDO
165
166 IF (cg2dNormaliseRHS) THEN
167 C - Normalise RHS :
168 #ifdef LETS_MAKE_JAM
169 C _GLOBAL_MAX_RL( rhsMax, myThid )
170 rhsMaxGlobal=1.
171 #else
172 #ifdef ALLOW_CONST_RHSMAX
173 rhsMaxGlobal=1.
174 #else
175 rhsMaxGlobal=rhsMax
176 _GLOBAL_MAX_RL( rhsMaxGlobal, myThid )
177 #endif /* ALLOW_CONST_RHSMAX */
178 #endif
179 #ifdef ALLOW_AUTODIFF_TAMC
180 CADJ STORE rhsNorm = comlev1_cg2d, key = ikey, byte = isbyte
181 #endif /* ALLOW_AUTODIFF_TAMC */
182 IF ( rhsMaxGlobal .NE. 0. ) THEN
183 rhsNorm = 1. _d 0 / rhsMaxGlobal
184 ELSE
185 rhsNorm = 1. _d 0
186 ENDIF
187
188 #ifdef ALLOW_AUTODIFF_TAMC
189 CADJ STORE cg2d_b = comlev1_cg2d, key = ikey, byte = isbyte
190 CADJ STORE cg2d_x = comlev1_cg2d, key = ikey, byte = isbyte
191 #endif /* ALLOW_AUTODIFF_TAMC */
192 DO bj=myByLo(myThid),myByHi(myThid)
193 DO bi=myBxLo(myThid),myBxHi(myThid)
194 DO J=1,sNy
195 DO I=1,sNx
196 cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
197 cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
198 ENDDO
199 ENDDO
200 ENDDO
201 ENDDO
202 C- end Normalise RHS
203 ENDIF
204
205 C-- Update overlaps
206 c CALL EXCH_XY_RL( cg2d_b, myThid )
207 CALL EXCH_XY_RL( cg2d_x, myThid )
208 CcnhDebugStarts
209 C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
210 C suff = 'normalised'
211 C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
212 CcnhDebugEnds
213
214 C-- Initial residual calculation
215 err = 0. _d 0
216 err_sq = 0. _d 0
217 sumRHS = 0. _d 0
218 #ifdef ALLOW_AUTODIFF_TAMC
219 CADJ STORE cg2d_b = comlev1_cg2d, key = ikey, byte = isbyte
220 CADJ STORE cg2d_x = comlev1_cg2d, key = ikey, byte = isbyte
221 #endif /* ALLOW_AUTODIFF_TAMC */
222 DO bj=myByLo(myThid),myByHi(myThid)
223 DO bi=myBxLo(myThid),myBxHi(myThid)
224 DO J=1-1,sNy+1
225 DO I=1-1,sNx+1
226 cg2d_s(I,J,bi,bj) = 0.
227 ENDDO
228 ENDDO
229 DO J=1,sNy
230 DO I=1,sNx
231 cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
232 & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
233 & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
234 & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
235 & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
236 & +aC2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
237 & )
238 c & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
239 c & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
240 c & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
241 c & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
242 c & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
243 c & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTfreesurf*cg2dNorm
244 cML & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
245 c & )
246 err_sq = err_sq +
247 & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
248 sumRHS = sumRHS +
249 & cg2d_b(I,J,bi,bj)
250 ENDDO
251 ENDDO
252 ENDDO
253 ENDDO
254
255 c CALL EXCH_S3D_RL( cg2d_r, 1, myThid )
256 CALL EXCH_XY_RL ( cg2d_r, myThid )
257 _GLOBAL_SUM_RL( sumRHS, myThid )
258 _GLOBAL_SUM_RL( err_sq, myThid )
259 IF ( err_sq .NE. 0. ) THEN
260 err = SQRT(err_sq)
261 ELSE
262 err = 0.
263 ENDIF
264 actualIts = 0
265 actualResidual = err
266
267 IF ( debugLevel .GE. debLevZero ) THEN
268 _BEGIN_MASTER( myThid )
269 WRITE(standardmessageunit,'(A,1P2E22.14)')
270 & ' cg2d: Sum(rhs),rhsMax = ', sumRHS,rhsMaxGlobal
271 _END_MASTER( myThid )
272 ENDIF
273 C _BARRIER
274 c _BEGIN_MASTER( myThid )
275 c WRITE(*,'(A,I6,1PE30.14)') ' CG2D_NSA iters, err = ',
276 c & actualIts, actualResidual
277 c _END_MASTER( myThid )
278 firstResidual=actualResidual
279 cg2dTolerance_sq = cg2dTolerance**2
280
281 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
282 Cml begin main solver loop
283 #if ((defined ALLOW_AUTODIFF_TAMC) && (defined ALLOW_LOOP_DIRECTIVE))
284 CADJ LOOP = iteration, cg2d_x = comlev_cg2d_iter
285 #endif /* ALLOW_AUTODIFF_TAMC and ALLOW_LOOP_DIRECTIVE */
286 DO it2d=1, numIters
287 #ifdef ALLOW_LOOP_DIRECTIVE
288 CML it2d = 0
289 CML DO WHILE ( err_sq .GT. cg2dTolerance_sq .and. it2d .LT. numIters )
290 CML it2d = it2d+1
291 #endif /* ALLOW_LOOP_DIRECTIVE */
292
293 #ifdef ALLOW_AUTODIFF_TAMC
294 icg2dkey = (ikey-1)*numItersMax + it2d
295 CMLCADJ STORE err = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
296 CADJ STORE err_sq = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
297 #endif /* ALLOW_AUTODIFF_TAMC */
298 CML IF ( err .LT. cg2dTolerance ) THEN
299 IF ( err_sq .LT. cg2dTolerance_sq ) THEN
300 Cml DO NOTHING
301 ELSE
302
303 CcnhDebugStarts
304 C WRITE(*,*) ' CG2D_NSA: Iteration ',it2d-1,' residual = ',
305 C & actualResidual
306 CcnhDebugEnds
307 C-- Solve preconditioning equation and update
308 C-- conjugate direction vector "s".
309 eta_qrN = 0. _d 0
310 #ifdef ALLOW_AUTODIFF_TAMC
311 CADJ STORE cg2d_r = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
312 CADJ STORE cg2d_s = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
313 #endif /* ALLOW_AUTODIFF_TAMC */
314 DO bj=myByLo(myThid),myByHi(myThid)
315 DO bi=myBxLo(myThid),myBxHi(myThid)
316 DO J=1,sNy
317 DO I=1,sNx
318 cg2d_z(I,J,bi,bj) =
319 & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
320 & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
321 & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
322 & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
323 & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
324 CcnhDebugStarts
325 C cg2d_z(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
326 CcnhDebugEnds
327 eta_qrN = eta_qrN
328 & +cg2d_z(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
329 ENDDO
330 ENDDO
331 ENDDO
332 ENDDO
333
334 _GLOBAL_SUM_RL(eta_qrN, myThid)
335 CcnhDebugStarts
336 C WRITE(*,*) ' CG2D_NSA: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
337 CcnhDebugEnds
338 #ifdef ALLOW_AUTODIFF_TAMC
339 CMLCADJ STORE eta_qrNM1 = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
340 CADJ STORE recip_eta_qrNM1 = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
341 #endif /* ALLOW_AUTODIFF_TAMC */
342 CML cgBeta = eta_qrN/eta_qrNM1
343 cgBeta = eta_qrN*recip_eta_qrNM1
344 CcnhDebugStarts
345 C WRITE(*,*) ' CG2D_NSA: Iteration ',it2d-1,' beta = ',cgBeta
346 CcnhDebugEnds
347 Cml store normalisation factor for the next interation
348 Cml (in case there is one).
349 CML store the inverse of the normalization factor for higher precision
350 CML eta_qrNM1 = eta_qrN
351 recip_eta_qrNM1 = 1./eta_qrN
352
353 DO bj=myByLo(myThid),myByHi(myThid)
354 DO bi=myBxLo(myThid),myBxHi(myThid)
355 DO J=1,sNy
356 DO I=1,sNx
357 cg2d_s(I,J,bi,bj) = cg2d_z(I,J,bi,bj)
358 & + cgBeta*cg2d_s(I,J,bi,bj)
359 ENDDO
360 ENDDO
361 ENDDO
362 ENDDO
363
364 C-- Do exchanges that require messages i.e. between
365 C-- processes.
366 c CALL EXCH_S3D_RL( cg2d_s, 1, myThid )
367 CALL EXCH_XY_RL ( cg2d_s, myThid )
368
369 C== Evaluate laplace operator on conjugate gradient vector
370 C== q = A.s
371 alpha = 0. _d 0
372 alpha_aux = 0. _d 0
373 #ifdef ALLOW_AUTODIFF_TAMC
374 #ifndef ALLOW_LOOP_DIRECTIVE
375 CADJ STORE cg2d_s = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
376 #endif /* not ALLOW_LOOP_DIRECTIVE */
377 #endif /* ALLOW_AUTODIFF_TAMC */
378 DO bj=myByLo(myThid),myByHi(myThid)
379 DO bi=myBxLo(myThid),myBxHi(myThid)
380 DO J=1,sNy
381 DO I=1,sNx
382 cg2d_q(I,J,bi,bj) =
383 & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
384 & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
385 & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
386 & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
387 & +aC2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
388 c & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
389 c & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
390 c & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
391 c & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
392 c & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
393 c & cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTfreesurf*cg2dNorm
394 cML & cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
395 alpha_aux = alpha_aux+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
396 ENDDO
397 ENDDO
398 ENDDO
399 ENDDO
400 _GLOBAL_SUM_RL(alpha_aux,myThid)
401 CcnhDebugStarts
402 C WRITE(*,*) ' CG2D_NSA: Iteration ',it2d-1,' SUM(s*q)= ',alpha_aux
403 CcnhDebugEnds
404 alpha = eta_qrN/alpha_aux
405 CcnhDebugStarts
406 C WRITE(*,*) ' CG2D_NSA: Iteration ',it2d-1,' alpha= ',alpha
407 CcnhDebugEnds
408
409 C== Update solution and residual vectors
410 C Now compute "interior" points.
411 err = 0. _d 0
412 err_sq = 0. _d 0
413 #ifdef ALLOW_AUTODIFF_TAMC
414 #ifndef ALLOW_LOOP_DIRECTIVE
415 CADJ STORE cg2d_r = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
416 #endif /* ALLOW_LOOP_DIRECTIVE */
417 #endif /* ALLOW_AUTODIFF_TAMC */
418 DO bj=myByLo(myThid),myByHi(myThid)
419 DO bi=myBxLo(myThid),myBxHi(myThid)
420 DO J=1,sNy
421 DO I=1,sNx
422 cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
423 cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
424 err_sq = err_sq+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
425 ENDDO
426 ENDDO
427 ENDDO
428 ENDDO
429
430 _GLOBAL_SUM_RL( err_sq , myThid )
431 if ( err_sq .ne. 0. ) then
432 err = SQRT(err_sq)
433 else
434 err = 0.
435 end if
436 actualIts = it2d
437 actualResidual = err
438
439 c CALL EXCH_S3D_RL( cg2d_r, 1, myThid )
440 CALL EXCH_XY_RL ( cg2d_r, myThid )
441
442 Cml end of IF ( err .LT. cg2dTolerance ) THEN; ELSE
443 ENDIF
444 Cml end main solver loop
445 ENDDO
446
447 IF (cg2dNormaliseRHS) THEN
448 #ifdef ALLOW_AUTODIFF_TAMC
449 CADJ STORE rhsNorm = comlev1_cg2d, key = ikey, byte = isbyte
450 CADJ STORE cg2d_x = comlev1_cg2d, key = ikey, byte = isbyte
451 #endif /* ALLOW_AUTODIFF_TAMC */
452 C-- Un-normalise the answer
453 DO bj=myByLo(myThid),myByHi(myThid)
454 DO bi=myBxLo(myThid),myBxHi(myThid)
455 DO J=1,sNy
456 DO I=1,sNx
457 cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
458 ENDDO
459 ENDDO
460 ENDDO
461 ENDDO
462 ENDIF
463
464 C The following exchange was moved up to solve_for_pressure
465 C for compatibility with TAMC.
466 C _EXCH_XY_RL(cg2d_x, myThid )
467 c _BEGIN_MASTER( myThid )
468 c WRITE(*,'(A,I6,1PE30.14)') ' CG2D_NSA iters, err = ',
469 c & actualIts, actualResidual
470 c _END_MASTER( )
471
472 C-- Return parameters to caller
473 lastResidual=actualResidual
474 numIters=actualIts
475
476 #endif /* ALLOW_CG2D_NSA */
477 RETURN
478 END
479
480 # if ((defined ALLOW_AUTODIFF_TAMC) && (defined ALLOW_LOOP_DIRECTIVE))
481 C
482 C These routines are routinely part of the TAMC/TAF library that is
483 C not included in the MITcgm, therefore they are mimicked here.
484 C
485 subroutine adstore(chardum,int1,idow,int2,int3,icount)
486
487 implicit none
488
489 #include "SIZE.h"
490 #include "tamc.h"
491
492 character*(*) chardum
493 integer int1, int2, int3, idow, icount
494
495 C the length of this vector must be greater or equal
496 C twice the number of timesteps
497 integer nidow
498 #ifdef ALLOW_TAMC_CHECKPOINTING
499 parameter ( nidow = 2*nchklev_1*nchklev_2*nchklev_3 )
500 #else
501 parameter ( nidow = 1000000 )
502 #endif /* ALLOW_TAMC_CHECKPOINTING */
503 integer istoreidow(nidow)
504 common /istorecommon/ istoreidow
505
506 print *, 'adstore: ', chardum, int1, idow, int2, int3, icount
507
508 if ( icount .gt. nidow ) then
509 print *, 'adstore: error: icount > nidow = ', nidow
510 stop 'ABNORMAL STOP in adstore'
511 endif
512
513 istoreidow(icount) = idow
514
515 return
516 end
517
518 subroutine adresto(chardum,int1,idow,int2,int3,icount)
519
520 implicit none
521
522 #include "SIZE.h"
523 #include "tamc.h"
524
525 character*(*) chardum
526 integer int1, int2, int3, idow, icount
527
528
529 C the length of this vector must be greater or equal
530 C twice the number of timesteps
531 integer nidow
532 #ifdef ALLOW_TAMC_CHECKPOINTING
533 parameter ( nidow = 2*nchklev_1*nchklev_2*nchklev_3 )
534 #else
535 parameter ( nidow = 1000000 )
536 #endif /* ALLOW_TAMC_CHECKPOINTING */
537 integer istoreidow(nidow)
538 common /istorecommon/ istoreidow
539
540 print *, 'adresto: ', chardum, int1, idow, int2, int3, icount
541
542 if ( icount .gt. nidow ) then
543 print *, 'adstore: error: icount > nidow = ', nidow
544 stop 'ABNORMAL STOP in adstore'
545 endif
546
547 idow = istoreidow(icount)
548
549 return
550 end
551 # endif /* ALLOW_AUTODIFF_TAMC and ALLOW_LOOP_DIRECTIVE */
552

  ViewVC Help
Powered by ViewVC 1.1.22