/[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.145 - (hide annotations) (download)
Wed Jan 20 03:50:56 2010 UTC (14 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62c, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62b
Changes since 1.144: +24 -1 lines
add pressure gradient diagnostics (without surface pressure contribution)

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

  ViewVC Help
Powered by ViewVC 1.1.22