/[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.176 - (hide annotations) (download)
Wed Dec 24 19:09:33 2014 UTC (9 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint65o
Changes since 1.175: +19 -19 lines
- add one more level to vertical viscosity local arrays (Nr+1, previously Nr)
  since no-slip bottom BC uses viscosity @ k+1 to update velocity @ level k
- for now and until vertical mixing scheme are updated to fill up Nr+1 level,
  just copy Nr value to Nr+1

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

  ViewVC Help
Powered by ViewVC 1.1.22