/[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.166 - (hide annotations) (download)
Sun Sep 15 14:28:31 2013 UTC (10 years, 9 months ago) by m_bates
Branch: MAIN
CVS Tags: checkpoint64p, checkpoint64o
Changes since 1.165: +7 -1 lines
Changes to implement a residual model.  Also, calculation of the mean velocity from the residual and bolus.

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

  ViewVC Help
Powered by ViewVC 1.1.22