/[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.167 - (hide annotations) (download)
Tue Nov 5 13:34:31 2013 UTC (10 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64q, checkpoint64r
Changes since 1.166: +76 -7 lines
- move to pkg/mom_common and model/src (previously in tutorial_deep_convection
  code) 2nd version of isotropic 3-D Smagorinsky code interface: strain and
  viscosity are locally declared in dynmics.F and pass as argument to CALC_GW;

1 jmc 1.167 C $Header: /u/gcmpack/MITgcm/verification/tutorial_deep_convection/code/dynamics.F,v 1.1 2013/09/30 18:19:41 jmc Exp $
2 heimbach 1.78 C $Name: $
3 cnh 1.1
4 edhill 1.100 #include "PACKAGES_CONFIG.h"
5 adcroft 1.24 #include "CPP_OPTIONS.h"
6 jmc 1.167 #ifdef ALLOW_MOM_COMMON
7     # include "MOM_COMMON_OPTIONS.h"
8     #endif
9 heimbach 1.131 #ifdef ALLOW_OBCS
10     # include "OBCS_OPTIONS.h"
11     #endif
12    
13 jmc 1.125 #undef DYNAMICS_GUGV_EXCH_CHECK
14 cnh 1.1
15 cnh 1.82 CBOP
16     C !ROUTINE: DYNAMICS
17     C !INTERFACE:
18 cnh 1.8 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
19 cnh 1.82 C !DESCRIPTION: \bv
20     C *==========================================================*
21 jmc 1.144 C | SUBROUTINE DYNAMICS
22     C | o Controlling routine for the explicit part of the model
23     C | dynamics.
24 cnh 1.82 C *==========================================================*
25 jmc 1.144 C | This routine evaluates the "dynamics" terms for each
26     C | block of ocean in turn. Because the blocks of ocean have
27     C | overlap regions they are independent of one another.
28     C | If terms involving lateral integrals are needed in this
29     C | routine care will be needed. Similarly finite-difference
30     C | operations with stencils wider than the overlap region
31     C | require special consideration.
32 cnh 1.82 C | The algorithm...
33     C |
34     C | "Correction Step"
35     C | =================
36     C | Here we update the horizontal velocities with the surface
37     C | pressure such that the resulting flow is either consistent
38     C | with the free-surface evolution or the rigid-lid:
39     C | U[n] = U* + dt x d/dx P
40     C | V[n] = V* + dt x d/dy P
41 jmc 1.122 C | W[n] = W* + dt x d/dz P (NH mode)
42 cnh 1.82 C |
43     C | "Calculation of Gs"
44     C | ===================
45     C | This is where all the accelerations and tendencies (ie.
46     C | physics, parameterizations etc...) are calculated
47     C | rho = rho ( theta[n], salt[n] )
48     C | b = b(rho, theta)
49     C | K31 = K31 ( rho )
50     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
51     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
52     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
53     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
54     C |
55     C | "Time-stepping" or "Prediction"
56     C | ================================
57     C | The models variables are stepped forward with the appropriate
58     C | time-stepping scheme (currently we use Adams-Bashforth II)
59     C | - For momentum, the result is always *only* a "prediction"
60     C | in that the flow may be divergent and will be "corrected"
61     C | later with a surface pressure gradient.
62     C | - Normally for tracers the result is the new field at time
63     C | level [n+1} *BUT* in the case of implicit diffusion the result
64     C | is also *only* a prediction.
65     C | - We denote "predictors" with an asterisk (*).
66     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
67     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
68     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
70     C | With implicit diffusion:
71     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
72     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
73     C | (1 + dt * K * d_zz) theta[n] = theta*
74     C | (1 + dt * K * d_zz) salt[n] = salt*
75     C |
76     C *==========================================================*
77     C \ev
78     C !USES:
79 adcroft 1.40 IMPLICIT NONE
80 cnh 1.1 C == Global variables ===
81     #include "SIZE.h"
82     #include "EEPARAMS.h"
83 adcroft 1.6 #include "PARAMS.h"
84 jmc 1.165 #include "GRID.h"
85 adcroft 1.3 #include "DYNVARS.h"
86 jmc 1.167 #ifdef ALLOW_MOM_COMMON
87     # include "MOM_VISC.h"
88     #endif
89 edhill 1.103 #ifdef ALLOW_CD_CODE
90 jmc 1.165 # include "CD_CODE_VARS.h"
91 edhill 1.103 #endif
92 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
93 heimbach 1.53 # include "tamc.h"
94     # include "tamc_keys.h"
95 heimbach 1.67 # include "FFIELDS.h"
96 heimbach 1.91 # include "EOS.h"
97 heimbach 1.67 # ifdef ALLOW_KPP
98     # include "KPP.h"
99     # endif
100 heimbach 1.131 # ifdef ALLOW_PTRACERS
101     # include "PTRACERS_SIZE.h"
102 jmc 1.139 # include "PTRACERS_FIELDS.h"
103 heimbach 1.131 # endif
104     # ifdef ALLOW_OBCS
105 jmc 1.165 # include "OBCS_PARAMS.h"
106 jmc 1.157 # include "OBCS_FIELDS.h"
107 heimbach 1.131 # ifdef ALLOW_PTRACERS
108     # include "OBCS_PTRACERS.h"
109     # endif
110     # endif
111 heimbach 1.133 # ifdef ALLOW_MOM_FLUXFORM
112     # include "MOM_FLUXFORM.h"
113     # endif
114 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
115 jmc 1.62
116 cnh 1.82 C !CALLING SEQUENCE:
117     C DYNAMICS()
118     C |
119 jmc 1.122 C |-- CALC_EP_FORCING
120     C |
121 cnh 1.82 C |-- CALC_GRAD_PHI_SURF
122     C |
123     C |-- CALC_VISCOSITY
124     C |
125 jmc 1.167 C |-- MOM_CALC_3D_STRAIN
126     C |
127 m_bates 1.166 C |-- CALC_EDDY_STRESS
128     C |
129 jmc 1.136 C |-- CALC_PHI_HYD
130 cnh 1.82 C |
131 jmc 1.136 C |-- MOM_FLUXFORM
132 cnh 1.82 C |
133 jmc 1.136 C |-- MOM_VECINV
134 cnh 1.82 C |
135 jmc 1.167 C |-- MOM_CALC_SMAG_3D
136     C |-- MOM_UV_SMAG_3D
137     C |
138 jmc 1.136 C |-- TIMESTEP
139 cnh 1.82 C |
140 jmc 1.136 C |-- MOM_U_IMPLICIT_R
141     C |-- MOM_V_IMPLICIT_R
142 jmc 1.122 C |
143 jmc 1.136 C |-- IMPLDIFF
144 cnh 1.82 C |
145     C |-- OBCS_APPLY_UV
146     C |
147 jmc 1.122 C |-- CALC_GW
148     C |
149     C |-- DIAGNOSTICS_FILL
150     C |-- DEBUG_STATS_RL
151 cnh 1.82
152     C !INPUT/OUTPUT PARAMETERS:
153 cnh 1.1 C == Routine arguments ==
154 jmc 1.140 C myTime :: Current time in simulation
155     C myIter :: Current iteration number in simulation
156     C myThid :: Thread number for this instance of the routine.
157 cnh 1.8 _RL myTime
158     INTEGER myIter
159 adcroft 1.47 INTEGER myThid
160 cnh 1.1
161 jmc 1.145 C !FUNCTIONS:
162     #ifdef ALLOW_DIAGNOSTICS
163     LOGICAL DIAGNOSTICS_IS_ON
164     EXTERNAL DIAGNOSTICS_IS_ON
165     #endif
166    
167 cnh 1.82 C !LOCAL VARIABLES:
168 cnh 1.1 C == Local variables
169 jmc 1.113 C fVer[UV] o fVer: Vertical flux term - note fVer
170     C is "pipelined" in the vertical
171     C so we need an fVer for each
172     C variable.
173 jmc 1.94 C phiHydC :: hydrostatic potential anomaly at cell center
174     C In z coords phiHyd is the hydrostatic potential
175     C (=pressure/rho0) anomaly
176     C In p coords phiHyd is the geopotential height anomaly.
177     C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
178     C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
179     C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
180 jmc 1.92 C phiSurfY or geopotential (atmos) in X and Y direction
181 jmc 1.110 C guDissip :: dissipation tendency (all explicit terms), u component
182     C gvDissip :: dissipation tendency (all explicit terms), v component
183 jmc 1.165 C KappaRU :: vertical viscosity for velocity U-component
184     C KappaRV :: vertical viscosity for velocity V-component
185 jmc 1.162 C iMin, iMax :: Ranges and sub-block indices on which calculations
186     C jMin, jMax are applied.
187     C bi, bj :: tile indices
188     C k :: current level index
189     C km1, kp1 :: index of level above (k-1) and below (k+1)
190     C kUp, kDown :: Index for interface above and below. kUp and kDown are
191     C are switched with k to be the appropriate index into fVerU,V
192 cnh 1.30 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
193     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
194 jmc 1.94 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
195     _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
196 jmc 1.161 _RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
197     _RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
198 jmc 1.63 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
199     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
200 jmc 1.110 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
201     _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
202 jmc 1.161 _RL KappaRU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
203     _RL KappaRV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
204 jmc 1.167 #ifdef ALLOW_SMAG_3D
205     C str11 :: strain component Vxx @ grid-cell center
206     C str22 :: strain component Vyy @ grid-cell center
207     C str33 :: strain component Vzz @ grid-cell center
208     C str12 :: strain component Vxy @ grid-cell corner
209     C str13 :: strain component Vxz @ above uVel
210     C str23 :: strain component Vyz @ above vVel
211     C viscAh3d_00 :: Smagorinsky viscosity @ grid-cell center
212     C viscAh3d_12 :: Smagorinsky viscosity @ grid-cell corner
213     C viscAh3d_13 :: Smagorinsky viscosity @ above uVel
214     C viscAh3d_23 :: Smagorinsky viscosity @ above vVel
215     C addDissU :: zonal momentum tendency from 3-D Smag. viscosity
216     C addDissV :: merid momentum tendency from 3-D Smag. viscosity
217     _RL str11(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
218     _RL str22(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
219     _RL str33(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
220     _RL str12(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
221     _RL str13(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
222     _RL str23(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
223     _RL viscAh3d_00(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
224     _RL viscAh3d_12(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
225     _RL viscAh3d_13(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
226     _RL viscAh3d_23(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
227     _RL addDissU(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
228     _RL addDissV(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
229     #elif ( defined ALLOW_NONHYDROSTATIC )
230     _RL str13(1), str23(1), str33(1)
231     _RL viscAh3d_00(1), viscAh3d_13(1), viscAh3d_23(1)
232     #endif
233 adcroft 1.12
234 cnh 1.1 INTEGER bi, bj
235     INTEGER i, j
236 jmc 1.162 INTEGER k, km1, kp1, kUp, kDown
237 jmc 1.167 INTEGER iMin, iMax
238     INTEGER jMin, jMax
239     PARAMETER( iMin = 0 , iMax = sNx+1 )
240     PARAMETER( jMin = 0 , jMax = sNy+1 )
241 cnh 1.1
242 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
243 jmc 1.145 LOGICAL dPhiHydDiagIsOn
244 jmc 1.120 _RL tmpFac
245 jmc 1.113 #endif /* ALLOW_DIAGNOSTICS */
246    
247 adcroft 1.11 C--- The algorithm...
248     C
249     C "Correction Step"
250     C =================
251     C Here we update the horizontal velocities with the surface
252     C pressure such that the resulting flow is either consistent
253     C with the free-surface evolution or the rigid-lid:
254     C U[n] = U* + dt x d/dx P
255     C V[n] = V* + dt x d/dy P
256     C
257     C "Calculation of Gs"
258     C ===================
259     C This is where all the accelerations and tendencies (ie.
260 heimbach 1.53 C physics, parameterizations etc...) are calculated
261 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
262 cnh 1.27 C b = b(rho, theta)
263 adcroft 1.11 C K31 = K31 ( rho )
264 jmc 1.61 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
265     C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
266     C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
267     C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
268 adcroft 1.11 C
269 adcroft 1.12 C "Time-stepping" or "Prediction"
270 adcroft 1.11 C ================================
271     C The models variables are stepped forward with the appropriate
272     C time-stepping scheme (currently we use Adams-Bashforth II)
273     C - For momentum, the result is always *only* a "prediction"
274     C in that the flow may be divergent and will be "corrected"
275     C later with a surface pressure gradient.
276     C - Normally for tracers the result is the new field at time
277     C level [n+1} *BUT* in the case of implicit diffusion the result
278     C is also *only* a prediction.
279     C - We denote "predictors" with an asterisk (*).
280     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
281     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
282     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
283     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
284 adcroft 1.12 C With implicit diffusion:
285 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
286     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
287 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
288     C (1 + dt * K * d_zz) salt[n] = salt*
289 adcroft 1.11 C---
290 cnh 1.82 CEOP
291 adcroft 1.11
292 jmc 1.123 #ifdef ALLOW_DEBUG
293 jmc 1.153 IF (debugMode) CALL DEBUG_ENTER( 'DYNAMICS', myThid )
294 jmc 1.123 #endif
295    
296 jmc 1.145 #ifdef ALLOW_DIAGNOSTICS
297     dPhiHydDiagIsOn = .FALSE.
298     IF ( useDiagnostics )
299     & dPhiHydDiagIsOn = DIAGNOSTICS_IS_ON( 'Um_dPHdx', myThid )
300     & .OR. DIAGNOSTICS_IS_ON( 'Vm_dPHdy', myThid )
301     #endif
302    
303 jmc 1.165 C-- Call to routine for calculation of Eliassen-Palm-flux-forced
304     C U-tendency, if desired:
305 heimbach 1.88 #ifdef INCLUDE_EP_FORCING_CODE
306     CALL CALC_EP_FORCING(myThid)
307     #endif
308    
309 heimbach 1.154 #ifdef ALLOW_AUTODIFF_MONITOR_DIAG
310 jmc 1.161 CALL DUMMY_IN_DYNAMICS( myTime, myIter, myThid )
311 heimbach 1.154 #endif
312    
313 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
314     C-- HPF directive to help TAMC
315     CHPF$ INDEPENDENT
316     #endif /* ALLOW_AUTODIFF_TAMC */
317    
318 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
319 heimbach 1.76
320     #ifdef ALLOW_AUTODIFF_TAMC
321     C-- HPF directive to help TAMC
322     CHPF$ INDEPENDENT, NEW (fVerU,fVerV
323 jmc 1.94 CHPF$& ,phiHydF
324 heimbach 1.76 CHPF$& ,KappaRU,KappaRV
325     CHPF$& )
326     #endif /* ALLOW_AUTODIFF_TAMC */
327    
328 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
329 heimbach 1.76
330     #ifdef ALLOW_AUTODIFF_TAMC
331     act1 = bi - myBxLo(myThid)
332     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
333     act2 = bj - myByLo(myThid)
334     max2 = myByHi(myThid) - myByLo(myThid) + 1
335     act3 = myThid - 1
336     max3 = nTx*nTy
337     act4 = ikey_dynamics - 1
338 heimbach 1.91 idynkey = (act1 + 1) + act2*max1
339 heimbach 1.76 & + act3*max1*max2
340     & + act4*max1*max2*max3
341     #endif /* ALLOW_AUTODIFF_TAMC */
342    
343 heimbach 1.97 C-- Set up work arrays with valid (i.e. not NaN) values
344 jmc 1.161 C These initial values do not alter the numerical results. They
345 heimbach 1.97 C just ensure that all memory references are to valid floating
346     C point numbers. This prevents spurious hardware signals due to
347     C uninitialised but inert locations.
348    
349 jmc 1.140 #ifdef ALLOW_AUTODIFF_TAMC
350 jmc 1.94 DO k=1,Nr
351     DO j=1-OLy,sNy+OLy
352     DO i=1-OLx,sNx+OLx
353 heimbach 1.97 cph(
354     c-- need some re-initialisation here to break dependencies
355     cph)
356 jmc 1.122 gU(i,j,k,bi,bj) = 0. _d 0
357     gV(i,j,k,bi,bj) = 0. _d 0
358 heimbach 1.87 ENDDO
359 jmc 1.94 ENDDO
360     ENDDO
361 jmc 1.140 #endif /* ALLOW_AUTODIFF_TAMC */
362 jmc 1.94 DO j=1-OLy,sNy+OLy
363     DO i=1-OLx,sNx+OLx
364 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
365     fVerU (i,j,2) = 0. _d 0
366     fVerV (i,j,1) = 0. _d 0
367     fVerV (i,j,2) = 0. _d 0
368 jmc 1.136 phiHydF (i,j) = 0. _d 0
369     phiHydC (i,j) = 0. _d 0
370 jmc 1.146 #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
371 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
372 jmc 1.136 dPhiHydY(i,j) = 0. _d 0
373 jmc 1.146 #endif
374 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
375     phiSurfY(i,j) = 0. _d 0
376 jmc 1.110 guDissip(i,j) = 0. _d 0
377     gvDissip(i,j) = 0. _d 0
378 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
379 gforget 1.143 phiHydLow(i,j,bi,bj) = 0. _d 0
380 jmc 1.150 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
381 heimbach 1.138 # ifndef DISABLE_RSTAR_CODE
382 jahn 1.164 # ifndef ALLOW_AUTODIFF_OPENAD
383 heimbach 1.138 dWtransC(i,j,bi,bj) = 0. _d 0
384     dWtransU(i,j,bi,bj) = 0. _d 0
385     dWtransV(i,j,bi,bj) = 0. _d 0
386 jahn 1.164 # endif
387 heimbach 1.138 # endif
388     # endif
389     #endif
390 heimbach 1.76 ENDDO
391     ENDDO
392 heimbach 1.49
393 jmc 1.63 C-- Start computation of dynamics
394    
395 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
396 jmc 1.161 CADJ STORE wVel (:,:,:,bi,bj) =
397 heimbach 1.141 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
398 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
399    
400 jmc 1.161 C-- Explicit part of the Surface Potential Gradient (add in TIMESTEP)
401 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
402     IF (implicSurfPress.NE.1.) THEN
403 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
404     I bi,bj,iMin,iMax,jMin,jMax,
405     I etaN,
406     O phiSurfX,phiSurfY,
407 jmc 1.136 I myThid )
408 jmc 1.63 ENDIF
409 heimbach 1.83
410     #ifdef ALLOW_AUTODIFF_TAMC
411 jmc 1.161 CADJ STORE uVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
412     CADJ STORE vVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
413 heimbach 1.83 #ifdef ALLOW_KPP
414 jmc 1.150 CADJ STORE KPPviscAz (:,:,:,bi,bj)
415 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
416 heimbach 1.83 #endif /* ALLOW_KPP */
417     #endif /* ALLOW_AUTODIFF_TAMC */
418 adcroft 1.58
419 jmc 1.165 #if (defined INCLUDE_CALC_DIFFUSIVITY_CALL) && !(defined ALLOW_AUTODIFF)
420     IF ( .NOT.momViscosity ) THEN
421     #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL and not ALLOW_AUTODIFF */
422     DO k=1,Nr
423     DO j=1-OLy,sNy+OLy
424     DO i=1-OLx,sNx+OLx
425     KappaRU(i,j,k) = 0. _d 0
426     KappaRV(i,j,k) = 0. _d 0
427     ENDDO
428     ENDDO
429     ENDDO
430     #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
431 jmc 1.140 C-- Calculate the total vertical viscosity
432 jmc 1.165 #ifdef ALLOW_AUTODIFF
433     IF ( momViscosity ) THEN
434     #else
435     ELSE
436     #endif
437     CALL CALC_VISCOSITY(
438 jmc 1.140 I bi,bj, iMin,iMax,jMin,jMax,
439     O KappaRU, KappaRV,
440     I myThid )
441 jmc 1.165 ENDIF
442     #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */
443 heimbach 1.77
444 jmc 1.167 #ifdef ALLOW_SMAG_3D
445     IF ( useSmag3D ) THEN
446     CALL MOM_CALC_3D_STRAIN(
447     O str11, str22, str33, str12, str13, str23,
448     I bi, bj, myThid )
449     ENDIF
450     #endif /* ALLOW_SMAG_3D */
451    
452 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
453 jmc 1.150 CADJ STORE KappaRU(:,:,:)
454 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
455 jmc 1.150 CADJ STORE KappaRV(:,:,:)
456 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
457 heimbach 1.101 #endif /* ALLOW_AUTODIFF_TAMC */
458    
459 mlosch 1.159 #ifdef ALLOW_OBCS
460     C-- For Stevens boundary conditions velocities need to be extrapolated
461     C (copied) to a narrow strip outside the domain
462 jmc 1.165 IF ( useOBCS ) THEN
463 jmc 1.160 CALL OBCS_COPY_UV_N(
464 jmc 1.161 U uVel(1-OLx,1-OLy,1,bi,bj),
465     U vVel(1-OLx,1-OLy,1,bi,bj),
466 mlosch 1.159 I Nr, bi, bj, myThid )
467 jmc 1.165 ENDIF
468 mlosch 1.159 #endif /* ALLOW_OBCS */
469    
470 m_bates 1.166 #ifdef ALLOW_EDDYPSI
471     CALL CALC_EDDY_STRESS(bi,bj,myThid)
472     #endif
473    
474 adcroft 1.58 C-- Start of dynamics loop
475     DO k=1,Nr
476    
477     C-- km1 Points to level above k (=k-1)
478     C-- kup Cycles through 1,2 to point to layer above
479     C-- kDown Cycles through 2,1 to point to current layer
480    
481     km1 = MAX(1,k-1)
482 heimbach 1.77 kp1 = MIN(k+1,Nr)
483 adcroft 1.58 kup = 1+MOD(k+1,2)
484     kDown= 1+MOD(k,2)
485    
486 jmc 1.144 #ifdef ALLOW_AUTODIFF_TAMC
487 heimbach 1.91 kkey = (idynkey-1)*Nr + k
488 heimbach 1.99 c
489 jmc 1.161 CADJ STORE totPhiHyd (:,:,k,bi,bj)
490 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
491 jmc 1.161 CADJ STORE phiHydLow (:,:,bi,bj)
492 gforget 1.143 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
493 jmc 1.150 CADJ STORE theta (:,:,k,bi,bj)
494 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
495 jmc 1.150 CADJ STORE salt (:,:,k,bi,bj)
496 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
497 jmc 1.161 CADJ STORE gT(:,:,k,bi,bj)
498 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
499 jmc 1.161 CADJ STORE gS(:,:,k,bi,bj)
500 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
501 heimbach 1.126 # ifdef NONLIN_FRSURF
502     cph-test
503 jmc 1.150 CADJ STORE phiHydC (:,:)
504 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
505 jmc 1.150 CADJ STORE phiHydF (:,:)
506 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
507 jmc 1.161 CADJ STORE guDissip (:,:)
508 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
509 jmc 1.161 CADJ STORE gvDissip (:,:)
510 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
511 jmc 1.150 CADJ STORE fVerU (:,:,:)
512 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
513 jmc 1.150 CADJ STORE fVerV (:,:,:)
514 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
515 jmc 1.161 CADJ STORE gU(:,:,k,bi,bj)
516 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
517 jmc 1.161 CADJ STORE gV(:,:,k,bi,bj)
518 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
519 heimbach 1.148 # ifndef ALLOW_ADAMSBASHFORTH_3
520 jmc 1.161 CADJ STORE guNm1(:,:,k,bi,bj)
521 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
522 jmc 1.161 CADJ STORE gvNm1(:,:,k,bi,bj)
523 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
524 heimbach 1.148 # else
525 jmc 1.161 CADJ STORE guNm(:,:,k,bi,bj,1)
526 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
527 jmc 1.161 CADJ STORE guNm(:,:,k,bi,bj,2)
528 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
529 jmc 1.161 CADJ STORE gvNm(:,:,k,bi,bj,1)
530 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
531 jmc 1.161 CADJ STORE gvNm(:,:,k,bi,bj,2)
532 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
533     # endif
534 heimbach 1.126 # ifdef ALLOW_CD_CODE
535 jmc 1.161 CADJ STORE uNM1(:,:,k,bi,bj)
536 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
537 jmc 1.161 CADJ STORE vNM1(:,:,k,bi,bj)
538 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
539 jmc 1.150 CADJ STORE uVelD(:,:,k,bi,bj)
540 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
541 jmc 1.150 CADJ STORE vVelD(:,:,k,bi,bj)
542 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
543     # endif
544     # endif
545 heimbach 1.134 # ifdef ALLOW_DEPTH_CONTROL
546 jmc 1.150 CADJ STORE fVerU (:,:,:)
547 heimbach 1.134 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
548 jmc 1.150 CADJ STORE fVerV (:,:,:)
549 heimbach 1.134 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
550     # endif
551 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
552    
553 jmc 1.165 C-- Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
554 jmc 1.128 IF ( implicitIntGravWave ) THEN
555     CALL CALC_PHI_HYD(
556     I bi,bj,iMin,iMax,jMin,jMax,k,
557     I gT, gS,
558     U phiHydF,
559     O phiHydC, dPhiHydX, dPhiHydY,
560     I myTime, myIter, myThid )
561     ELSE
562     CALL CALC_PHI_HYD(
563 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
564     I theta, salt,
565 jmc 1.94 U phiHydF,
566     O phiHydC, dPhiHydX, dPhiHydY,
567 jmc 1.92 I myTime, myIter, myThid )
568 jmc 1.128 ENDIF
569 jmc 1.145 #ifdef ALLOW_DIAGNOSTICS
570     IF ( dPhiHydDiagIsOn ) THEN
571     tmpFac = -1. _d 0
572     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
573     & 'Um_dPHdx', k, 1, 2, bi, bj, myThid )
574     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
575     & 'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
576     ENDIF
577     #endif /* ALLOW_DIAGNOSTICS */
578 mlosch 1.89
579 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
580 jmc 1.96 C and step forward storing the result in gU, gV, etc...
581 adcroft 1.58 IF ( momStepping ) THEN
582 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
583 jmc 1.162 # ifdef NONLIN_FRSURF
584     # if (defined ALLOW_MOM_FLUXFORM) && !(defined DISABLE_RSTAR_CODE)
585 jmc 1.150 CADJ STORE dWtransC(:,:,bi,bj)
586 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
587 jmc 1.150 CADJ STORE dWtransU(:,:,bi,bj)
588 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
589 jmc 1.150 CADJ STORE dWtransV(:,:,bi,bj)
590 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
591     # endif
592 jmc 1.162 CADJ STORE fVerU(:,:,:)
593     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
594     CADJ STORE fVerV(:,:,:)
595     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
596     # endif /* NONLIN_FRSURF */
597     #endif /* ALLOW_AUTODIFF_TAMC */
598 heimbach 1.132 IF (.NOT. vectorInvariantMomentum) THEN
599 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
600 heimbach 1.132 CALL MOM_FLUXFORM(
601 jmc 1.162 I bi,bj,k,iMin,iMax,jMin,jMax,
602 jmc 1.121 I KappaRU, KappaRV,
603 jmc 1.162 U fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
604     O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
605 jmc 1.121 O guDissip, gvDissip,
606 adcroft 1.80 I myTime, myIter, myThid)
607 adcroft 1.79 #endif
608 heimbach 1.132 ELSE
609 edhill 1.105 #ifdef ALLOW_MOM_VECINV
610 heimbach 1.126 CALL MOM_VECINV(
611 jmc 1.162 I bi,bj,k,iMin,iMax,jMin,jMax,
612 jmc 1.121 I KappaRU, KappaRV,
613 jmc 1.162 I fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
614     O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
615 jmc 1.110 O guDissip, gvDissip,
616 adcroft 1.80 I myTime, myIter, myThid)
617 heimbach 1.132 #endif
618 heimbach 1.126 ENDIF
619 jmc 1.165
620 jmc 1.167 #ifdef ALLOW_SMAG_3D
621     IF ( useSmag3D ) THEN
622     CALL MOM_CALC_SMAG_3D(
623     I str11, str22, str33, str12, str13, str23,
624     O viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
625     I smag3D_hLsC, smag3D_hLsW, smag3D_hLsS, smag3D_hLsZ,
626     I k, bi, bj, myThid )
627     CALL MOM_UV_SMAG_3D(
628     I str11, str22, str12, str13, str23,
629     I viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
630     O addDissU, addDissV,
631     I iMin,iMax,jMin,jMax, k, bi, bj, myThid )
632     DO j= jMin,jMax
633     DO i= iMin,iMax
634     guDissip(i,j) = guDissip(i,j) + addDissU(i,j)
635     gvDissip(i,j) = gvDissip(i,j) + addDissV(i,j)
636     ENDDO
637     ENDDO
638     ENDIF
639     #endif /* ALLOW_SMAG_3D */
640    
641 adcroft 1.58 CALL TIMESTEP(
642 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
643 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
644 jmc 1.110 I guDissip, gvDissip,
645 jmc 1.96 I myTime, myIter, myThid)
646 adcroft 1.58
647     ENDIF
648    
649     C-- end of dynamics k loop (1:Nr)
650     ENDDO
651    
652 jmc 1.106 C-- Implicit Vertical advection & viscosity
653 gforget 1.147 #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
654     defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
655 jmc 1.106 IF ( momImplVertAdv ) THEN
656 jmc 1.136 CALL MOM_U_IMPLICIT_R( kappaRU,
657 jmc 1.106 I bi, bj, myTime, myIter, myThid )
658     CALL MOM_V_IMPLICIT_R( kappaRV,
659     I bi, bj, myTime, myIter, myThid )
660     ELSEIF ( implicitViscosity ) THEN
661     #else /* INCLUDE_IMPLVERTADV_CODE */
662     IF ( implicitViscosity ) THEN
663     #endif /* INCLUDE_IMPLVERTADV_CODE */
664 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
665 heimbach 1.101 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
666 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
667 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
668 adcroft 1.42 CALL IMPLDIFF(
669     I bi, bj, iMin, iMax, jMin, jMax,
670 jmc 1.160 I -1, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
671 jmc 1.96 U gU,
672 adcroft 1.42 I myThid )
673 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
674 heimbach 1.101 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
675 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
676 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
677 adcroft 1.42 CALL IMPLDIFF(
678     I bi, bj, iMin, iMax, jMin, jMax,
679 jmc 1.160 I -2, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
680 jmc 1.96 U gV,
681 adcroft 1.42 I myThid )
682 jmc 1.106 ENDIF
683 heimbach 1.49
684 mlosch 1.159 #ifdef ALLOW_OBCS
685 adcroft 1.58 C-- Apply open boundary conditions
686 jmc 1.151 IF ( useOBCS ) THEN
687 jmc 1.160 C-- but first save intermediate velocities to be used in the
688 mlosch 1.159 C next time step for the Stevens boundary conditions
689 jmc 1.160 CALL OBCS_SAVE_UV_N(
690     I bi, bj, iMin, iMax, jMin, jMax, 0,
691 mlosch 1.159 I gU, gV, myThid )
692 jmc 1.151 CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
693 jmc 1.106 ENDIF
694 mlosch 1.159 #endif /* ALLOW_OBCS */
695 heimbach 1.49
696 edhill 1.102 #ifdef ALLOW_CD_CODE
697 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
698 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
699 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
700 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
701 adcroft 1.42 CALL IMPLDIFF(
702     I bi, bj, iMin, iMax, jMin, jMax,
703 jmc 1.160 I 0, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
704 adcroft 1.42 U vVelD,
705     I myThid )
706 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
707 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
708 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
709 adcroft 1.42 CALL IMPLDIFF(
710     I bi, bj, iMin, iMax, jMin, jMax,
711 jmc 1.160 I 0, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
712 adcroft 1.42 U uVelD,
713     I myThid )
714 jmc 1.106 ENDIF
715 edhill 1.102 #endif /* ALLOW_CD_CODE */
716 jmc 1.106 C-- End implicit Vertical advection & viscosity
717 heimbach 1.109
718 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
719    
720 jmc 1.122 #ifdef ALLOW_NONHYDROSTATIC
721     C-- Step forward W field in N-H algorithm
722 jmc 1.136 IF ( nonHydrostatic ) THEN
723 jmc 1.122 #ifdef ALLOW_DEBUG
724 jmc 1.153 IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
725 jmc 1.122 #endif
726     CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
727 baylor 1.135 CALL CALC_GW(
728 jmc 1.136 I bi,bj, KappaRU, KappaRV,
729 jmc 1.167 I str13, str23, str33,
730     I viscAh3d_00, viscAh3d_13, viscAh3d_23,
731 jmc 1.136 I myTime, myIter, myThid )
732     ENDIF
733     IF ( nonHydrostatic.OR.implicitIntGravWave )
734     & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
735     IF ( nonHydrostatic )
736 jmc 1.128 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
737 jmc 1.122 #endif
738    
739     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
740    
741 jmc 1.136 C- end of bi,bj loops
742     ENDDO
743     ENDDO
744    
745     #ifdef ALLOW_OBCS
746 jmc 1.152 IF (useOBCS) THEN
747 jmc 1.155 CALL OBCS_EXCHANGES( myThid )
748 jmc 1.152 ENDIF
749 jmc 1.136 #endif
750    
751 mlosch 1.90 Cml(
752     C In order to compare the variance of phiHydLow of a p/z-coordinate
753     C run with etaH of a z/p-coordinate run the drift of phiHydLow
754     C has to be removed by something like the following subroutine:
755 jmc 1.144 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
756     C & 'phiHydLow', myTime, myThid )
757 mlosch 1.90 Cml)
758 adcroft 1.69
759 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
760 jmc 1.130 IF ( useDiagnostics ) THEN
761 jmc 1.113
762     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
763 jmc 1.120 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
764 molod 1.116
765 jmc 1.120 tmpFac = 1. _d 0
766     CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
767     & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
768 molod 1.116
769 jmc 1.120 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
770     & 'PHIBOTSQ',0, 1,0,1,1,myThid)
771 jmc 1.113
772     ENDIF
773     #endif /* ALLOW_DIAGNOSTICS */
774 jmc 1.136
775 edhill 1.104 #ifdef ALLOW_DEBUG
776 jmc 1.158 IF ( debugLevel .GE. debLevD ) THEN
777 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
778 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
779 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
780     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
781     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
782     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
783 jmc 1.115 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
784     CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
785     CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
786     CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
787     #ifndef ALLOW_ADAMSBASHFORTH_3
788     CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
789     CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
790     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
791     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
792     #endif
793 adcroft 1.70 ENDIF
794 adcroft 1.69 #endif
795 cnh 1.1
796 jmc 1.125 #ifdef DYNAMICS_GUGV_EXCH_CHECK
797 jmc 1.144 C- jmc: For safety checking only: This Exchange here should not change
798     C the solution. If solution changes, it means something is wrong,
799 jmc 1.125 C but it does not mean that it is less wrong with this exchange.
800 jmc 1.158 IF ( debugLevel .GE. debLevE ) THEN
801 jmc 1.125 CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
802     ENDIF
803     #endif
804    
805 jmc 1.123 #ifdef ALLOW_DEBUG
806 jmc 1.153 IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
807 jmc 1.123 #endif
808    
809 cnh 1.1 RETURN
810     END

  ViewVC Help
Powered by ViewVC 1.1.22