/[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.95 - (hide annotations) (download)
Fri Feb 28 02:20:52 2003 UTC (21 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint48i_post, checkpoint50, checkpoint50b_pre, checkpoint50a_post, checkpoint49
Changes since 1.94: +3 -3 lines
Changes to restore differentiability of code w.r.t. previous tag
(mostly adding new routines to make list and replacing
pressure by totPhiHyd).

1 heimbach 1.95 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.94 2003/02/18 15:25:09 jmc 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.76 C-- Set up work arrays with valid (i.e. not NaN) values
208     C These inital values do not alter the numerical results. They
209     C just ensure that all memory references are to valid floating
210     C point numbers. This prevents spurious hardware signals due to
211     C uninitialised but inert locations.
212     DO j=1-OLy,sNy+OLy
213     DO i=1-OLx,sNx+OLx
214     phiSurfX(i,j) = 0. _d 0
215     phiSurfY(i,j) = 0. _d 0
216     ENDDO
217     ENDDO
218    
219 heimbach 1.88 C-- Call to routine for calculation of
220     C Eliassen-Palm-flux-forced U-tendency,
221     C if desired:
222     #ifdef INCLUDE_EP_FORCING_CODE
223     CALL CALC_EP_FORCING(myThid)
224     #endif
225    
226 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
227     C-- HPF directive to help TAMC
228     CHPF$ INDEPENDENT
229     #endif /* ALLOW_AUTODIFF_TAMC */
230    
231 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
232 heimbach 1.76
233     #ifdef ALLOW_AUTODIFF_TAMC
234     C-- HPF directive to help TAMC
235     CHPF$ INDEPENDENT, NEW (fVerU,fVerV
236 jmc 1.94 CHPF$& ,phiHydF
237 heimbach 1.76 CHPF$& ,KappaRU,KappaRV
238     CHPF$& )
239     #endif /* ALLOW_AUTODIFF_TAMC */
240    
241 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
242 heimbach 1.76
243     #ifdef ALLOW_AUTODIFF_TAMC
244     act1 = bi - myBxLo(myThid)
245     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
246     act2 = bj - myByLo(myThid)
247     max2 = myByHi(myThid) - myByLo(myThid) + 1
248     act3 = myThid - 1
249     max3 = nTx*nTy
250     act4 = ikey_dynamics - 1
251 heimbach 1.91 idynkey = (act1 + 1) + act2*max1
252 heimbach 1.76 & + act3*max1*max2
253     & + act4*max1*max2*max3
254     #endif /* ALLOW_AUTODIFF_TAMC */
255    
256     C-- Set up work arrays that need valid initial values
257 jmc 1.94 DO k=1,Nr
258     DO j=1-OLy,sNy+OLy
259     DO i=1-OLx,sNx+OLx
260 heimbach 1.87 KappaRU(i,j,k) = 0. _d 0
261     KappaRV(i,j,k) = 0. _d 0
262     ENDDO
263 jmc 1.94 ENDDO
264     ENDDO
265     DO j=1-OLy,sNy+OLy
266     DO i=1-OLx,sNx+OLx
267 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
268     fVerU (i,j,2) = 0. _d 0
269     fVerV (i,j,1) = 0. _d 0
270     fVerV (i,j,2) = 0. _d 0
271 jmc 1.94 phiHydF (i,j) = 0. _d 0
272     phiHydC (i,j) = 0. _d 0
273 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
274     dPhiHydY(i,j) = 0. _d 0
275 heimbach 1.76 ENDDO
276     ENDDO
277 heimbach 1.49
278 jmc 1.63 C-- Start computation of dynamics
279 jmc 1.93 iMin = 0
280     iMax = sNx+1
281     jMin = 0
282     jMax = sNy+1
283 jmc 1.63
284 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
285 heimbach 1.91 CADJ STORE wvel (:,:,:,bi,bj) =
286     CADJ & comlev1_bibj, key = idynkey, byte = isbyte
287 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
288    
289 jmc 1.65 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
290 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
291     IF (implicSurfPress.NE.1.) THEN
292 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
293     I bi,bj,iMin,iMax,jMin,jMax,
294     I etaN,
295     O phiSurfX,phiSurfY,
296     I myThid )
297 jmc 1.63 ENDIF
298 heimbach 1.83
299     #ifdef ALLOW_AUTODIFF_TAMC
300 heimbach 1.91 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
301     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
302 heimbach 1.83 #ifdef ALLOW_KPP
303     CADJ STORE KPPviscAz (:,:,:,bi,bj)
304 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
305 heimbach 1.83 #endif /* ALLOW_KPP */
306     #endif /* ALLOW_AUTODIFF_TAMC */
307 adcroft 1.58
308 heimbach 1.77 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
309     C-- Calculate the total vertical diffusivity
310     DO k=1,Nr
311     CALL CALC_VISCOSITY(
312     I bi,bj,iMin,iMax,jMin,jMax,k,
313     O KappaRU,KappaRV,
314     I myThid)
315     ENDDO
316     #endif
317    
318 adcroft 1.58 C-- Start of dynamics loop
319     DO k=1,Nr
320    
321     C-- km1 Points to level above k (=k-1)
322     C-- kup Cycles through 1,2 to point to layer above
323     C-- kDown Cycles through 2,1 to point to current layer
324    
325     km1 = MAX(1,k-1)
326 heimbach 1.77 kp1 = MIN(k+1,Nr)
327 adcroft 1.58 kup = 1+MOD(k+1,2)
328     kDown= 1+MOD(k,2)
329    
330 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
331 heimbach 1.91 kkey = (idynkey-1)*Nr + k
332 heimbach 1.95 CADJ STORE totphihyd (:,:,k,bi,bj)
333     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
334 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
335    
336 adcroft 1.58 C-- Integrate hydrostatic balance for phiHyd with BC of
337     C phiHyd(z=0)=0
338     C distinguishe between Stagger and Non Stagger time stepping
339     IF (staggerTimeStep) THEN
340     CALL CALC_PHI_HYD(
341     I bi,bj,iMin,iMax,jMin,jMax,k,
342 adcroft 1.81 I gT, gS,
343 jmc 1.94 U phiHydF,
344     O phiHydC, dPhiHydX, dPhiHydY,
345 jmc 1.92 I myTime, myIter, myThid )
346 adcroft 1.58 ELSE
347     CALL CALC_PHI_HYD(
348     I bi,bj,iMin,iMax,jMin,jMax,k,
349     I theta, salt,
350 jmc 1.94 U phiHydF,
351     O phiHydC, dPhiHydX, dPhiHydY,
352 jmc 1.92 I myTime, myIter, myThid )
353 adcroft 1.58 ENDIF
354 mlosch 1.89
355 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
356     C and step forward storing the result in gUnm1, gVnm1, etc...
357     IF ( momStepping ) THEN
358 adcroft 1.79 #ifndef DISABLE_MOM_FLUXFORM
359     IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
360 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
361 jmc 1.94 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
362 adcroft 1.58 U fVerU, fVerV,
363 adcroft 1.80 I myTime, myIter, myThid)
364 adcroft 1.79 #endif
365     #ifndef DISABLE_MOM_VECINV
366     IF (vectorInvariantMomentum) CALL MOM_VECINV(
367     I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
368 jmc 1.92 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
369 adcroft 1.79 U fVerU, fVerV,
370 adcroft 1.80 I myTime, myIter, myThid)
371 adcroft 1.79 #endif
372 adcroft 1.58 CALL TIMESTEP(
373 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
374 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
375 adcroft 1.58 I myIter, myThid)
376    
377     #ifdef ALLOW_OBCS
378     C-- Apply open boundary conditions
379     IF (useOBCS) THEN
380     CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
381     END IF
382     #endif /* ALLOW_OBCS */
383    
384     #ifdef ALLOW_AUTODIFF_TAMC
385     #ifdef INCLUDE_CD_CODE
386     ELSE
387     DO j=1-OLy,sNy+OLy
388     DO i=1-OLx,sNx+OLx
389     guCD(i,j,k,bi,bj) = 0.0
390     gvCD(i,j,k,bi,bj) = 0.0
391     END DO
392     END DO
393     #endif /* INCLUDE_CD_CODE */
394     #endif /* ALLOW_AUTODIFF_TAMC */
395     ENDIF
396    
397    
398     C-- end of dynamics k loop (1:Nr)
399     ENDDO
400    
401 adcroft 1.44 C-- Implicit viscosity
402 adcroft 1.58 IF (implicitViscosity.AND.momStepping) THEN
403     #ifdef ALLOW_AUTODIFF_TAMC
404 heimbach 1.91 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
405 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
406 adcroft 1.42 CALL IMPLDIFF(
407     I bi, bj, iMin, iMax, jMin, jMax,
408     I deltaTmom, KappaRU,recip_HFacW,
409     U gUNm1,
410     I myThid )
411 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
412 heimbach 1.91 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
413 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
414 adcroft 1.42 CALL IMPLDIFF(
415     I bi, bj, iMin, iMax, jMin, jMax,
416     I deltaTmom, KappaRV,recip_HFacS,
417     U gVNm1,
418     I myThid )
419 heimbach 1.49
420 adcroft 1.58 #ifdef ALLOW_OBCS
421     C-- Apply open boundary conditions
422     IF (useOBCS) THEN
423     DO K=1,Nr
424     CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
425     ENDDO
426     END IF
427     #endif /* ALLOW_OBCS */
428 heimbach 1.49
429 adcroft 1.58 #ifdef INCLUDE_CD_CODE
430     #ifdef ALLOW_AUTODIFF_TAMC
431 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
432 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
433 adcroft 1.42 CALL IMPLDIFF(
434     I bi, bj, iMin, iMax, jMin, jMax,
435     I deltaTmom, KappaRU,recip_HFacW,
436     U vVelD,
437     I myThid )
438 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
439 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
440 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
441 adcroft 1.42 CALL IMPLDIFF(
442     I bi, bj, iMin, iMax, jMin, jMax,
443     I deltaTmom, KappaRV,recip_HFacS,
444     U uVelD,
445     I myThid )
446 adcroft 1.58 #endif /* INCLUDE_CD_CODE */
447     C-- End If implicitViscosity.AND.momStepping
448 heimbach 1.53 ENDIF
449 cnh 1.1
450     ENDDO
451     ENDDO
452 mlosch 1.90
453     Cml(
454     C In order to compare the variance of phiHydLow of a p/z-coordinate
455     C run with etaH of a z/p-coordinate run the drift of phiHydLow
456     C has to be removed by something like the following subroutine:
457     C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
458     C & 'phiHydLow', myThid )
459     Cml)
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