/[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.168 - (hide annotations) (download)
Fri Jan 3 16:19:04 2014 UTC (10 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64s, checkpoint64t
Changes since 1.167: +26 -42 lines
remove some unnecessary TAF storage directives (note: not always removed,
 e.g., double storage of kappaRU,kappaRV in former version of dynamics.F)

1 jmc 1.168 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.167 2013/11/05 13:34:31 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.168 #ifdef ALLOW_AUTODIFF
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 c-- need some re-initialisation here to break dependencies
354 jmc 1.122 gU(i,j,k,bi,bj) = 0. _d 0
355     gV(i,j,k,bi,bj) = 0. _d 0
356 heimbach 1.87 ENDDO
357 jmc 1.94 ENDDO
358     ENDDO
359 jmc 1.168 #endif /* ALLOW_AUTODIFF */
360 jmc 1.94 DO j=1-OLy,sNy+OLy
361     DO i=1-OLx,sNx+OLx
362 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
363     fVerU (i,j,2) = 0. _d 0
364     fVerV (i,j,1) = 0. _d 0
365     fVerV (i,j,2) = 0. _d 0
366 jmc 1.136 phiHydF (i,j) = 0. _d 0
367     phiHydC (i,j) = 0. _d 0
368 jmc 1.146 #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
369 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
370 jmc 1.136 dPhiHydY(i,j) = 0. _d 0
371 jmc 1.146 #endif
372 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
373     phiSurfY(i,j) = 0. _d 0
374 jmc 1.110 guDissip(i,j) = 0. _d 0
375     gvDissip(i,j) = 0. _d 0
376 jmc 1.168 #ifdef ALLOW_AUTODIFF
377 gforget 1.143 phiHydLow(i,j,bi,bj) = 0. _d 0
378 jmc 1.150 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
379 heimbach 1.138 # ifndef DISABLE_RSTAR_CODE
380     dWtransC(i,j,bi,bj) = 0. _d 0
381     dWtransU(i,j,bi,bj) = 0. _d 0
382     dWtransV(i,j,bi,bj) = 0. _d 0
383     # endif
384     # endif
385 jmc 1.168 #endif /* ALLOW_AUTODIFF */
386 heimbach 1.76 ENDDO
387     ENDDO
388 heimbach 1.49
389 jmc 1.63 C-- Start computation of dynamics
390    
391 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
392 jmc 1.161 CADJ STORE wVel (:,:,:,bi,bj) =
393 heimbach 1.141 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
394 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
395    
396 jmc 1.161 C-- Explicit part of the Surface Potential Gradient (add in TIMESTEP)
397 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
398     IF (implicSurfPress.NE.1.) THEN
399 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
400     I bi,bj,iMin,iMax,jMin,jMax,
401     I etaN,
402     O phiSurfX,phiSurfY,
403 jmc 1.136 I myThid )
404 jmc 1.63 ENDIF
405 heimbach 1.83
406     #ifdef ALLOW_AUTODIFF_TAMC
407 jmc 1.161 CADJ STORE uVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
408     CADJ STORE vVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
409 heimbach 1.83 #ifdef ALLOW_KPP
410 jmc 1.150 CADJ STORE KPPviscAz (:,:,:,bi,bj)
411 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
412 heimbach 1.83 #endif /* ALLOW_KPP */
413     #endif /* ALLOW_AUTODIFF_TAMC */
414 adcroft 1.58
415 jmc 1.168 #ifndef ALLOW_AUTODIFF
416 jmc 1.165 IF ( .NOT.momViscosity ) THEN
417 jmc 1.168 #endif
418 jmc 1.165 DO k=1,Nr
419     DO j=1-OLy,sNy+OLy
420     DO i=1-OLx,sNx+OLx
421     KappaRU(i,j,k) = 0. _d 0
422     KappaRV(i,j,k) = 0. _d 0
423     ENDDO
424     ENDDO
425     ENDDO
426 jmc 1.168 #ifndef ALLOW_AUTODIFF
427     ENDIF
428     #endif
429 jmc 1.165 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
430 jmc 1.140 C-- Calculate the total vertical viscosity
431 jmc 1.165 IF ( momViscosity ) THEN
432     CALL CALC_VISCOSITY(
433 jmc 1.140 I bi,bj, iMin,iMax,jMin,jMax,
434     O KappaRU, KappaRV,
435     I myThid )
436 jmc 1.165 ENDIF
437     #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */
438 heimbach 1.77
439 jmc 1.167 #ifdef ALLOW_SMAG_3D
440     IF ( useSmag3D ) THEN
441     CALL MOM_CALC_3D_STRAIN(
442     O str11, str22, str33, str12, str13, str23,
443     I bi, bj, myThid )
444     ENDIF
445     #endif /* ALLOW_SMAG_3D */
446    
447 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
448 jmc 1.150 CADJ STORE KappaRU(:,:,:)
449 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
450 jmc 1.150 CADJ STORE KappaRV(:,:,:)
451 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
452 heimbach 1.101 #endif /* ALLOW_AUTODIFF_TAMC */
453    
454 mlosch 1.159 #ifdef ALLOW_OBCS
455     C-- For Stevens boundary conditions velocities need to be extrapolated
456     C (copied) to a narrow strip outside the domain
457 jmc 1.165 IF ( useOBCS ) THEN
458 jmc 1.160 CALL OBCS_COPY_UV_N(
459 jmc 1.161 U uVel(1-OLx,1-OLy,1,bi,bj),
460     U vVel(1-OLx,1-OLy,1,bi,bj),
461 mlosch 1.159 I Nr, bi, bj, myThid )
462 jmc 1.165 ENDIF
463 mlosch 1.159 #endif /* ALLOW_OBCS */
464    
465 m_bates 1.166 #ifdef ALLOW_EDDYPSI
466     CALL CALC_EDDY_STRESS(bi,bj,myThid)
467     #endif
468    
469 adcroft 1.58 C-- Start of dynamics loop
470     DO k=1,Nr
471    
472     C-- km1 Points to level above k (=k-1)
473     C-- kup Cycles through 1,2 to point to layer above
474     C-- kDown Cycles through 2,1 to point to current layer
475    
476     km1 = MAX(1,k-1)
477 heimbach 1.77 kp1 = MIN(k+1,Nr)
478 adcroft 1.58 kup = 1+MOD(k+1,2)
479     kDown= 1+MOD(k,2)
480    
481 jmc 1.144 #ifdef ALLOW_AUTODIFF_TAMC
482 heimbach 1.91 kkey = (idynkey-1)*Nr + k
483 jmc 1.161 CADJ STORE totPhiHyd (:,:,k,bi,bj)
484 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
485 jmc 1.150 CADJ STORE theta (:,:,k,bi,bj)
486 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
487 jmc 1.150 CADJ STORE salt (:,:,k,bi,bj)
488 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
489 jmc 1.161 CADJ STORE gT(:,:,k,bi,bj)
490 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
491 jmc 1.161 CADJ STORE gS(:,:,k,bi,bj)
492 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
493 heimbach 1.126 # ifdef NONLIN_FRSURF
494     cph-test
495 jmc 1.150 CADJ STORE phiHydC (:,:)
496 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
497 jmc 1.150 CADJ STORE phiHydF (:,:)
498 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
499 jmc 1.161 CADJ STORE gU(:,:,k,bi,bj)
500 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
501 jmc 1.161 CADJ STORE gV(:,:,k,bi,bj)
502 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
503 heimbach 1.148 # ifndef ALLOW_ADAMSBASHFORTH_3
504 jmc 1.161 CADJ STORE guNm1(:,:,k,bi,bj)
505 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
506 jmc 1.161 CADJ STORE gvNm1(:,:,k,bi,bj)
507 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
508 heimbach 1.148 # else
509 jmc 1.161 CADJ STORE guNm(:,:,k,bi,bj,1)
510 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
511 jmc 1.161 CADJ STORE guNm(:,:,k,bi,bj,2)
512 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
513 jmc 1.161 CADJ STORE gvNm(:,:,k,bi,bj,1)
514 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
515 jmc 1.161 CADJ STORE gvNm(:,:,k,bi,bj,2)
516 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
517     # endif
518 heimbach 1.126 # ifdef ALLOW_CD_CODE
519 jmc 1.161 CADJ STORE uNM1(:,:,k,bi,bj)
520 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
521 jmc 1.161 CADJ STORE vNM1(:,:,k,bi,bj)
522 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
523 jmc 1.150 CADJ STORE uVelD(:,:,k,bi,bj)
524 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
525 jmc 1.150 CADJ STORE vVelD(:,:,k,bi,bj)
526 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
527     # endif
528 jmc 1.168 # endif /* NONLIN_FRSURF */
529 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
530    
531 jmc 1.165 C-- Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
532 jmc 1.128 IF ( implicitIntGravWave ) THEN
533     CALL CALC_PHI_HYD(
534     I bi,bj,iMin,iMax,jMin,jMax,k,
535     I gT, gS,
536     U phiHydF,
537     O phiHydC, dPhiHydX, dPhiHydY,
538     I myTime, myIter, myThid )
539     ELSE
540     CALL CALC_PHI_HYD(
541 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
542     I theta, salt,
543 jmc 1.94 U phiHydF,
544     O phiHydC, dPhiHydX, dPhiHydY,
545 jmc 1.92 I myTime, myIter, myThid )
546 jmc 1.128 ENDIF
547 jmc 1.145 #ifdef ALLOW_DIAGNOSTICS
548     IF ( dPhiHydDiagIsOn ) THEN
549     tmpFac = -1. _d 0
550     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
551     & 'Um_dPHdx', k, 1, 2, bi, bj, myThid )
552     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
553     & 'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
554     ENDIF
555     #endif /* ALLOW_DIAGNOSTICS */
556 mlosch 1.89
557 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
558 jmc 1.96 C and step forward storing the result in gU, gV, etc...
559 adcroft 1.58 IF ( momStepping ) THEN
560 jmc 1.168 #ifdef ALLOW_AUTODIFF
561     DO j=1-OLy,sNy+OLy
562     DO i=1-OLx,sNx+OLx
563     guDissip(i,j) = 0. _d 0
564     gvDissip(i,j) = 0. _d 0
565     ENDDO
566     ENDDO
567     #endif /* ALLOW_AUTODIFF */
568 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
569 jmc 1.168 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
570     # ifndef DISABLE_RSTAR_CODE
571 jmc 1.150 CADJ STORE dWtransC(:,:,bi,bj)
572 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
573 jmc 1.150 CADJ STORE dWtransU(:,:,bi,bj)
574 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
575 jmc 1.150 CADJ STORE dWtransV(:,:,bi,bj)
576 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
577     # endif
578 jmc 1.168 # endif /* NONLIN_FRSURF and ALLOW_MOM_FLUXFORM */
579     # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
580     CADJ STORE fVerU(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
581     CADJ STORE fVerV(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
582     # endif
583 jmc 1.162 #endif /* ALLOW_AUTODIFF_TAMC */
584 heimbach 1.132 IF (.NOT. vectorInvariantMomentum) THEN
585 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
586 heimbach 1.132 CALL MOM_FLUXFORM(
587 jmc 1.162 I bi,bj,k,iMin,iMax,jMin,jMax,
588 jmc 1.121 I KappaRU, KappaRV,
589 jmc 1.162 U fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
590     O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
591 jmc 1.121 O guDissip, gvDissip,
592 adcroft 1.80 I myTime, myIter, myThid)
593 adcroft 1.79 #endif
594 heimbach 1.132 ELSE
595 edhill 1.105 #ifdef ALLOW_MOM_VECINV
596 heimbach 1.126 CALL MOM_VECINV(
597 jmc 1.162 I bi,bj,k,iMin,iMax,jMin,jMax,
598 jmc 1.121 I KappaRU, KappaRV,
599 jmc 1.162 I fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
600     O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
601 jmc 1.110 O guDissip, gvDissip,
602 adcroft 1.80 I myTime, myIter, myThid)
603 heimbach 1.132 #endif
604 heimbach 1.126 ENDIF
605 jmc 1.165
606 jmc 1.167 #ifdef ALLOW_SMAG_3D
607     IF ( useSmag3D ) THEN
608     CALL MOM_CALC_SMAG_3D(
609     I str11, str22, str33, str12, str13, str23,
610     O viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
611     I smag3D_hLsC, smag3D_hLsW, smag3D_hLsS, smag3D_hLsZ,
612     I k, bi, bj, myThid )
613     CALL MOM_UV_SMAG_3D(
614     I str11, str22, str12, str13, str23,
615     I viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
616     O addDissU, addDissV,
617     I iMin,iMax,jMin,jMax, k, bi, bj, myThid )
618     DO j= jMin,jMax
619     DO i= iMin,iMax
620     guDissip(i,j) = guDissip(i,j) + addDissU(i,j)
621     gvDissip(i,j) = gvDissip(i,j) + addDissV(i,j)
622     ENDDO
623     ENDDO
624     ENDIF
625     #endif /* ALLOW_SMAG_3D */
626    
627 adcroft 1.58 CALL TIMESTEP(
628 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
629 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
630 jmc 1.110 I guDissip, gvDissip,
631 jmc 1.96 I myTime, myIter, myThid)
632 adcroft 1.58
633     ENDIF
634    
635     C-- end of dynamics k loop (1:Nr)
636     ENDDO
637    
638 jmc 1.106 C-- Implicit Vertical advection & viscosity
639 gforget 1.147 #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
640     defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
641 jmc 1.106 IF ( momImplVertAdv ) THEN
642 jmc 1.136 CALL MOM_U_IMPLICIT_R( kappaRU,
643 jmc 1.106 I bi, bj, myTime, myIter, myThid )
644     CALL MOM_V_IMPLICIT_R( kappaRV,
645     I bi, bj, myTime, myIter, myThid )
646     ELSEIF ( implicitViscosity ) THEN
647     #else /* INCLUDE_IMPLVERTADV_CODE */
648     IF ( implicitViscosity ) THEN
649     #endif /* INCLUDE_IMPLVERTADV_CODE */
650 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
651 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
652 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
653 adcroft 1.42 CALL IMPLDIFF(
654     I bi, bj, iMin, iMax, jMin, jMax,
655 jmc 1.160 I -1, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
656 jmc 1.96 U gU,
657 adcroft 1.42 I myThid )
658 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
659 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
660 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
661 adcroft 1.42 CALL IMPLDIFF(
662     I bi, bj, iMin, iMax, jMin, jMax,
663 jmc 1.160 I -2, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
664 jmc 1.96 U gV,
665 adcroft 1.42 I myThid )
666 jmc 1.106 ENDIF
667 heimbach 1.49
668 mlosch 1.159 #ifdef ALLOW_OBCS
669 adcroft 1.58 C-- Apply open boundary conditions
670 jmc 1.151 IF ( useOBCS ) THEN
671 jmc 1.160 C-- but first save intermediate velocities to be used in the
672 mlosch 1.159 C next time step for the Stevens boundary conditions
673 jmc 1.160 CALL OBCS_SAVE_UV_N(
674     I bi, bj, iMin, iMax, jMin, jMax, 0,
675 mlosch 1.159 I gU, gV, myThid )
676 jmc 1.151 CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
677 jmc 1.106 ENDIF
678 mlosch 1.159 #endif /* ALLOW_OBCS */
679 heimbach 1.49
680 edhill 1.102 #ifdef ALLOW_CD_CODE
681 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
682 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
683 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
684 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
685 adcroft 1.42 CALL IMPLDIFF(
686     I bi, bj, iMin, iMax, jMin, jMax,
687 jmc 1.160 I 0, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
688 adcroft 1.42 U vVelD,
689     I myThid )
690 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
691 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
692 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
693 adcroft 1.42 CALL IMPLDIFF(
694     I bi, bj, iMin, iMax, jMin, jMax,
695 jmc 1.160 I 0, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
696 adcroft 1.42 U uVelD,
697     I myThid )
698 jmc 1.106 ENDIF
699 edhill 1.102 #endif /* ALLOW_CD_CODE */
700 jmc 1.106 C-- End implicit Vertical advection & viscosity
701 heimbach 1.109
702 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
703    
704 jmc 1.122 #ifdef ALLOW_NONHYDROSTATIC
705     C-- Step forward W field in N-H algorithm
706 jmc 1.136 IF ( nonHydrostatic ) THEN
707 jmc 1.122 #ifdef ALLOW_DEBUG
708 jmc 1.153 IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
709 jmc 1.122 #endif
710     CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
711 baylor 1.135 CALL CALC_GW(
712 jmc 1.136 I bi,bj, KappaRU, KappaRV,
713 jmc 1.167 I str13, str23, str33,
714     I viscAh3d_00, viscAh3d_13, viscAh3d_23,
715 jmc 1.136 I myTime, myIter, myThid )
716     ENDIF
717     IF ( nonHydrostatic.OR.implicitIntGravWave )
718     & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
719     IF ( nonHydrostatic )
720 jmc 1.128 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
721 jmc 1.122 #endif
722    
723     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
724    
725 jmc 1.136 C- end of bi,bj loops
726     ENDDO
727     ENDDO
728    
729     #ifdef ALLOW_OBCS
730 jmc 1.152 IF (useOBCS) THEN
731 jmc 1.155 CALL OBCS_EXCHANGES( myThid )
732 jmc 1.152 ENDIF
733 jmc 1.136 #endif
734    
735 mlosch 1.90 Cml(
736     C In order to compare the variance of phiHydLow of a p/z-coordinate
737     C run with etaH of a z/p-coordinate run the drift of phiHydLow
738     C has to be removed by something like the following subroutine:
739 jmc 1.144 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
740     C & 'phiHydLow', myTime, myThid )
741 mlosch 1.90 Cml)
742 adcroft 1.69
743 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
744 jmc 1.130 IF ( useDiagnostics ) THEN
745 jmc 1.113
746     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
747 jmc 1.120 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
748 molod 1.116
749 jmc 1.120 tmpFac = 1. _d 0
750     CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
751     & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
752 molod 1.116
753 jmc 1.120 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
754     & 'PHIBOTSQ',0, 1,0,1,1,myThid)
755 jmc 1.113
756     ENDIF
757     #endif /* ALLOW_DIAGNOSTICS */
758 jmc 1.136
759 edhill 1.104 #ifdef ALLOW_DEBUG
760 jmc 1.158 IF ( debugLevel .GE. debLevD ) THEN
761 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
762 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
763 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
764     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
765     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
766     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
767 jmc 1.115 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
768     CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
769     CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
770     CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
771     #ifndef ALLOW_ADAMSBASHFORTH_3
772     CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
773     CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
774     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
775     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
776     #endif
777 adcroft 1.70 ENDIF
778 adcroft 1.69 #endif
779 cnh 1.1
780 jmc 1.125 #ifdef DYNAMICS_GUGV_EXCH_CHECK
781 jmc 1.144 C- jmc: For safety checking only: This Exchange here should not change
782     C the solution. If solution changes, it means something is wrong,
783 jmc 1.125 C but it does not mean that it is less wrong with this exchange.
784 jmc 1.158 IF ( debugLevel .GE. debLevE ) THEN
785 jmc 1.125 CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
786     ENDIF
787     #endif
788    
789 jmc 1.123 #ifdef ALLOW_DEBUG
790 jmc 1.153 IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
791 jmc 1.123 #endif
792    
793 cnh 1.1 RETURN
794     END

  ViewVC Help
Powered by ViewVC 1.1.22