/[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.150 - (hide annotations) (download)
Mon Oct 4 02:58:03 2010 UTC (13 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62m, checkpoint62l
Changes since 1.149: +37 -37 lines
fix previous modif (from Sep 19, cvs version 1.149, "Separate MOM_FLUXFORM
 store directives") which was causing extensive recomputations in AD test
experiment global_ocean.cs32x15.

1 jmc 1.150 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.149 2010/09/29 19:01:31 heimbach 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.146 #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
328 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
329 jmc 1.136 dPhiHydY(i,j) = 0. _d 0
330 jmc 1.146 #endif
331 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
332     phiSurfY(i,j) = 0. _d 0
333 jmc 1.110 guDissip(i,j) = 0. _d 0
334     gvDissip(i,j) = 0. _d 0
335 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
336 gforget 1.143 phiHydLow(i,j,bi,bj) = 0. _d 0
337 jmc 1.150 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
338 heimbach 1.138 # ifndef DISABLE_RSTAR_CODE
339     dWtransC(i,j,bi,bj) = 0. _d 0
340     dWtransU(i,j,bi,bj) = 0. _d 0
341     dWtransV(i,j,bi,bj) = 0. _d 0
342     # endif
343     # endif
344     #endif
345 heimbach 1.76 ENDDO
346     ENDDO
347 heimbach 1.49
348 jmc 1.63 C-- Start computation of dynamics
349 jmc 1.93 iMin = 0
350     iMax = sNx+1
351     jMin = 0
352     jMax = sNy+1
353 jmc 1.63
354 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
355 jmc 1.150 CADJ STORE wvel (:,:,:,bi,bj) =
356 heimbach 1.141 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
357 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
358    
359 jmc 1.65 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
360 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
361     IF (implicSurfPress.NE.1.) THEN
362 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
363     I bi,bj,iMin,iMax,jMin,jMax,
364     I etaN,
365     O phiSurfX,phiSurfY,
366 jmc 1.136 I myThid )
367 jmc 1.63 ENDIF
368 heimbach 1.83
369     #ifdef ALLOW_AUTODIFF_TAMC
370 heimbach 1.91 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
371     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
372 heimbach 1.83 #ifdef ALLOW_KPP
373 jmc 1.150 CADJ STORE KPPviscAz (:,:,:,bi,bj)
374 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
375 heimbach 1.83 #endif /* ALLOW_KPP */
376     #endif /* ALLOW_AUTODIFF_TAMC */
377 adcroft 1.58
378 heimbach 1.77 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
379 jmc 1.140 C-- Calculate the total vertical viscosity
380     CALL CALC_VISCOSITY(
381     I bi,bj, iMin,iMax,jMin,jMax,
382     O KappaRU, KappaRV,
383     I myThid )
384     #else
385 heimbach 1.77 DO k=1,Nr
386 jmc 1.140 DO j=1-OLy,sNy+OLy
387     DO i=1-OLx,sNx+OLx
388     KappaRU(i,j,k) = 0. _d 0
389     KappaRV(i,j,k) = 0. _d 0
390     ENDDO
391     ENDDO
392     ENDDO
393 heimbach 1.77 #endif
394    
395 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
396 jmc 1.150 CADJ STORE KappaRU(:,:,:)
397 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
398 jmc 1.150 CADJ STORE KappaRV(:,:,:)
399 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
400 heimbach 1.101 #endif /* ALLOW_AUTODIFF_TAMC */
401    
402 adcroft 1.58 C-- Start of dynamics loop
403     DO k=1,Nr
404    
405     C-- km1 Points to level above k (=k-1)
406     C-- kup Cycles through 1,2 to point to layer above
407     C-- kDown Cycles through 2,1 to point to current layer
408    
409     km1 = MAX(1,k-1)
410 heimbach 1.77 kp1 = MIN(k+1,Nr)
411 adcroft 1.58 kup = 1+MOD(k+1,2)
412     kDown= 1+MOD(k,2)
413    
414 jmc 1.144 #ifdef ALLOW_AUTODIFF_TAMC
415 heimbach 1.91 kkey = (idynkey-1)*Nr + k
416 heimbach 1.99 c
417 jmc 1.150 CADJ STORE totphihyd (:,:,k,bi,bj)
418 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
419 gforget 1.143 CADJ STORE phihydlow (:,:,bi,bj)
420     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
421 jmc 1.150 CADJ STORE theta (:,:,k,bi,bj)
422 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
423 jmc 1.150 CADJ STORE salt (:,:,k,bi,bj)
424 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
425 jmc 1.150 CADJ STORE gt(:,:,k,bi,bj)
426 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
427 jmc 1.150 CADJ STORE gs(:,:,k,bi,bj)
428 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
429 heimbach 1.126 # ifdef NONLIN_FRSURF
430     cph-test
431 jmc 1.150 CADJ STORE phiHydC (:,:)
432 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
433 jmc 1.150 CADJ STORE phiHydF (:,:)
434 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
435 jmc 1.150 CADJ STORE gudissip (:,:)
436 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
437 jmc 1.150 CADJ STORE gvdissip (:,:)
438 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
439 jmc 1.150 CADJ STORE fVerU (:,:,:)
440 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
441 jmc 1.150 CADJ STORE fVerV (:,:,:)
442 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
443 jmc 1.150 CADJ STORE gu(:,:,k,bi,bj)
444 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
445 jmc 1.150 CADJ STORE gv(:,:,k,bi,bj)
446 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
447 heimbach 1.148 # ifndef ALLOW_ADAMSBASHFORTH_3
448 jmc 1.150 CADJ STORE gunm1(:,:,k,bi,bj)
449 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
450 jmc 1.150 CADJ STORE gvnm1(:,:,k,bi,bj)
451 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
452 heimbach 1.148 # else
453 jmc 1.150 CADJ STORE gunm(:,:,k,bi,bj,1)
454 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
455 jmc 1.150 CADJ STORE gunm(:,:,k,bi,bj,2)
456 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
457 jmc 1.150 CADJ STORE gvnm(:,:,k,bi,bj,1)
458 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
459 jmc 1.150 CADJ STORE gvnm(:,:,k,bi,bj,2)
460 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
461     # endif
462 heimbach 1.126 # ifdef ALLOW_CD_CODE
463 jmc 1.150 CADJ STORE unm1(:,:,k,bi,bj)
464 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
465 jmc 1.150 CADJ STORE vnm1(:,:,k,bi,bj)
466 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
467 jmc 1.150 CADJ STORE uVelD(:,:,k,bi,bj)
468 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
469 jmc 1.150 CADJ STORE vVelD(:,:,k,bi,bj)
470 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
471     # endif
472     # endif
473 heimbach 1.134 # ifdef ALLOW_DEPTH_CONTROL
474 jmc 1.150 CADJ STORE fVerU (:,:,:)
475 heimbach 1.134 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
476 jmc 1.150 CADJ STORE fVerV (:,:,:)
477 heimbach 1.134 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
478     # endif
479 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
480    
481 jmc 1.144 C-- Integrate hydrostatic balance for phiHyd with BC of
482 adcroft 1.58 C phiHyd(z=0)=0
483 jmc 1.128 IF ( implicitIntGravWave ) THEN
484     CALL CALC_PHI_HYD(
485     I bi,bj,iMin,iMax,jMin,jMax,k,
486     I gT, gS,
487     U phiHydF,
488     O phiHydC, dPhiHydX, dPhiHydY,
489     I myTime, myIter, myThid )
490     ELSE
491     CALL CALC_PHI_HYD(
492 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
493     I theta, salt,
494 jmc 1.94 U phiHydF,
495     O phiHydC, dPhiHydX, dPhiHydY,
496 jmc 1.92 I myTime, myIter, myThid )
497 jmc 1.128 ENDIF
498 jmc 1.145 #ifdef ALLOW_DIAGNOSTICS
499     IF ( dPhiHydDiagIsOn ) THEN
500     tmpFac = -1. _d 0
501     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
502     & 'Um_dPHdx', k, 1, 2, bi, bj, myThid )
503     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
504     & 'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
505     ENDIF
506     #endif /* ALLOW_DIAGNOSTICS */
507 mlosch 1.89
508 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
509 jmc 1.96 C and step forward storing the result in gU, gV, etc...
510 adcroft 1.58 IF ( momStepping ) THEN
511 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
512 jmc 1.150 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
513 heimbach 1.138 # ifndef DISABLE_RSTAR_CODE
514 jmc 1.150 CADJ STORE dWtransC(:,:,bi,bj)
515 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
516 jmc 1.150 CADJ STORE dWtransU(:,:,bi,bj)
517 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
518 jmc 1.150 CADJ STORE dWtransV(:,:,bi,bj)
519 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
520     # endif
521     # endif
522     #endif
523 heimbach 1.132 IF (.NOT. vectorInvariantMomentum) THEN
524 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
525 heimbach 1.132 C
526     CALL MOM_FLUXFORM(
527 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
528 jmc 1.121 I KappaRU, KappaRV,
529 adcroft 1.58 U fVerU, fVerV,
530 jmc 1.121 O guDissip, gvDissip,
531 adcroft 1.80 I myTime, myIter, myThid)
532 adcroft 1.79 #endif
533 heimbach 1.132 ELSE
534 edhill 1.105 #ifdef ALLOW_MOM_VECINV
535 heimbach 1.126 C
536 jmc 1.144 # ifdef ALLOW_AUTODIFF_TAMC
537 heimbach 1.126 # ifdef NONLIN_FRSURF
538 jmc 1.150 CADJ STORE fVerU(:,:,:)
539 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
540 jmc 1.150 CADJ STORE fVerV(:,:,:)
541 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
542     # endif
543     # endif /* ALLOW_AUTODIFF_TAMC */
544     C
545     CALL MOM_VECINV(
546 adcroft 1.79 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
547 jmc 1.121 I KappaRU, KappaRV,
548 adcroft 1.79 U fVerU, fVerV,
549 jmc 1.110 O guDissip, gvDissip,
550 adcroft 1.80 I myTime, myIter, myThid)
551 heimbach 1.132 #endif
552 heimbach 1.126 ENDIF
553 heimbach 1.132 C
554 adcroft 1.58 CALL TIMESTEP(
555 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
556 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
557 jmc 1.110 I guDissip, gvDissip,
558 jmc 1.96 I myTime, myIter, myThid)
559 adcroft 1.58
560     #ifdef ALLOW_OBCS
561     C-- Apply open boundary conditions
562 jmc 1.96 IF (useOBCS) THEN
563     CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
564     ENDIF
565 adcroft 1.58 #endif /* ALLOW_OBCS */
566    
567     ENDIF
568    
569    
570     C-- end of dynamics k loop (1:Nr)
571     ENDDO
572    
573 jmc 1.106 C-- Implicit Vertical advection & viscosity
574 gforget 1.147 #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
575     defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
576 jmc 1.106 IF ( momImplVertAdv ) THEN
577 jmc 1.136 CALL MOM_U_IMPLICIT_R( kappaRU,
578 jmc 1.106 I bi, bj, myTime, myIter, myThid )
579     CALL MOM_V_IMPLICIT_R( kappaRV,
580     I bi, bj, myTime, myIter, myThid )
581     ELSEIF ( implicitViscosity ) THEN
582     #else /* INCLUDE_IMPLVERTADV_CODE */
583     IF ( implicitViscosity ) THEN
584     #endif /* INCLUDE_IMPLVERTADV_CODE */
585 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
586 heimbach 1.101 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
587 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
588 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
589 adcroft 1.42 CALL IMPLDIFF(
590     I bi, bj, iMin, iMax, jMin, jMax,
591 jmc 1.124 I -1, KappaRU,recip_HFacW,
592 jmc 1.96 U gU,
593 adcroft 1.42 I myThid )
594 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
595 heimbach 1.101 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
596 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
597 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
598 adcroft 1.42 CALL IMPLDIFF(
599     I bi, bj, iMin, iMax, jMin, jMax,
600 jmc 1.124 I -2, KappaRV,recip_HFacS,
601 jmc 1.96 U gV,
602 adcroft 1.42 I myThid )
603 jmc 1.106 ENDIF
604 heimbach 1.49
605 adcroft 1.58 #ifdef ALLOW_OBCS
606     C-- Apply open boundary conditions
607 jmc 1.106 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
608 adcroft 1.58 DO K=1,Nr
609 jmc 1.96 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
610 adcroft 1.58 ENDDO
611 jmc 1.106 ENDIF
612 adcroft 1.58 #endif /* ALLOW_OBCS */
613 heimbach 1.49
614 edhill 1.102 #ifdef ALLOW_CD_CODE
615 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
616 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
617 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
618 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
619 adcroft 1.42 CALL IMPLDIFF(
620     I bi, bj, iMin, iMax, jMin, jMax,
621 jmc 1.111 I 0, KappaRU,recip_HFacW,
622 adcroft 1.42 U vVelD,
623     I myThid )
624 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
625 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
626 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
627 adcroft 1.42 CALL IMPLDIFF(
628     I bi, bj, iMin, iMax, jMin, jMax,
629 jmc 1.111 I 0, KappaRV,recip_HFacS,
630 adcroft 1.42 U uVelD,
631     I myThid )
632 jmc 1.106 ENDIF
633 edhill 1.102 #endif /* ALLOW_CD_CODE */
634 jmc 1.106 C-- End implicit Vertical advection & viscosity
635 heimbach 1.109
636 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
637    
638 jmc 1.122 #ifdef ALLOW_NONHYDROSTATIC
639     C-- Step forward W field in N-H algorithm
640 jmc 1.136 IF ( nonHydrostatic ) THEN
641 jmc 1.122 #ifdef ALLOW_DEBUG
642 jmc 1.123 IF ( debugLevel .GE. debLevB )
643     & CALL DEBUG_CALL('CALC_GW', myThid )
644 jmc 1.122 #endif
645     CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
646 baylor 1.135 CALL CALC_GW(
647 jmc 1.136 I bi,bj, KappaRU, KappaRV,
648     I myTime, myIter, myThid )
649     ENDIF
650     IF ( nonHydrostatic.OR.implicitIntGravWave )
651     & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
652     IF ( nonHydrostatic )
653 jmc 1.128 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
654 jmc 1.122 #endif
655    
656     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
657    
658 jmc 1.136 C- end of bi,bj loops
659     ENDDO
660     ENDDO
661    
662     #ifdef ALLOW_OBCS
663     IF (useOBCS) THEN
664     CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
665     ENDIF
666     #endif
667    
668 mlosch 1.90 Cml(
669     C In order to compare the variance of phiHydLow of a p/z-coordinate
670     C run with etaH of a z/p-coordinate run the drift of phiHydLow
671     C has to be removed by something like the following subroutine:
672 jmc 1.144 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
673     C & 'phiHydLow', myTime, myThid )
674 mlosch 1.90 Cml)
675 adcroft 1.69
676 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
677 jmc 1.130 IF ( useDiagnostics ) THEN
678 jmc 1.113
679     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
680 jmc 1.120 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
681 molod 1.116
682 jmc 1.120 tmpFac = 1. _d 0
683     CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
684     & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
685 molod 1.116
686 jmc 1.120 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
687     & 'PHIBOTSQ',0, 1,0,1,1,myThid)
688 jmc 1.113
689     ENDIF
690     #endif /* ALLOW_DIAGNOSTICS */
691 jmc 1.136
692 edhill 1.104 #ifdef ALLOW_DEBUG
693 heimbach 1.98 If ( debugLevel .GE. debLevB ) THEN
694 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
695 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
696 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
697     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
698     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
699     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
700 jmc 1.115 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
701     CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
702     CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
703     CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
704     #ifndef ALLOW_ADAMSBASHFORTH_3
705     CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
706     CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
707     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
708     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
709     #endif
710 adcroft 1.70 ENDIF
711 adcroft 1.69 #endif
712 cnh 1.1
713 jmc 1.125 #ifdef DYNAMICS_GUGV_EXCH_CHECK
714 jmc 1.144 C- jmc: For safety checking only: This Exchange here should not change
715     C the solution. If solution changes, it means something is wrong,
716 jmc 1.125 C but it does not mean that it is less wrong with this exchange.
717     IF ( debugLevel .GT. debLevB ) THEN
718     CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
719     ENDIF
720     #endif
721    
722 jmc 1.123 #ifdef ALLOW_DEBUG
723     IF ( debugLevel .GE. debLevB )
724     & CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
725     #endif
726    
727 cnh 1.1 RETURN
728     END

  ViewVC Help
Powered by ViewVC 1.1.22