/[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.164 - (hide annotations) (download)
Thu Mar 21 18:15:44 2013 UTC (11 years, 2 months ago) by jahn
Branch: MAIN
CVS Tags: checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64g, checkpoint64f
Changes since 1.163: +3 -1 lines
revert to local variables for OpenAD

1 jahn 1.164 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.163 2012/11/15 15:55:42 dimitri 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 dimitri 1.163 #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 jmc 1.136 C |-- CALC_PHI_HYD
120 cnh 1.82 C |
121 jmc 1.136 C |-- MOM_FLUXFORM
122 cnh 1.82 C |
123 jmc 1.136 C |-- MOM_VECINV
124 cnh 1.82 C |
125 jmc 1.136 C |-- TIMESTEP
126 cnh 1.82 C |
127 jmc 1.136 C |-- MOM_U_IMPLICIT_R
128     C |-- MOM_V_IMPLICIT_R
129 jmc 1.122 C |
130 jmc 1.136 C |-- IMPLDIFF
131 cnh 1.82 C |
132     C |-- OBCS_APPLY_UV
133     C |
134 jmc 1.122 C |-- CALC_GW
135     C |
136     C |-- DIAGNOSTICS_FILL
137     C |-- DEBUG_STATS_RL
138 cnh 1.82
139     C !INPUT/OUTPUT PARAMETERS:
140 cnh 1.1 C == Routine arguments ==
141 jmc 1.140 C myTime :: Current time in simulation
142     C myIter :: Current iteration number in simulation
143     C myThid :: Thread number for this instance of the routine.
144 cnh 1.8 _RL myTime
145     INTEGER myIter
146 adcroft 1.47 INTEGER myThid
147 cnh 1.1
148 jmc 1.145 C !FUNCTIONS:
149     #ifdef ALLOW_DIAGNOSTICS
150     LOGICAL DIAGNOSTICS_IS_ON
151     EXTERNAL DIAGNOSTICS_IS_ON
152     #endif
153    
154 cnh 1.82 C !LOCAL VARIABLES:
155 cnh 1.1 C == Local variables
156 jmc 1.113 C fVer[UV] o fVer: Vertical flux term - note fVer
157     C is "pipelined" in the vertical
158     C so we need an fVer for each
159     C variable.
160 jmc 1.94 C phiHydC :: hydrostatic potential anomaly at cell center
161     C In z coords phiHyd is the hydrostatic potential
162     C (=pressure/rho0) anomaly
163     C In p coords phiHyd is the geopotential height anomaly.
164     C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
165     C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
166     C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
167 jmc 1.92 C phiSurfY or geopotential (atmos) in X and Y direction
168 jmc 1.110 C guDissip :: dissipation tendency (all explicit terms), u component
169     C gvDissip :: dissipation tendency (all explicit terms), v component
170 jmc 1.140 C KappaRU :: vertical viscosity
171     C KappaRV :: vertical viscosity
172 jmc 1.162 C iMin, iMax :: Ranges and sub-block indices on which calculations
173     C jMin, jMax are applied.
174     C bi, bj :: tile indices
175     C k :: current level index
176     C km1, kp1 :: index of level above (k-1) and below (k+1)
177     C kUp, kDown :: Index for interface above and below. kUp and kDown are
178     C are switched with k to be the appropriate index into fVerU,V
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.161 _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 jmc 1.161 _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 jmc 1.162 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 jmc 1.153 IF (debugMode) CALL DEBUG_ENTER( 'DYNAMICS', myThid )
251 jmc 1.123 #endif
252    
253 jmc 1.145 #ifdef ALLOW_DIAGNOSTICS
254     dPhiHydDiagIsOn = .FALSE.
255     IF ( useDiagnostics )
256     & dPhiHydDiagIsOn = DIAGNOSTICS_IS_ON( 'Um_dPHdx', myThid )
257     & .OR. DIAGNOSTICS_IS_ON( 'Vm_dPHdy', myThid )
258     #endif
259    
260 heimbach 1.88 C-- Call to routine for calculation of
261     C Eliassen-Palm-flux-forced U-tendency,
262     C if desired:
263     #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.87 KappaRU(i,j,k) = 0. _d 0
312     KappaRV(i,j,k) = 0. _d 0
313 heimbach 1.97 cph(
314     c-- need some re-initialisation here to break dependencies
315     cph)
316 jmc 1.122 gU(i,j,k,bi,bj) = 0. _d 0
317     gV(i,j,k,bi,bj) = 0. _d 0
318 heimbach 1.87 ENDDO
319 jmc 1.94 ENDDO
320     ENDDO
321 jmc 1.140 #endif /* ALLOW_AUTODIFF_TAMC */
322 jmc 1.94 DO j=1-OLy,sNy+OLy
323     DO i=1-OLx,sNx+OLx
324 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
325     fVerU (i,j,2) = 0. _d 0
326     fVerV (i,j,1) = 0. _d 0
327     fVerV (i,j,2) = 0. _d 0
328 jmc 1.136 phiHydF (i,j) = 0. _d 0
329     phiHydC (i,j) = 0. _d 0
330 jmc 1.146 #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
331 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
332 jmc 1.136 dPhiHydY(i,j) = 0. _d 0
333 jmc 1.146 #endif
334 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
335     phiSurfY(i,j) = 0. _d 0
336 jmc 1.110 guDissip(i,j) = 0. _d 0
337     gvDissip(i,j) = 0. _d 0
338 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
339 gforget 1.143 phiHydLow(i,j,bi,bj) = 0. _d 0
340 jmc 1.150 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
341 heimbach 1.138 # ifndef DISABLE_RSTAR_CODE
342 jahn 1.164 # ifndef ALLOW_AUTODIFF_OPENAD
343 heimbach 1.138 dWtransC(i,j,bi,bj) = 0. _d 0
344     dWtransU(i,j,bi,bj) = 0. _d 0
345     dWtransV(i,j,bi,bj) = 0. _d 0
346 jahn 1.164 # endif
347 heimbach 1.138 # endif
348     # endif
349     #endif
350 heimbach 1.76 ENDDO
351     ENDDO
352 heimbach 1.49
353 jmc 1.63 C-- Start computation of dynamics
354 jmc 1.93 iMin = 0
355     iMax = sNx+1
356     jMin = 0
357     jMax = sNy+1
358 jmc 1.63
359 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
360 jmc 1.161 CADJ STORE wVel (:,:,:,bi,bj) =
361 heimbach 1.141 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
362 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
363    
364 jmc 1.161 C-- Explicit part of the Surface Potential Gradient (add in TIMESTEP)
365 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
366     IF (implicSurfPress.NE.1.) THEN
367 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
368     I bi,bj,iMin,iMax,jMin,jMax,
369     I etaN,
370     O phiSurfX,phiSurfY,
371 jmc 1.136 I myThid )
372 jmc 1.63 ENDIF
373 heimbach 1.83
374     #ifdef ALLOW_AUTODIFF_TAMC
375 jmc 1.161 CADJ STORE uVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
376     CADJ STORE vVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
377 heimbach 1.83 #ifdef ALLOW_KPP
378 jmc 1.150 CADJ STORE KPPviscAz (:,:,:,bi,bj)
379 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
380 heimbach 1.83 #endif /* ALLOW_KPP */
381     #endif /* ALLOW_AUTODIFF_TAMC */
382 adcroft 1.58
383 heimbach 1.77 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
384 jmc 1.140 C-- Calculate the total vertical viscosity
385     CALL CALC_VISCOSITY(
386     I bi,bj, iMin,iMax,jMin,jMax,
387     O KappaRU, KappaRV,
388     I myThid )
389     #else
390 heimbach 1.77 DO k=1,Nr
391 jmc 1.140 DO j=1-OLy,sNy+OLy
392     DO i=1-OLx,sNx+OLx
393     KappaRU(i,j,k) = 0. _d 0
394     KappaRV(i,j,k) = 0. _d 0
395     ENDDO
396     ENDDO
397     ENDDO
398 heimbach 1.77 #endif
399    
400 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
401 jmc 1.150 CADJ STORE KappaRU(:,:,:)
402 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
403 jmc 1.150 CADJ STORE KappaRV(:,:,:)
404 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
405 heimbach 1.101 #endif /* ALLOW_AUTODIFF_TAMC */
406    
407 mlosch 1.159 #ifdef ALLOW_OBCS
408     C-- For Stevens boundary conditions velocities need to be extrapolated
409     C (copied) to a narrow strip outside the domain
410     IF ( useOBCS ) THEN
411 jmc 1.160 CALL OBCS_COPY_UV_N(
412 jmc 1.161 U uVel(1-OLx,1-OLy,1,bi,bj),
413     U vVel(1-OLx,1-OLy,1,bi,bj),
414 mlosch 1.159 I Nr, bi, bj, myThid )
415     ENDIF
416     #endif /* ALLOW_OBCS */
417    
418 adcroft 1.58 C-- Start of dynamics loop
419     DO k=1,Nr
420    
421     C-- km1 Points to level above k (=k-1)
422     C-- kup Cycles through 1,2 to point to layer above
423     C-- kDown Cycles through 2,1 to point to current layer
424    
425     km1 = MAX(1,k-1)
426 heimbach 1.77 kp1 = MIN(k+1,Nr)
427 adcroft 1.58 kup = 1+MOD(k+1,2)
428     kDown= 1+MOD(k,2)
429    
430 jmc 1.144 #ifdef ALLOW_AUTODIFF_TAMC
431 heimbach 1.91 kkey = (idynkey-1)*Nr + k
432 heimbach 1.99 c
433 jmc 1.161 CADJ STORE totPhiHyd (:,:,k,bi,bj)
434 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
435 jmc 1.161 CADJ STORE phiHydLow (:,:,bi,bj)
436 gforget 1.143 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
437 jmc 1.150 CADJ STORE theta (:,:,k,bi,bj)
438 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
439 jmc 1.150 CADJ STORE salt (:,:,k,bi,bj)
440 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
441 jmc 1.161 CADJ STORE gT(:,:,k,bi,bj)
442 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
443 jmc 1.161 CADJ STORE gS(:,:,k,bi,bj)
444 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
445 heimbach 1.126 # ifdef NONLIN_FRSURF
446     cph-test
447 jmc 1.150 CADJ STORE phiHydC (:,:)
448 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
449 jmc 1.150 CADJ STORE phiHydF (:,:)
450 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
451 jmc 1.161 CADJ STORE guDissip (:,:)
452 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
453 jmc 1.161 CADJ STORE gvDissip (:,:)
454 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
455 jmc 1.150 CADJ STORE fVerU (:,:,:)
456 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
457 jmc 1.150 CADJ STORE fVerV (:,:,:)
458 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
459 jmc 1.161 CADJ STORE gU(:,:,k,bi,bj)
460 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
461 jmc 1.161 CADJ STORE gV(:,:,k,bi,bj)
462 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
463 heimbach 1.148 # ifndef ALLOW_ADAMSBASHFORTH_3
464 jmc 1.161 CADJ STORE guNm1(:,:,k,bi,bj)
465 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
466 jmc 1.161 CADJ STORE gvNm1(:,:,k,bi,bj)
467 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
468 heimbach 1.148 # else
469 jmc 1.161 CADJ STORE guNm(:,:,k,bi,bj,1)
470 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
471 jmc 1.161 CADJ STORE guNm(:,:,k,bi,bj,2)
472 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
473 jmc 1.161 CADJ STORE gvNm(:,:,k,bi,bj,1)
474 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
475 jmc 1.161 CADJ STORE gvNm(:,:,k,bi,bj,2)
476 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
477     # endif
478 heimbach 1.126 # ifdef ALLOW_CD_CODE
479 jmc 1.161 CADJ STORE uNM1(:,:,k,bi,bj)
480 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
481 jmc 1.161 CADJ STORE vNM1(:,:,k,bi,bj)
482 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
483 jmc 1.150 CADJ STORE uVelD(:,:,k,bi,bj)
484 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
485 jmc 1.150 CADJ STORE vVelD(:,:,k,bi,bj)
486 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
487     # endif
488     # endif
489 heimbach 1.134 # ifdef ALLOW_DEPTH_CONTROL
490 jmc 1.150 CADJ STORE fVerU (:,:,:)
491 heimbach 1.134 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
492 jmc 1.150 CADJ STORE fVerV (:,:,:)
493 heimbach 1.134 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
494     # endif
495 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
496    
497 jmc 1.144 C-- Integrate hydrostatic balance for phiHyd with BC of
498 adcroft 1.58 C phiHyd(z=0)=0
499 jmc 1.128 IF ( implicitIntGravWave ) THEN
500     CALL CALC_PHI_HYD(
501     I bi,bj,iMin,iMax,jMin,jMax,k,
502     I gT, gS,
503     U phiHydF,
504     O phiHydC, dPhiHydX, dPhiHydY,
505     I myTime, myIter, myThid )
506     ELSE
507     CALL CALC_PHI_HYD(
508 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
509     I theta, salt,
510 jmc 1.94 U phiHydF,
511     O phiHydC, dPhiHydX, dPhiHydY,
512 jmc 1.92 I myTime, myIter, myThid )
513 jmc 1.128 ENDIF
514 jmc 1.145 #ifdef ALLOW_DIAGNOSTICS
515     IF ( dPhiHydDiagIsOn ) THEN
516     tmpFac = -1. _d 0
517     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
518     & 'Um_dPHdx', k, 1, 2, bi, bj, myThid )
519     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
520     & 'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
521     ENDIF
522     #endif /* ALLOW_DIAGNOSTICS */
523 mlosch 1.89
524 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
525 jmc 1.96 C and step forward storing the result in gU, gV, etc...
526 adcroft 1.58 IF ( momStepping ) THEN
527 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
528 jmc 1.162 # ifdef NONLIN_FRSURF
529     # if (defined ALLOW_MOM_FLUXFORM) && !(defined DISABLE_RSTAR_CODE)
530 jmc 1.150 CADJ STORE dWtransC(:,:,bi,bj)
531 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
532 jmc 1.150 CADJ STORE dWtransU(:,:,bi,bj)
533 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
534 jmc 1.150 CADJ STORE dWtransV(:,:,bi,bj)
535 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
536     # endif
537 jmc 1.162 CADJ STORE fVerU(:,:,:)
538     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
539     CADJ STORE fVerV(:,:,:)
540     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
541     # endif /* NONLIN_FRSURF */
542     #endif /* ALLOW_AUTODIFF_TAMC */
543 heimbach 1.132 IF (.NOT. vectorInvariantMomentum) THEN
544 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
545 heimbach 1.132 CALL MOM_FLUXFORM(
546 jmc 1.162 I bi,bj,k,iMin,iMax,jMin,jMax,
547 jmc 1.121 I KappaRU, KappaRV,
548 jmc 1.162 U fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
549     O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
550 jmc 1.121 O guDissip, gvDissip,
551 adcroft 1.80 I myTime, myIter, myThid)
552 adcroft 1.79 #endif
553 heimbach 1.132 ELSE
554 edhill 1.105 #ifdef ALLOW_MOM_VECINV
555 heimbach 1.126 CALL MOM_VECINV(
556 jmc 1.162 I bi,bj,k,iMin,iMax,jMin,jMax,
557 jmc 1.121 I KappaRU, KappaRV,
558 jmc 1.162 I fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
559     O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
560 jmc 1.110 O guDissip, gvDissip,
561 adcroft 1.80 I myTime, myIter, myThid)
562 heimbach 1.132 #endif
563 heimbach 1.126 ENDIF
564 heimbach 1.132 C
565 adcroft 1.58 CALL TIMESTEP(
566 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
567 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
568 jmc 1.110 I guDissip, gvDissip,
569 jmc 1.96 I myTime, myIter, myThid)
570 adcroft 1.58
571     ENDIF
572    
573     C-- end of dynamics k loop (1:Nr)
574     ENDDO
575    
576 jmc 1.106 C-- Implicit Vertical advection & viscosity
577 gforget 1.147 #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
578     defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
579 jmc 1.106 IF ( momImplVertAdv ) THEN
580 jmc 1.136 CALL MOM_U_IMPLICIT_R( kappaRU,
581 jmc 1.106 I bi, bj, myTime, myIter, myThid )
582     CALL MOM_V_IMPLICIT_R( kappaRV,
583     I bi, bj, myTime, myIter, myThid )
584     ELSEIF ( implicitViscosity ) THEN
585     #else /* INCLUDE_IMPLVERTADV_CODE */
586     IF ( implicitViscosity ) THEN
587     #endif /* INCLUDE_IMPLVERTADV_CODE */
588 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
589 heimbach 1.101 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
590 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
591 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
592 adcroft 1.42 CALL IMPLDIFF(
593     I bi, bj, iMin, iMax, jMin, jMax,
594 jmc 1.160 I -1, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
595 jmc 1.96 U gU,
596 adcroft 1.42 I myThid )
597 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
598 heimbach 1.101 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
599 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
600 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
601 adcroft 1.42 CALL IMPLDIFF(
602     I bi, bj, iMin, iMax, jMin, jMax,
603 jmc 1.160 I -2, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
604 jmc 1.96 U gV,
605 adcroft 1.42 I myThid )
606 jmc 1.106 ENDIF
607 heimbach 1.49
608 mlosch 1.159 #ifdef ALLOW_OBCS
609 adcroft 1.58 C-- Apply open boundary conditions
610 jmc 1.151 IF ( useOBCS ) THEN
611 jmc 1.160 C-- but first save intermediate velocities to be used in the
612 mlosch 1.159 C next time step for the Stevens boundary conditions
613 jmc 1.160 CALL OBCS_SAVE_UV_N(
614     I bi, bj, iMin, iMax, jMin, jMax, 0,
615 mlosch 1.159 I gU, gV, myThid )
616 jmc 1.151 CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
617 jmc 1.106 ENDIF
618 mlosch 1.159 #endif /* ALLOW_OBCS */
619 heimbach 1.49
620 edhill 1.102 #ifdef ALLOW_CD_CODE
621 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
622 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
623 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
624 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
625 adcroft 1.42 CALL IMPLDIFF(
626     I bi, bj, iMin, iMax, jMin, jMax,
627 jmc 1.160 I 0, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
628 adcroft 1.42 U vVelD,
629     I myThid )
630 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
631 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
632 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
633 adcroft 1.42 CALL IMPLDIFF(
634     I bi, bj, iMin, iMax, jMin, jMax,
635 jmc 1.160 I 0, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
636 adcroft 1.42 U uVelD,
637     I myThid )
638 jmc 1.106 ENDIF
639 edhill 1.102 #endif /* ALLOW_CD_CODE */
640 jmc 1.106 C-- End implicit Vertical advection & viscosity
641 heimbach 1.109
642 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
643    
644 jmc 1.122 #ifdef ALLOW_NONHYDROSTATIC
645     C-- Step forward W field in N-H algorithm
646 jmc 1.136 IF ( nonHydrostatic ) THEN
647 jmc 1.122 #ifdef ALLOW_DEBUG
648 jmc 1.153 IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
649 jmc 1.122 #endif
650     CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
651 baylor 1.135 CALL CALC_GW(
652 jmc 1.136 I bi,bj, KappaRU, KappaRV,
653     I myTime, myIter, myThid )
654     ENDIF
655     IF ( nonHydrostatic.OR.implicitIntGravWave )
656     & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
657     IF ( nonHydrostatic )
658 jmc 1.128 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
659 jmc 1.122 #endif
660    
661     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
662    
663 jmc 1.136 C- end of bi,bj loops
664     ENDDO
665     ENDDO
666    
667     #ifdef ALLOW_OBCS
668 jmc 1.152 IF (useOBCS) THEN
669 jmc 1.155 CALL OBCS_EXCHANGES( myThid )
670 jmc 1.152 ENDIF
671 jmc 1.136 #endif
672    
673 mlosch 1.90 Cml(
674     C In order to compare the variance of phiHydLow of a p/z-coordinate
675     C run with etaH of a z/p-coordinate run the drift of phiHydLow
676     C has to be removed by something like the following subroutine:
677 jmc 1.144 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
678     C & 'phiHydLow', myTime, myThid )
679 mlosch 1.90 Cml)
680 adcroft 1.69
681 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
682 jmc 1.130 IF ( useDiagnostics ) THEN
683 jmc 1.113
684     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
685 jmc 1.120 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
686 molod 1.116
687 jmc 1.120 tmpFac = 1. _d 0
688     CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
689     & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
690 molod 1.116
691 jmc 1.120 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
692     & 'PHIBOTSQ',0, 1,0,1,1,myThid)
693 jmc 1.113
694     ENDIF
695     #endif /* ALLOW_DIAGNOSTICS */
696 jmc 1.136
697 edhill 1.104 #ifdef ALLOW_DEBUG
698 jmc 1.158 IF ( debugLevel .GE. debLevD ) THEN
699 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
700 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
701 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
702     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
703     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
704     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
705 jmc 1.115 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
706     CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
707     CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
708     CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
709     #ifndef ALLOW_ADAMSBASHFORTH_3
710     CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
711     CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
712     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
713     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
714     #endif
715 adcroft 1.70 ENDIF
716 adcroft 1.69 #endif
717 cnh 1.1
718 jmc 1.125 #ifdef DYNAMICS_GUGV_EXCH_CHECK
719 jmc 1.144 C- jmc: For safety checking only: This Exchange here should not change
720     C the solution. If solution changes, it means something is wrong,
721 jmc 1.125 C but it does not mean that it is less wrong with this exchange.
722 jmc 1.158 IF ( debugLevel .GE. debLevE ) THEN
723 jmc 1.125 CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
724     ENDIF
725     #endif
726    
727 jmc 1.123 #ifdef ALLOW_DEBUG
728 jmc 1.153 IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
729 jmc 1.123 #endif
730    
731 cnh 1.1 RETURN
732     END

  ViewVC Help
Powered by ViewVC 1.1.22