/[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.163 - (hide annotations) (download)
Thu Nov 15 15:55:42 2012 UTC (11 years, 6 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d
Changes since 1.162: +2 -1 lines
adding tidal velocity forcing capability to obcs
 Modified Files:
  model/src/dynamics.F forward_step.F the_main_loop.F
  pkg/obcs/OBCS_FIELDS.h OBCS_OPTIONS.h OBCS_PARAMS.h
  OBCS_SEAICE.h obcs_apply_eta.F obcs_apply_r_star.F
  obcs_apply_surf_dr.F obcs_apply_ts.F obcs_apply_w.F
  obcs_calc.F obcs_check.F obcs_init_variables.F obcs_readparms.F
  verification/seaice_obcs/code/OBCS_OPTIONS.h
 Added Files:
  pkg/obcs/obcs_add_tides.F
  verification/seaice_obcs/input.tides/*
  verification/seaice_obcs/results/output.tides.txt

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

  ViewVC Help
Powered by ViewVC 1.1.22