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

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

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


Revision 1.1 - (hide annotations) (download)
Wed Jun 7 01:45:42 2006 UTC (17 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58l_post, mitgcm_mapl_00, checkpoint58u_post, checkpoint58w_post, checkpoint60, checkpoint61, checkpoint58r_post, checkpoint58n_post, checkpoint58x_post, checkpoint58t_post, checkpoint58h_post, checkpoint58q_post, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint58j_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j, checkpoint59, checkpoint58i_post, checkpoint58g_post, checkpoint58o_post, checkpoint58y_post, checkpoint58k_post, checkpoint58v_post, checkpoint58s_post, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint58p_post, checkpoint61a, checkpoint61l, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint58m_post
New routines for bottom topog. control.

1 heimbach 1.1 C $Header: $
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 "GRID.h"
58     #include "CG2D.h"
59     #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 myThid - Thread on which I am working.
68     C cg2d_b - The source term or "right hand side"
69     C cg2d_x - The solution
70     C firstResidual - the initial residual before any iterations
71     C lastResidual - the actual residual reached
72     C numIters - Entry: the maximum number of iterations allowed
73     C Exit: the actual number of iterations used
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 - Block index in X and Y.
87     C bj
88     C eta_qrN - Used in computing search directions
89     C eta_qrNM1 suffix N and NM1 denote current and
90     C cgBeta previous iterations respectively.
91     C recip_eta_qrNM1 reciprocal of eta_qrNM1
92     C alpha
93     C alpha_aux - to avoid the statement: alpha = 1./alpha (for TAMC/TAF)
94     C sumRHS - Sum of right-hand-side. Sometimes this is a
95     C useful debuggin/trouble shooting diagnostic.
96     C For neumann problems sumRHS needs to be ~0.
97     C or they converge at a non-zero residual.
98     C err - Measure of residual of Ax - b, usually the norm.
99     C err_sq - square of err (for TAMC/TAF)
100     C I, J, N - Loop counters ( N counts CG iterations )
101     INTEGER actualIts
102     _RL actualResidual
103     INTEGER bi, bj
104     INTEGER I, J, it2d
105    
106     _RL err
107     _RL err_sq
108     _RL eta_qrN
109     _RL eta_qrNM1
110     _RL recip_eta_qrNM1
111     _RL cgBeta
112     _RL alpha
113     _RL alpha_aux
114     _RL sumRHS
115     _RL rhsMax, rhsMaxGlobal
116     _RL rhsNorm
117     _RL cg2dTolerance_sq
118     CEOP
119    
120     #ifdef ALLOW_AUTODIFF_TAMC
121     IF ( numIters .GT. numItersMax ) THEN
122     WRITE(standardMessageUnit,'(A,I10)')
123     & 'CG2D_NSA: numIters > numItersMax = ', numItersMax
124     STOP 'NON-NORMAL in CG2D_NSA'
125     ENDIF
126     #endif /* ALLOW_AUTODIFF_TAMC */
127    
128     CcnhDebugStarts
129     C CHARACTER*(MAX_LEN_FNAM) suff
130     CcnhDebugEnds
131    
132     #ifdef ALLOW_AUTODIFF_TAMC
133     act1 = myThid - 1
134     max1 = nTx*nTy
135     act2 = ikey_dynamics - 1
136     ikey = (act1 + 1) + act2*max1
137     #endif /* ALLOW_AUTODIFF_TAMC */
138    
139     C-- Initialise inverter
140     eta_qrNM1 = 1. _d 0
141     recip_eta_qrNM1 = 1./eta_qrNM1
142    
143     CcnhDebugStarts
144     C _EXCH_XY_R8( cg2d_b, myThid )
145     C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid )
146     C suff = 'unnormalised'
147     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
148     C STOP
149     CcnhDebugEnds
150    
151     C-- Normalise RHS
152     #ifdef ALLOW_AUTODIFF_TAMC
153     CADJ STORE cg2d_b = comlev1_cg2d, key = ikey, byte = isbyte
154     #endif /* ALLOW_AUTODIFF_TAMC */
155     rhsMax = 0. _d 0
156     DO bj=myByLo(myThid),myByHi(myThid)
157     DO bi=myBxLo(myThid),myBxHi(myThid)
158     DO J=1,sNy
159     DO I=1,sNx
160     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm
161     rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax)
162     ENDDO
163     ENDDO
164     ENDDO
165     ENDDO
166    
167     IF (cg2dNormaliseRHS) THEN
168     C - Normalise RHS :
169     #ifdef LETS_MAKE_JAM
170     C _GLOBAL_MAX_R8( rhsMax, myThid )
171     rhsMaxGlobal=1.
172     #else
173     #ifdef ALLOW_CONST_RHSMAX
174     rhsMaxGlobal=1.
175     #else
176     rhsMaxGlobal=rhsMax
177     _GLOBAL_MAX_R8( rhsMaxGlobal, myThid )
178     #endif /* ALLOW_CONST_RHSMAX */
179     #endif
180     #ifdef ALLOW_AUTODIFF_TAMC
181     CADJ STORE rhsNorm = comlev1_cg2d, key = ikey, byte = isbyte
182     #endif /* ALLOW_AUTODIFF_TAMC */
183     IF ( rhsMaxGlobal .NE. 0. ) THEN
184     rhsNorm = 1. _d 0 / rhsMaxGlobal
185     ELSE
186     rhsNorm = 1. _d 0
187     ENDIF
188    
189     #ifdef ALLOW_AUTODIFF_TAMC
190     CADJ STORE cg2d_b = comlev1_cg2d, key = ikey, byte = isbyte
191     CADJ STORE cg2d_x = comlev1_cg2d, key = ikey, byte = isbyte
192     #endif /* ALLOW_AUTODIFF_TAMC */
193     DO bj=myByLo(myThid),myByHi(myThid)
194     DO bi=myBxLo(myThid),myBxHi(myThid)
195     DO J=1,sNy
196     DO I=1,sNx
197     cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
198     cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
199     ENDDO
200     ENDDO
201     ENDDO
202     ENDDO
203     C- end Normalise RHS
204     ENDIF
205    
206     C-- Update overlaps
207     _EXCH_XY_R8( cg2d_b, myThid )
208     _EXCH_XY_R8( cg2d_x, myThid )
209     CcnhDebugStarts
210     C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
211     C suff = 'normalised'
212     C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid)
213     CcnhDebugEnds
214    
215     C-- Initial residual calculation
216     err = 0. _d 0
217     err_sq = 0. _d 0
218     sumRHS = 0. _d 0
219     #ifdef ALLOW_AUTODIFF_TAMC
220     CADJ STORE cg2d_b = comlev1_cg2d, key = ikey, byte = isbyte
221     CADJ STORE cg2d_x = comlev1_cg2d, key = ikey, byte = isbyte
222     #endif /* ALLOW_AUTODIFF_TAMC */
223     DO bj=myByLo(myThid),myByHi(myThid)
224     DO bi=myBxLo(myThid),myBxHi(myThid)
225     DO J=1,sNy
226     DO I=1,sNx
227     cg2d_s(I,J,bi,bj) = 0.
228     cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
229     & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
230     & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
231     & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
232     & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
233     & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
234     & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
235     & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
236     & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)
237     & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
238     & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTfreesurf*cg2dNorm
239     CML & cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm
240     & )
241     err_sq = err_sq +
242     & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
243     sumRHS = sumRHS +
244     & cg2d_b(I,J,bi,bj)
245     ENDDO
246     ENDDO
247     ENDDO
248     ENDDO
249    
250     #ifdef LETS_MAKE_JAM
251     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
252     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
253     #else
254     _EXCH_XY_R8( cg2d_r, myThid )
255     _EXCH_XY_R8( cg2d_s, myThid )
256     #endif
257     _GLOBAL_SUM_R8( sumRHS, myThid )
258     _GLOBAL_SUM_R8( err_sq, myThid )
259     if ( err_sq .ne. 0. ) then
260     err = SQRT(err_sq)
261     else
262     err = 0.
263     end if
264     actualIts = 0
265     actualResidual = err
266    
267     _BEGIN_MASTER( myThid )
268     write(standardMessageUnit,'(A,1P2E22.14)')
269     & ' cg2d: Sum(rhs),rhsMax = ', sumRHS,rhsMaxGlobal
270     _END_MASTER( )
271     C _BARRIER
272     c _BEGIN_MASTER( myThid )
273     c WRITE(*,'(A,I6,1PE30.14)') ' CG2D_NSA iters, err = ',
274     c & actualIts, actualResidual
275     c _END_MASTER( )
276     firstResidual=actualResidual
277     cg2dTolerance_sq = cg2dTolerance**2
278    
279     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
280     Cml begin main solver loop
281     #if ((defined ALLOW_AUTODIFF_TAMC) && (defined ALLOW_LOOP_DIRECTIVE))
282     CADJ LOOP = iteration, cg2d_x = comlev_cg2d_iter
283     #endif /* ALLOW_AUTODIFF_TAMC and ALLOW_LOOP_DIRECTIVE */
284     DO it2d=1, numIters
285     #ifdef ALLOW_LOOP_DIRECTIVE
286     CML it2d = 0
287     CML DO WHILE ( err_sq .GT. cg2dTolerance_sq .and. it2d .LT. numIters )
288     CML it2d = it2d+1
289     #endif /* ALLOW_LOOP_DIRECTIVE */
290    
291     #ifdef ALLOW_AUTODIFF_TAMC
292     icg2dkey = (ikey-1)*numItersMax + it2d
293     CMLCADJ STORE err = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
294     CADJ STORE err_sq = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
295     #endif /* ALLOW_AUTODIFF_TAMC */
296     CML IF ( err .LT. cg2dTolerance ) THEN
297     IF ( err_sq .LT. cg2dTolerance_sq ) THEN
298     Cml DO NOTHING
299     ELSE
300    
301     CcnhDebugStarts
302     C WRITE(*,*) ' CG2D_NSA: Iteration ',it2d-1,' residual = ',
303     C & actualResidual
304     CcnhDebugEnds
305     C-- Solve preconditioning equation and update
306     C-- conjugate direction vector "s".
307     eta_qrN = 0. _d 0
308     #ifdef ALLOW_AUTODIFF_TAMC
309     CADJ STORE cg2d_r = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
310     CADJ STORE cg2d_s = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
311     #endif /* ALLOW_AUTODIFF_TAMC */
312     DO bj=myByLo(myThid),myByHi(myThid)
313     DO bi=myBxLo(myThid),myBxHi(myThid)
314     DO J=1,sNy
315     DO I=1,sNx
316     cg2d_z(I,J,bi,bj) =
317     & pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
318     & +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
319     & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
320     & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
321     & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
322     CcnhDebugStarts
323     C cg2d_z(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj)
324     CcnhDebugEnds
325     eta_qrN = eta_qrN
326     & +cg2d_z(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
327     ENDDO
328     ENDDO
329     ENDDO
330     ENDDO
331    
332     _GLOBAL_SUM_R8(eta_qrN, myThid)
333     CcnhDebugStarts
334     C WRITE(*,*) ' CG2D_NSA: Iteration ',it2d-1,' eta_qrN = ',eta_qrN
335     CcnhDebugEnds
336     #ifdef ALLOW_AUTODIFF_TAMC
337     CMLCADJ STORE eta_qrNM1 = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
338     CADJ STORE recip_eta_qrNM1 = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
339     #endif /* ALLOW_AUTODIFF_TAMC */
340     CML cgBeta = eta_qrN/eta_qrNM1
341     cgBeta = eta_qrN*recip_eta_qrNM1
342     CcnhDebugStarts
343     C WRITE(*,*) ' CG2D_NSA: Iteration ',it2d-1,' beta = ',cgBeta
344     CcnhDebugEnds
345     Cml store normalisation factor for the next interation
346     Cml (in case there is one).
347     CML store the inverse of the normalization factor for higher precision
348     CML eta_qrNM1 = eta_qrN
349     recip_eta_qrNM1 = 1./eta_qrN
350    
351     DO bj=myByLo(myThid),myByHi(myThid)
352     DO bi=myBxLo(myThid),myBxHi(myThid)
353     DO J=1,sNy
354     DO I=1,sNx
355     cg2d_s(I,J,bi,bj) = cg2d_z(I,J,bi,bj)
356     & + cgBeta*cg2d_s(I,J,bi,bj)
357     ENDDO
358     ENDDO
359     ENDDO
360     ENDDO
361    
362     C-- Do exchanges that require messages i.e. between
363     C-- processes.
364     #ifdef LETS_MAKE_JAM
365     CALL EXCH_XY_O1_R8_JAM( cg2d_s )
366     #else
367     _EXCH_XY_R8( cg2d_s, myThid )
368     #endif
369    
370     C== Evaluate laplace operator on conjugate gradient vector
371     C== q = A.s
372     alpha = 0. _d 0
373     alpha_aux = 0. _d 0
374     #ifdef ALLOW_AUTODIFF_TAMC
375     #ifndef ALLOW_LOOP_DIRECTIVE
376     CADJ STORE cg2d_s = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
377     #endif /* not ALLOW_LOOP_DIRECTIVE */
378     #endif /* ALLOW_AUTODIFF_TAMC */
379     DO bj=myByLo(myThid),myByHi(myThid)
380     DO bi=myBxLo(myThid),myBxHi(myThid)
381     DO J=1,sNy
382     DO I=1,sNx
383     cg2d_q(I,J,bi,bj) =
384     & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
385     & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
386     & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
387     & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
388     & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
389     & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
390     & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
391     & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
392     & -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)*
393     & 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_R8(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_R8( 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     #ifdef LETS_MAKE_JAM
440     CALL EXCH_XY_O1_R8_JAM( cg2d_r )
441     #else
442     _EXCH_XY_R8( cg2d_r, myThid )
443     _EXCH_XY_R8( cg2d_x, myThid )
444     #endif
445    
446     Cml end of IF ( err .LT. cg2dTolerance ) THEN; ELSE
447     ENDIF
448     Cml end main solver loop
449     ENDDO
450    
451     IF (cg2dNormaliseRHS) THEN
452     #ifdef ALLOW_AUTODIFF_TAMC
453     CADJ STORE rhsNorm = comlev1_cg2d, key = ikey, byte = isbyte
454     CADJ STORE cg2d_x = comlev1_cg2d, key = ikey, byte = isbyte
455     #endif /* ALLOW_AUTODIFF_TAMC */
456     C-- Un-normalise the answer
457     DO bj=myByLo(myThid),myByHi(myThid)
458     DO bi=myBxLo(myThid),myBxHi(myThid)
459     DO J=1,sNy
460     DO I=1,sNx
461     cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
462     ENDDO
463     ENDDO
464     ENDDO
465     ENDDO
466     ENDIF
467    
468     C The following exchange was moved up to solve_for_pressure
469     C for compatibility with TAMC.
470     C _EXCH_XY_R8(cg2d_x, myThid )
471     c _BEGIN_MASTER( myThid )
472     c WRITE(*,'(A,I6,1PE30.14)') ' CG2D_NSA iters, err = ',
473     c & actualIts, actualResidual
474     c _END_MASTER( )
475    
476     C-- Return parameters to caller
477     lastResidual=actualResidual
478     numIters=actualIts
479    
480     #endif /* ALLOW_CG2D_NSA */
481     RETURN
482     END
483    
484     # if ((defined ALLOW_AUTODIFF_TAMC) && (defined ALLOW_LOOP_DIRECTIVE))
485     C
486     C These routines are routinely part of the TAMC/TAF library that is
487     C not included in the MITcgm, therefore they are mimicked here.
488     C
489     subroutine adstore(chardum,int1,idow,int2,int3,icount)
490    
491     implicit none
492    
493     #include "SIZE.h"
494     #include "tamc.h"
495    
496     character*(*) chardum
497     integer int1, int2, int3, idow, icount
498    
499     C the length of this vector must be greater or equal
500     C twice the number of timesteps
501     integer nidow
502     #ifdef ALLOW_TAMC_CHECKPOINTING
503     parameter ( nidow = 2*nchklev_1*nchklev_2*nchklev_3 )
504     #else
505     parameter ( nidow = 1000000 )
506     #endif /* ALLOW_TAMC_CHECKPOINTING */
507     integer istoreidow(nidow)
508     common /istorecommon/ istoreidow
509    
510     print *, 'adstore: ', chardum, int1, idow, int2, int3, icount
511    
512     if ( icount .gt. nidow ) then
513     print *, 'adstore: error: icount > nidow = ', nidow
514     stop 'ABNORMAL STOP in adstore'
515     endif
516    
517     istoreidow(icount) = idow
518    
519     return
520     end
521    
522     subroutine adresto(chardum,int1,idow,int2,int3,icount)
523    
524     implicit none
525    
526     #include "SIZE.h"
527     #include "tamc.h"
528    
529     character*(*) chardum
530     integer int1, int2, int3, idow, icount
531    
532    
533     C the length of this vector must be greater or equal
534     C twice the number of timesteps
535     integer nidow
536     #ifdef ALLOW_TAMC_CHECKPOINTING
537     parameter ( nidow = 2*nchklev_1*nchklev_2*nchklev_3 )
538     #else
539     parameter ( nidow = 1000000 )
540     #endif /* ALLOW_TAMC_CHECKPOINTING */
541     integer istoreidow(nidow)
542     common /istorecommon/ istoreidow
543    
544     print *, 'adresto: ', chardum, int1, idow, int2, int3, icount
545    
546     if ( icount .gt. nidow ) then
547     print *, 'adstore: error: icount > nidow = ', nidow
548     stop 'ABNORMAL STOP in adstore'
549     endif
550    
551     idow = istoreidow(icount)
552    
553     return
554     end
555     # endif /* ALLOW_AUTODIFF_TAMC and ALLOW_LOOP_DIRECTIVE */
556    
557    

  ViewVC Help
Powered by ViewVC 1.1.22