/[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.98 - (hide annotations) (download)
Tue Jul 8 15:00:26 2003 UTC (20 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51f_post, checkpoint51d_post, branchpoint-genmake2, checkpoint51c_post, checkpoint51e_post, checkpoint51f_pre
Branch point for: branch-genmake2
Changes since 1.97: +2 -2 lines
o introducing integer flag debugLevel
o introducing pathname variable mdsioLocalDir for mdsio

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

  ViewVC Help
Powered by ViewVC 1.1.22