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

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

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


Revision 1.86 - (hide annotations) (download)
Sun Mar 24 02:36:39 2002 UTC (22 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint45a_post, checkpoint44h_post, checkpoint45
Changes since 1.85: +1 -4 lines
o Modified initialisations to break adjoint dependencies
o removed some store directives
o added options files for KPP, GMREDI

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

  ViewVC Help
Powered by ViewVC 1.1.22