/[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.161 - (hide annotations) (download)
Mon Mar 5 18:21:12 2012 UTC (12 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63k
Changes since 1.160: +29 -31 lines
update comments (calling tree)

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

  ViewVC Help
Powered by ViewVC 1.1.22