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

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

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


Revision 1.91 - (show annotations) (download)
Fri Nov 15 03:01:21 2002 UTC (21 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint47e_post, checkpoint47c_post, checkpoint48b_post, checkpoint48c_pre, checkpoint47d_pre, checkpoint47a_post, checkpoint47i_post, checkpoint47d_post, checkpoint47g_post, checkpoint48a_post, checkpoint47j_post, branch-exfmods-tag, checkpoint48c_post, checkpoint47b_post, checkpoint47f_post, checkpoint47, checkpoint48, checkpoint47h_post
Branch point for: branch-exfmods-curt
Changes since 1.90: +15 -12 lines
differentiable version of checkpoint46n_post
o external_fields_load now part of differentiation list
o pressure needs multiple storing;
  would be nice to have store_pressure at beginning or
  end of forward_step, e.g. by having phiHyd global (5-dim.)
  (NB: pressure is needed for certain cases in find_rho,
  which is also invoked through convective_adjustment).
o recomputations in find_rho for cases
 'JMD95'/'UNESCO' or 'MDJWF' are OK.
o #define ATMOSPHERIC_LOADING should be differentiable
o ini_forcing shifted to begining of initialise_varia

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.90 2002/09/18 16:38:01 mlosch Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: DYNAMICS
8 C !INTERFACE:
9 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
10 C !DESCRIPTION: \bv
11 C *==========================================================*
12 C | SUBROUTINE DYNAMICS
13 C | o Controlling routine for the explicit part of the model
14 C | dynamics.
15 C *==========================================================*
16 C | This routine evaluates the "dynamics" terms for each
17 C | block of ocean in turn. Because the blocks of ocean have
18 C | overlap regions they are independent of one another.
19 C | If terms involving lateral integrals are needed in this
20 C | routine care will be needed. Similarly finite-difference
21 C | operations with stencils wider than the overlap region
22 C | require special consideration.
23 C | The algorithm...
24 C |
25 C | "Correction Step"
26 C | =================
27 C | Here we update the horizontal velocities with the surface
28 C | pressure such that the resulting flow is either consistent
29 C | with the free-surface evolution or the rigid-lid:
30 C | U[n] = U* + dt x d/dx P
31 C | V[n] = V* + dt x d/dy P
32 C |
33 C | "Calculation of Gs"
34 C | ===================
35 C | This is where all the accelerations and tendencies (ie.
36 C | physics, parameterizations etc...) are calculated
37 C | rho = rho ( theta[n], salt[n] )
38 C | b = b(rho, theta)
39 C | K31 = K31 ( rho )
40 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
41 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
42 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
43 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
44 C |
45 C | "Time-stepping" or "Prediction"
46 C | ================================
47 C | The models variables are stepped forward with the appropriate
48 C | time-stepping scheme (currently we use Adams-Bashforth II)
49 C | - For momentum, the result is always *only* a "prediction"
50 C | in that the flow may be divergent and will be "corrected"
51 C | later with a surface pressure gradient.
52 C | - Normally for tracers the result is the new field at time
53 C | level [n+1} *BUT* in the case of implicit diffusion the result
54 C | is also *only* a prediction.
55 C | - We denote "predictors" with an asterisk (*).
56 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
57 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
58 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
59 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
60 C | With implicit diffusion:
61 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
62 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
63 C | (1 + dt * K * d_zz) theta[n] = theta*
64 C | (1 + dt * K * d_zz) salt[n] = salt*
65 C |
66 C *==========================================================*
67 C \ev
68 C !USES:
69 IMPLICIT NONE
70 C == Global variables ===
71 #include "SIZE.h"
72 #include "EEPARAMS.h"
73 #include "PARAMS.h"
74 #include "DYNVARS.h"
75 #include "GRID.h"
76 #ifdef ALLOW_PASSIVE_TRACER
77 #include "TR1.h"
78 #endif
79 #ifdef ALLOW_AUTODIFF_TAMC
80 # include "tamc.h"
81 # include "tamc_keys.h"
82 # include "FFIELDS.h"
83 # include "EOS.h"
84 # ifdef ALLOW_KPP
85 # include "KPP.h"
86 # endif
87 #endif /* ALLOW_AUTODIFF_TAMC */
88 #ifdef ALLOW_TIMEAVE
89 #include "TIMEAVE_STATV.h"
90 #endif
91
92 C !CALLING SEQUENCE:
93 C DYNAMICS()
94 C |
95 C |-- CALC_GRAD_PHI_SURF
96 C |
97 C |-- CALC_VISCOSITY
98 C |
99 C |-- CALC_PHI_HYD
100 C |
101 C |-- STORE_PRESSURE
102 C |
103 C |-- MOM_FLUXFORM
104 C |
105 C |-- MOM_VECINV
106 C |
107 C |-- TIMESTEP
108 C |
109 C |-- OBCS_APPLY_UV
110 C |
111 C |-- IMPLDIFF
112 C |
113 C |-- OBCS_APPLY_UV
114 C |
115 C |-- CALL TIMEAVE_CUMUL_1T
116 C |-- CALL DEBUG_STATS_RL
117
118 C !INPUT/OUTPUT PARAMETERS:
119 C == Routine arguments ==
120 C myTime - Current time in simulation
121 C myIter - Current iteration number in simulation
122 C myThid - Thread number for this instance of the routine.
123 _RL myTime
124 INTEGER myIter
125 INTEGER myThid
126
127 C !LOCAL VARIABLES:
128 C == Local variables
129 C fVer[STUV] o fVer: Vertical flux term - note fVer
130 C is "pipelined" in the vertical
131 C so we need an fVer for each
132 C variable.
133 C rhoK, rhoKM1 - Density at current level, and level above
134 C phiHyd - Hydrostatic part of the potential phiHydi.
135 C In z coords phiHydiHyd is the hydrostatic
136 C Potential (=pressure/rho0) anomaly
137 C In p coords phiHydiHyd is the geopotential
138 C surface height anomaly.
139 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
140 C phiSurfY or geopotentiel (atmos) in X and Y direction
141 C iMin, iMax - Ranges and sub-block indices on which calculations
142 C jMin, jMax are applied.
143 C bi, bj
144 C k, kup, - Index for layer above and below. kup and kDown
145 C kDown, km1 are switched with layer to be the appropriate
146 C index into fVerTerm.
147 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
148 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
149 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
150 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
155 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
156
157 INTEGER iMin, iMax
158 INTEGER jMin, jMax
159 INTEGER bi, bj
160 INTEGER i, j
161 INTEGER k, km1, kp1, kup, kDown
162
163 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
164 c CHARACTER*(MAX_LEN_MBUF) suff
165 c LOGICAL DIFFERENT_MULTIPLE
166 c EXTERNAL DIFFERENT_MULTIPLE
167 Cjmc(end)
168
169 C--- The algorithm...
170 C
171 C "Correction Step"
172 C =================
173 C Here we update the horizontal velocities with the surface
174 C pressure such that the resulting flow is either consistent
175 C with the free-surface evolution or the rigid-lid:
176 C U[n] = U* + dt x d/dx P
177 C V[n] = V* + dt x d/dy P
178 C
179 C "Calculation of Gs"
180 C ===================
181 C This is where all the accelerations and tendencies (ie.
182 C physics, parameterizations etc...) are calculated
183 C rho = rho ( theta[n], salt[n] )
184 C b = b(rho, theta)
185 C K31 = K31 ( rho )
186 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
187 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
188 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
189 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
190 C
191 C "Time-stepping" or "Prediction"
192 C ================================
193 C The models variables are stepped forward with the appropriate
194 C time-stepping scheme (currently we use Adams-Bashforth II)
195 C - For momentum, the result is always *only* a "prediction"
196 C in that the flow may be divergent and will be "corrected"
197 C later with a surface pressure gradient.
198 C - Normally for tracers the result is the new field at time
199 C level [n+1} *BUT* in the case of implicit diffusion the result
200 C is also *only* a prediction.
201 C - We denote "predictors" with an asterisk (*).
202 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
203 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
204 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
205 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
206 C With implicit diffusion:
207 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
208 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
209 C (1 + dt * K * d_zz) theta[n] = theta*
210 C (1 + dt * K * d_zz) salt[n] = salt*
211 C---
212 CEOP
213
214 C-- Set up work arrays with valid (i.e. not NaN) values
215 C These inital values do not alter the numerical results. They
216 C just ensure that all memory references are to valid floating
217 C point numbers. This prevents spurious hardware signals due to
218 C uninitialised but inert locations.
219 DO j=1-OLy,sNy+OLy
220 DO i=1-OLx,sNx+OLx
221 rhoKM1 (i,j) = 0. _d 0
222 rhok (i,j) = 0. _d 0
223 phiSurfX(i,j) = 0. _d 0
224 phiSurfY(i,j) = 0. _d 0
225 ENDDO
226 ENDDO
227
228 C-- Call to routine for calculation of
229 C Eliassen-Palm-flux-forced U-tendency,
230 C if desired:
231 #ifdef INCLUDE_EP_FORCING_CODE
232 CALL CALC_EP_FORCING(myThid)
233 #endif
234
235 #ifdef ALLOW_AUTODIFF_TAMC
236 C-- HPF directive to help TAMC
237 CHPF$ INDEPENDENT
238 #endif /* ALLOW_AUTODIFF_TAMC */
239
240 DO bj=myByLo(myThid),myByHi(myThid)
241
242 #ifdef ALLOW_AUTODIFF_TAMC
243 C-- HPF directive to help TAMC
244 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
245 CHPF$& ,phiHyd
246 CHPF$& ,KappaRU,KappaRV
247 CHPF$& )
248 #endif /* ALLOW_AUTODIFF_TAMC */
249
250 DO bi=myBxLo(myThid),myBxHi(myThid)
251
252 #ifdef ALLOW_AUTODIFF_TAMC
253 act1 = bi - myBxLo(myThid)
254 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
255 act2 = bj - myByLo(myThid)
256 max2 = myByHi(myThid) - myByLo(myThid) + 1
257 act3 = myThid - 1
258 max3 = nTx*nTy
259 act4 = ikey_dynamics - 1
260 idynkey = (act1 + 1) + act2*max1
261 & + act3*max1*max2
262 & + act4*max1*max2*max3
263 #endif /* ALLOW_AUTODIFF_TAMC */
264
265 C-- Set up work arrays that need valid initial values
266 DO j=1-OLy,sNy+OLy
267 DO i=1-OLx,sNx+OLx
268 DO k=1,Nr
269 phiHyd(i,j,k) = 0. _d 0
270 KappaRU(i,j,k) = 0. _d 0
271 KappaRV(i,j,k) = 0. _d 0
272 ENDDO
273 fVerU (i,j,1) = 0. _d 0
274 fVerU (i,j,2) = 0. _d 0
275 fVerV (i,j,1) = 0. _d 0
276 fVerV (i,j,2) = 0. _d 0
277 ENDDO
278 ENDDO
279
280 C-- Start computation of dynamics
281 iMin = 1-OLx+2
282 iMax = sNx+OLx-1
283 jMin = 1-OLy+2
284 jMax = sNy+OLy-1
285
286 #ifdef ALLOW_AUTODIFF_TAMC
287 CADJ STORE wvel (:,:,:,bi,bj) =
288 CADJ & comlev1_bibj, key = idynkey, byte = isbyte
289 #endif /* ALLOW_AUTODIFF_TAMC */
290
291 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
292 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
293 IF (implicSurfPress.NE.1.) THEN
294 CALL CALC_GRAD_PHI_SURF(
295 I bi,bj,iMin,iMax,jMin,jMax,
296 I etaN,
297 O phiSurfX,phiSurfY,
298 I myThid )
299 ENDIF
300
301 #ifdef ALLOW_AUTODIFF_TAMC
302 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
303 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
304 #ifdef ALLOW_KPP
305 CADJ STORE KPPviscAz (:,:,:,bi,bj)
306 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
307 #endif /* ALLOW_KPP */
308 #endif /* ALLOW_AUTODIFF_TAMC */
309
310 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
311 C-- Calculate the total vertical diffusivity
312 DO k=1,Nr
313 CALL CALC_VISCOSITY(
314 I bi,bj,iMin,iMax,jMin,jMax,k,
315 O KappaRU,KappaRV,
316 I myThid)
317 ENDDO
318 #endif
319
320 C-- Start of dynamics loop
321 DO k=1,Nr
322
323 C-- km1 Points to level above k (=k-1)
324 C-- kup Cycles through 1,2 to point to layer above
325 C-- kDown Cycles through 2,1 to point to current layer
326
327 km1 = MAX(1,k-1)
328 kp1 = MIN(k+1,Nr)
329 kup = 1+MOD(k+1,2)
330 kDown= 1+MOD(k,2)
331
332 #ifdef ALLOW_AUTODIFF_TAMC
333 kkey = (idynkey-1)*Nr + k
334 CADJ STORE pressure(:,:,k,bi,bj) = comlev1_bibj_k ,
335 CADJ & key=kkey , byte=isbyte
336 #endif /* ALLOW_AUTODIFF_TAMC */
337
338 C-- Integrate hydrostatic balance for phiHyd with BC of
339 C phiHyd(z=0)=0
340 C distinguishe between Stagger and Non Stagger time stepping
341 IF (staggerTimeStep) THEN
342 CALL CALC_PHI_HYD(
343 I bi,bj,iMin,iMax,jMin,jMax,k,
344 I gT, gS,
345 U phiHyd,
346 I myThid )
347 ELSE
348 CALL CALC_PHI_HYD(
349 I bi,bj,iMin,iMax,jMin,jMax,k,
350 I theta, salt,
351 U phiHyd,
352 I myThid )
353 ENDIF
354
355 C calculate pressure from phiHyd and store it on common block
356 C variable pressure
357 CALL STORE_PRESSURE( bi, bj, k, phiHyd, myThid )
358
359 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
360 C and step forward storing the result in gUnm1, gVnm1, etc...
361 IF ( momStepping ) THEN
362 #ifndef DISABLE_MOM_FLUXFORM
363 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
364 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
365 I phiHyd,KappaRU,KappaRV,
366 U fVerU, fVerV,
367 I myTime, myIter, myThid)
368 #endif
369 #ifndef DISABLE_MOM_VECINV
370 IF (vectorInvariantMomentum) CALL MOM_VECINV(
371 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
372 I phiHyd,KappaRU,KappaRV,
373 U fVerU, fVerV,
374 I myTime, myIter, myThid)
375 #endif
376 CALL TIMESTEP(
377 I bi,bj,iMin,iMax,jMin,jMax,k,
378 I phiHyd, phiSurfX, phiSurfY,
379 I myIter, myThid)
380
381 #ifdef ALLOW_OBCS
382 C-- Apply open boundary conditions
383 IF (useOBCS) THEN
384 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
385 END IF
386 #endif /* ALLOW_OBCS */
387
388 #ifdef ALLOW_AUTODIFF_TAMC
389 #ifdef INCLUDE_CD_CODE
390 ELSE
391 DO j=1-OLy,sNy+OLy
392 DO i=1-OLx,sNx+OLx
393 guCD(i,j,k,bi,bj) = 0.0
394 gvCD(i,j,k,bi,bj) = 0.0
395 END DO
396 END DO
397 #endif /* INCLUDE_CD_CODE */
398 #endif /* ALLOW_AUTODIFF_TAMC */
399 ENDIF
400
401
402 C-- end of dynamics k loop (1:Nr)
403 ENDDO
404
405 C-- Implicit viscosity
406 IF (implicitViscosity.AND.momStepping) THEN
407 #ifdef ALLOW_AUTODIFF_TAMC
408 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
409 #endif /* ALLOW_AUTODIFF_TAMC */
410 CALL IMPLDIFF(
411 I bi, bj, iMin, iMax, jMin, jMax,
412 I deltaTmom, KappaRU,recip_HFacW,
413 U gUNm1,
414 I myThid )
415 #ifdef ALLOW_AUTODIFF_TAMC
416 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
417 #endif /* ALLOW_AUTODIFF_TAMC */
418 CALL IMPLDIFF(
419 I bi, bj, iMin, iMax, jMin, jMax,
420 I deltaTmom, KappaRV,recip_HFacS,
421 U gVNm1,
422 I myThid )
423
424 #ifdef ALLOW_OBCS
425 C-- Apply open boundary conditions
426 IF (useOBCS) THEN
427 DO K=1,Nr
428 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
429 ENDDO
430 END IF
431 #endif /* ALLOW_OBCS */
432
433 #ifdef INCLUDE_CD_CODE
434 #ifdef ALLOW_AUTODIFF_TAMC
435 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
436 #endif /* ALLOW_AUTODIFF_TAMC */
437 CALL IMPLDIFF(
438 I bi, bj, iMin, iMax, jMin, jMax,
439 I deltaTmom, KappaRU,recip_HFacW,
440 U vVelD,
441 I myThid )
442 #ifdef ALLOW_AUTODIFF_TAMC
443 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
444 #endif /* ALLOW_AUTODIFF_TAMC */
445 CALL IMPLDIFF(
446 I bi, bj, iMin, iMax, jMin, jMax,
447 I deltaTmom, KappaRV,recip_HFacS,
448 U uVelD,
449 I myThid )
450 #endif /* INCLUDE_CD_CODE */
451 C-- End If implicitViscosity.AND.momStepping
452 ENDIF
453
454 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
455 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
456 c & .AND. buoyancyRelation .ne. 'OCEANIC' ) THEN
457 c WRITE(suff,'(I10.10)') myIter+1
458 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
459 c ENDIF
460 Cjmc(end)
461
462 #ifdef ALLOW_TIMEAVE
463 IF (taveFreq.GT.0.) THEN
464 CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
465 I deltaTclock, bi, bj, myThid)
466 ENDIF
467 #endif /* ALLOW_TIMEAVE */
468
469 ENDDO
470 ENDDO
471
472 Cml(
473 C In order to compare the variance of phiHydLow of a p/z-coordinate
474 C run with etaH of a z/p-coordinate run the drift of phiHydLow
475 C has to be removed by something like the following subroutine:
476 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
477 C & 'phiHydLow', myThid )
478 Cml)
479
480 #ifndef DISABLE_DEBUGMODE
481 If (debugMode) THEN
482 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
483 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
484 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
485 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
486 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
487 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
488 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
489 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
490 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
491 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
492 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
493 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
494 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
495 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
496 ENDIF
497 #endif
498
499 RETURN
500 END

  ViewVC Help
Powered by ViewVC 1.1.22