/[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.170 - (hide annotations) (download)
Fri Apr 4 20:54:11 2014 UTC (10 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65a, checkpoint65, checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64w, checkpoint64v
Changes since 1.169: +7 -4 lines
- Start to include explicitly AUTODIFF_OPTIONS.h, COST_OPTIONS.h,
  and CTRL_OPTIONS.h in src files (to enable to skip the ECCO_CPPOPTIONS.h)
  For now, only in pkgs used in verification/hs94.1x64x5.
- Replace ALLOW_AUTODIFF_TAMC by ALLOW_AUTODIFF (except for tape/storage
  which are specific to TAF/TAMC).

1 jmc 1.170 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.169 2014/02/06 23:16:46 gforget 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 jmc 1.170 #ifdef ALLOW_AUTODIFF
7     # include "AUTODIFF_OPTIONS.h"
8     #endif
9 jmc 1.167 #ifdef ALLOW_MOM_COMMON
10     # include "MOM_COMMON_OPTIONS.h"
11     #endif
12 heimbach 1.131 #ifdef ALLOW_OBCS
13     # include "OBCS_OPTIONS.h"
14     #endif
15    
16 jmc 1.125 #undef DYNAMICS_GUGV_EXCH_CHECK
17 cnh 1.1
18 cnh 1.82 CBOP
19     C !ROUTINE: DYNAMICS
20     C !INTERFACE:
21 cnh 1.8 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
22 cnh 1.82 C !DESCRIPTION: \bv
23     C *==========================================================*
24 jmc 1.144 C | SUBROUTINE DYNAMICS
25     C | o Controlling routine for the explicit part of the model
26     C | dynamics.
27 cnh 1.82 C *==========================================================*
28 jmc 1.144 C | This routine evaluates the "dynamics" terms for each
29     C | block of ocean in turn. Because the blocks of ocean have
30     C | overlap regions they are independent of one another.
31     C | If terms involving lateral integrals are needed in this
32     C | routine care will be needed. Similarly finite-difference
33     C | operations with stencils wider than the overlap region
34     C | require special consideration.
35 cnh 1.82 C | The algorithm...
36     C |
37     C | "Correction Step"
38     C | =================
39     C | Here we update the horizontal velocities with the surface
40     C | pressure such that the resulting flow is either consistent
41     C | with the free-surface evolution or the rigid-lid:
42     C | U[n] = U* + dt x d/dx P
43     C | V[n] = V* + dt x d/dy P
44 jmc 1.122 C | W[n] = W* + dt x d/dz P (NH mode)
45 cnh 1.82 C |
46     C | "Calculation of Gs"
47     C | ===================
48     C | This is where all the accelerations and tendencies (ie.
49     C | physics, parameterizations etc...) are calculated
50     C | rho = rho ( theta[n], salt[n] )
51     C | b = b(rho, theta)
52     C | K31 = K31 ( rho )
53     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
54     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
55     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
56     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
57     C |
58     C | "Time-stepping" or "Prediction"
59     C | ================================
60     C | The models variables are stepped forward with the appropriate
61     C | time-stepping scheme (currently we use Adams-Bashforth II)
62     C | - For momentum, the result is always *only* a "prediction"
63     C | in that the flow may be divergent and will be "corrected"
64     C | later with a surface pressure gradient.
65     C | - Normally for tracers the result is the new field at time
66     C | level [n+1} *BUT* in the case of implicit diffusion the result
67     C | is also *only* a prediction.
68     C | - We denote "predictors" with an asterisk (*).
69     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
70     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
71     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
72     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
73     C | With implicit diffusion:
74     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
75     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
76     C | (1 + dt * K * d_zz) theta[n] = theta*
77     C | (1 + dt * K * d_zz) salt[n] = salt*
78     C |
79     C *==========================================================*
80     C \ev
81     C !USES:
82 adcroft 1.40 IMPLICIT NONE
83 cnh 1.1 C == Global variables ===
84     #include "SIZE.h"
85     #include "EEPARAMS.h"
86 adcroft 1.6 #include "PARAMS.h"
87 jmc 1.165 #include "GRID.h"
88 adcroft 1.3 #include "DYNVARS.h"
89 jmc 1.167 #ifdef ALLOW_MOM_COMMON
90     # include "MOM_VISC.h"
91     #endif
92 edhill 1.103 #ifdef ALLOW_CD_CODE
93 jmc 1.165 # include "CD_CODE_VARS.h"
94 edhill 1.103 #endif
95 jmc 1.170 #ifdef ALLOW_AUTODIFF
96 heimbach 1.53 # include "tamc.h"
97     # include "tamc_keys.h"
98 heimbach 1.67 # include "FFIELDS.h"
99 heimbach 1.91 # include "EOS.h"
100 heimbach 1.67 # ifdef ALLOW_KPP
101     # include "KPP.h"
102     # endif
103 heimbach 1.131 # ifdef ALLOW_PTRACERS
104     # include "PTRACERS_SIZE.h"
105 jmc 1.139 # include "PTRACERS_FIELDS.h"
106 heimbach 1.131 # endif
107     # ifdef ALLOW_OBCS
108 jmc 1.165 # include "OBCS_PARAMS.h"
109 jmc 1.157 # include "OBCS_FIELDS.h"
110 heimbach 1.131 # ifdef ALLOW_PTRACERS
111     # include "OBCS_PTRACERS.h"
112     # endif
113     # endif
114 heimbach 1.133 # ifdef ALLOW_MOM_FLUXFORM
115     # include "MOM_FLUXFORM.h"
116     # endif
117 jmc 1.170 #endif /* ALLOW_AUTODIFF */
118 jmc 1.62
119 cnh 1.82 C !CALLING SEQUENCE:
120     C DYNAMICS()
121     C |
122 jmc 1.122 C |-- CALC_EP_FORCING
123     C |
124 cnh 1.82 C |-- CALC_GRAD_PHI_SURF
125     C |
126     C |-- CALC_VISCOSITY
127     C |
128 jmc 1.167 C |-- MOM_CALC_3D_STRAIN
129     C |
130 m_bates 1.166 C |-- CALC_EDDY_STRESS
131     C |
132 jmc 1.136 C |-- CALC_PHI_HYD
133 cnh 1.82 C |
134 jmc 1.136 C |-- MOM_FLUXFORM
135 cnh 1.82 C |
136 jmc 1.136 C |-- MOM_VECINV
137 cnh 1.82 C |
138 jmc 1.167 C |-- MOM_CALC_SMAG_3D
139     C |-- MOM_UV_SMAG_3D
140     C |
141 jmc 1.136 C |-- TIMESTEP
142 cnh 1.82 C |
143 jmc 1.136 C |-- MOM_U_IMPLICIT_R
144     C |-- MOM_V_IMPLICIT_R
145 jmc 1.122 C |
146 jmc 1.136 C |-- IMPLDIFF
147 cnh 1.82 C |
148     C |-- OBCS_APPLY_UV
149     C |
150 jmc 1.122 C |-- CALC_GW
151     C |
152     C |-- DIAGNOSTICS_FILL
153     C |-- DEBUG_STATS_RL
154 cnh 1.82
155     C !INPUT/OUTPUT PARAMETERS:
156 cnh 1.1 C == Routine arguments ==
157 jmc 1.140 C myTime :: Current time in simulation
158     C myIter :: Current iteration number in simulation
159     C myThid :: Thread number for this instance of the routine.
160 cnh 1.8 _RL myTime
161     INTEGER myIter
162 adcroft 1.47 INTEGER myThid
163 cnh 1.1
164 jmc 1.145 C !FUNCTIONS:
165     #ifdef ALLOW_DIAGNOSTICS
166     LOGICAL DIAGNOSTICS_IS_ON
167     EXTERNAL DIAGNOSTICS_IS_ON
168     #endif
169    
170 cnh 1.82 C !LOCAL VARIABLES:
171 cnh 1.1 C == Local variables
172 jmc 1.113 C fVer[UV] o fVer: Vertical flux term - note fVer
173     C is "pipelined" in the vertical
174     C so we need an fVer for each
175     C variable.
176 jmc 1.94 C phiHydC :: hydrostatic potential anomaly at cell center
177     C In z coords phiHyd is the hydrostatic potential
178     C (=pressure/rho0) anomaly
179     C In p coords phiHyd is the geopotential height anomaly.
180     C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
181     C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
182     C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
183 jmc 1.92 C phiSurfY or geopotential (atmos) in X and Y direction
184 jmc 1.110 C guDissip :: dissipation tendency (all explicit terms), u component
185     C gvDissip :: dissipation tendency (all explicit terms), v component
186 jmc 1.165 C KappaRU :: vertical viscosity for velocity U-component
187     C KappaRV :: vertical viscosity for velocity V-component
188 jmc 1.162 C iMin, iMax :: Ranges and sub-block indices on which calculations
189     C jMin, jMax are applied.
190     C bi, bj :: tile indices
191     C k :: current level index
192     C km1, kp1 :: index of level above (k-1) and below (k+1)
193     C kUp, kDown :: Index for interface above and below. kUp and kDown are
194     C are switched with k to be the appropriate index into fVerU,V
195 cnh 1.30 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
196     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
197 jmc 1.94 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
198     _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
199 jmc 1.161 _RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
200     _RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
201 jmc 1.63 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
202     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
203 jmc 1.110 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
204     _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
205 jmc 1.161 _RL KappaRU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
206     _RL KappaRV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
207 jmc 1.167 #ifdef ALLOW_SMAG_3D
208     C str11 :: strain component Vxx @ grid-cell center
209     C str22 :: strain component Vyy @ grid-cell center
210     C str33 :: strain component Vzz @ grid-cell center
211     C str12 :: strain component Vxy @ grid-cell corner
212     C str13 :: strain component Vxz @ above uVel
213     C str23 :: strain component Vyz @ above vVel
214     C viscAh3d_00 :: Smagorinsky viscosity @ grid-cell center
215     C viscAh3d_12 :: Smagorinsky viscosity @ grid-cell corner
216     C viscAh3d_13 :: Smagorinsky viscosity @ above uVel
217     C viscAh3d_23 :: Smagorinsky viscosity @ above vVel
218     C addDissU :: zonal momentum tendency from 3-D Smag. viscosity
219     C addDissV :: merid momentum tendency from 3-D Smag. viscosity
220     _RL str11(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
221     _RL str22(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
222     _RL str33(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
223     _RL str12(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
224     _RL str13(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
225     _RL str23(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
226     _RL viscAh3d_00(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
227     _RL viscAh3d_12(1-OLx:sNx+OLx,1-OLy:sNy+OLy, Nr )
228     _RL viscAh3d_13(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
229     _RL viscAh3d_23(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
230     _RL addDissU(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
231     _RL addDissV(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
232     #elif ( defined ALLOW_NONHYDROSTATIC )
233     _RL str13(1), str23(1), str33(1)
234     _RL viscAh3d_00(1), viscAh3d_13(1), viscAh3d_23(1)
235     #endif
236 adcroft 1.12
237 cnh 1.1 INTEGER bi, bj
238     INTEGER i, j
239 jmc 1.162 INTEGER k, km1, kp1, kUp, kDown
240 jmc 1.167 INTEGER iMin, iMax
241     INTEGER jMin, jMax
242     PARAMETER( iMin = 0 , iMax = sNx+1 )
243     PARAMETER( jMin = 0 , jMax = sNy+1 )
244 cnh 1.1
245 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
246 jmc 1.145 LOGICAL dPhiHydDiagIsOn
247 jmc 1.120 _RL tmpFac
248 jmc 1.113 #endif /* ALLOW_DIAGNOSTICS */
249    
250 adcroft 1.11 C--- The algorithm...
251     C
252     C "Correction Step"
253     C =================
254     C Here we update the horizontal velocities with the surface
255     C pressure such that the resulting flow is either consistent
256     C with the free-surface evolution or the rigid-lid:
257     C U[n] = U* + dt x d/dx P
258     C V[n] = V* + dt x d/dy P
259     C
260     C "Calculation of Gs"
261     C ===================
262     C This is where all the accelerations and tendencies (ie.
263 heimbach 1.53 C physics, parameterizations etc...) are calculated
264 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
265 cnh 1.27 C b = b(rho, theta)
266 adcroft 1.11 C K31 = K31 ( rho )
267 jmc 1.61 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
268     C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
269     C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
270     C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
271 adcroft 1.11 C
272 adcroft 1.12 C "Time-stepping" or "Prediction"
273 adcroft 1.11 C ================================
274     C The models variables are stepped forward with the appropriate
275     C time-stepping scheme (currently we use Adams-Bashforth II)
276     C - For momentum, the result is always *only* a "prediction"
277     C in that the flow may be divergent and will be "corrected"
278     C later with a surface pressure gradient.
279     C - Normally for tracers the result is the new field at time
280     C level [n+1} *BUT* in the case of implicit diffusion the result
281     C is also *only* a prediction.
282     C - We denote "predictors" with an asterisk (*).
283     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
284     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
285     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
286     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
287 adcroft 1.12 C With implicit diffusion:
288 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
289     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
290 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
291     C (1 + dt * K * d_zz) salt[n] = salt*
292 adcroft 1.11 C---
293 cnh 1.82 CEOP
294 adcroft 1.11
295 jmc 1.123 #ifdef ALLOW_DEBUG
296 jmc 1.153 IF (debugMode) CALL DEBUG_ENTER( 'DYNAMICS', myThid )
297 jmc 1.123 #endif
298    
299 jmc 1.145 #ifdef ALLOW_DIAGNOSTICS
300     dPhiHydDiagIsOn = .FALSE.
301     IF ( useDiagnostics )
302     & dPhiHydDiagIsOn = DIAGNOSTICS_IS_ON( 'Um_dPHdx', myThid )
303     & .OR. DIAGNOSTICS_IS_ON( 'Vm_dPHdy', myThid )
304     #endif
305    
306 jmc 1.165 C-- Call to routine for calculation of Eliassen-Palm-flux-forced
307     C U-tendency, if desired:
308 heimbach 1.88 #ifdef INCLUDE_EP_FORCING_CODE
309     CALL CALC_EP_FORCING(myThid)
310     #endif
311    
312 heimbach 1.154 #ifdef ALLOW_AUTODIFF_MONITOR_DIAG
313 jmc 1.161 CALL DUMMY_IN_DYNAMICS( myTime, myIter, myThid )
314 heimbach 1.154 #endif
315    
316 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
317     C-- HPF directive to help TAMC
318     CHPF$ INDEPENDENT
319     #endif /* ALLOW_AUTODIFF_TAMC */
320    
321 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
322 heimbach 1.76
323     #ifdef ALLOW_AUTODIFF_TAMC
324     C-- HPF directive to help TAMC
325     CHPF$ INDEPENDENT, NEW (fVerU,fVerV
326 jmc 1.94 CHPF$& ,phiHydF
327 heimbach 1.76 CHPF$& ,KappaRU,KappaRV
328     CHPF$& )
329     #endif /* ALLOW_AUTODIFF_TAMC */
330    
331 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
332 heimbach 1.76
333     #ifdef ALLOW_AUTODIFF_TAMC
334     act1 = bi - myBxLo(myThid)
335     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
336     act2 = bj - myByLo(myThid)
337     max2 = myByHi(myThid) - myByLo(myThid) + 1
338     act3 = myThid - 1
339     max3 = nTx*nTy
340     act4 = ikey_dynamics - 1
341 heimbach 1.91 idynkey = (act1 + 1) + act2*max1
342 heimbach 1.76 & + act3*max1*max2
343     & + act4*max1*max2*max3
344     #endif /* ALLOW_AUTODIFF_TAMC */
345    
346 heimbach 1.97 C-- Set up work arrays with valid (i.e. not NaN) values
347 jmc 1.161 C These initial values do not alter the numerical results. They
348 heimbach 1.97 C just ensure that all memory references are to valid floating
349     C point numbers. This prevents spurious hardware signals due to
350     C uninitialised but inert locations.
351    
352 jmc 1.168 #ifdef ALLOW_AUTODIFF
353 jmc 1.94 DO k=1,Nr
354     DO j=1-OLy,sNy+OLy
355     DO i=1-OLx,sNx+OLx
356 heimbach 1.97 c-- need some re-initialisation here to break dependencies
357 jmc 1.122 gU(i,j,k,bi,bj) = 0. _d 0
358     gV(i,j,k,bi,bj) = 0. _d 0
359 heimbach 1.87 ENDDO
360 jmc 1.94 ENDDO
361     ENDDO
362 jmc 1.168 #endif /* ALLOW_AUTODIFF */
363 jmc 1.94 DO j=1-OLy,sNy+OLy
364     DO i=1-OLx,sNx+OLx
365 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
366     fVerU (i,j,2) = 0. _d 0
367     fVerV (i,j,1) = 0. _d 0
368     fVerV (i,j,2) = 0. _d 0
369 jmc 1.136 phiHydF (i,j) = 0. _d 0
370     phiHydC (i,j) = 0. _d 0
371 jmc 1.146 #ifndef INCLUDE_PHIHYD_CALCULATION_CODE
372 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
373 jmc 1.136 dPhiHydY(i,j) = 0. _d 0
374 jmc 1.146 #endif
375 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
376     phiSurfY(i,j) = 0. _d 0
377 jmc 1.110 guDissip(i,j) = 0. _d 0
378     gvDissip(i,j) = 0. _d 0
379 jmc 1.168 #ifdef ALLOW_AUTODIFF
380 gforget 1.143 phiHydLow(i,j,bi,bj) = 0. _d 0
381 jmc 1.150 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
382 heimbach 1.138 # ifndef DISABLE_RSTAR_CODE
383     dWtransC(i,j,bi,bj) = 0. _d 0
384     dWtransU(i,j,bi,bj) = 0. _d 0
385     dWtransV(i,j,bi,bj) = 0. _d 0
386     # endif
387     # endif
388 jmc 1.168 #endif /* ALLOW_AUTODIFF */
389 heimbach 1.76 ENDDO
390     ENDDO
391 heimbach 1.49
392 jmc 1.63 C-- Start computation of dynamics
393    
394 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
395 jmc 1.161 CADJ STORE wVel (:,:,:,bi,bj) =
396 heimbach 1.141 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
397 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
398    
399 jmc 1.161 C-- Explicit part of the Surface Potential Gradient (add in TIMESTEP)
400 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
401     IF (implicSurfPress.NE.1.) THEN
402 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
403     I bi,bj,iMin,iMax,jMin,jMax,
404     I etaN,
405     O phiSurfX,phiSurfY,
406 jmc 1.136 I myThid )
407 jmc 1.63 ENDIF
408 heimbach 1.83
409     #ifdef ALLOW_AUTODIFF_TAMC
410 jmc 1.161 CADJ STORE uVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
411     CADJ STORE vVel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
412 heimbach 1.83 #ifdef ALLOW_KPP
413 jmc 1.150 CADJ STORE KPPviscAz (:,:,:,bi,bj)
414 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
415 heimbach 1.83 #endif /* ALLOW_KPP */
416     #endif /* ALLOW_AUTODIFF_TAMC */
417 adcroft 1.58
418 jmc 1.168 #ifndef ALLOW_AUTODIFF
419 jmc 1.165 IF ( .NOT.momViscosity ) THEN
420 jmc 1.168 #endif
421 jmc 1.165 DO k=1,Nr
422     DO j=1-OLy,sNy+OLy
423     DO i=1-OLx,sNx+OLx
424     KappaRU(i,j,k) = 0. _d 0
425     KappaRV(i,j,k) = 0. _d 0
426     ENDDO
427     ENDDO
428     ENDDO
429 jmc 1.168 #ifndef ALLOW_AUTODIFF
430     ENDIF
431     #endif
432 jmc 1.165 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
433 jmc 1.140 C-- Calculate the total vertical viscosity
434 jmc 1.165 IF ( momViscosity ) THEN
435     CALL CALC_VISCOSITY(
436 jmc 1.140 I bi,bj, iMin,iMax,jMin,jMax,
437     O KappaRU, KappaRV,
438     I myThid )
439 jmc 1.165 ENDIF
440     #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */
441 heimbach 1.77
442 jmc 1.167 #ifdef ALLOW_SMAG_3D
443     IF ( useSmag3D ) THEN
444     CALL MOM_CALC_3D_STRAIN(
445     O str11, str22, str33, str12, str13, str23,
446     I bi, bj, myThid )
447     ENDIF
448     #endif /* ALLOW_SMAG_3D */
449    
450 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
451 jmc 1.150 CADJ STORE KappaRU(:,:,:)
452 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
453 jmc 1.150 CADJ STORE KappaRV(:,:,:)
454 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
455 heimbach 1.101 #endif /* ALLOW_AUTODIFF_TAMC */
456    
457 mlosch 1.159 #ifdef ALLOW_OBCS
458     C-- For Stevens boundary conditions velocities need to be extrapolated
459     C (copied) to a narrow strip outside the domain
460 jmc 1.165 IF ( useOBCS ) THEN
461 jmc 1.160 CALL OBCS_COPY_UV_N(
462 jmc 1.161 U uVel(1-OLx,1-OLy,1,bi,bj),
463     U vVel(1-OLx,1-OLy,1,bi,bj),
464 mlosch 1.159 I Nr, bi, bj, myThid )
465 jmc 1.165 ENDIF
466 mlosch 1.159 #endif /* ALLOW_OBCS */
467    
468 m_bates 1.166 #ifdef ALLOW_EDDYPSI
469     CALL CALC_EDDY_STRESS(bi,bj,myThid)
470     #endif
471    
472 adcroft 1.58 C-- Start of dynamics loop
473     DO k=1,Nr
474    
475     C-- km1 Points to level above k (=k-1)
476     C-- kup Cycles through 1,2 to point to layer above
477     C-- kDown Cycles through 2,1 to point to current layer
478    
479     km1 = MAX(1,k-1)
480 heimbach 1.77 kp1 = MIN(k+1,Nr)
481 adcroft 1.58 kup = 1+MOD(k+1,2)
482     kDown= 1+MOD(k,2)
483    
484 jmc 1.144 #ifdef ALLOW_AUTODIFF_TAMC
485 heimbach 1.91 kkey = (idynkey-1)*Nr + k
486 jmc 1.161 CADJ STORE totPhiHyd (:,:,k,bi,bj)
487 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
488 gforget 1.169 CADJ STORE phiHydLow (:,:,bi,bj)
489 jmc 1.170 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
490 jmc 1.150 CADJ STORE theta (:,:,k,bi,bj)
491 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
492 jmc 1.150 CADJ STORE salt (:,:,k,bi,bj)
493 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
494 jmc 1.161 CADJ STORE gT(:,:,k,bi,bj)
495 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
496 jmc 1.161 CADJ STORE gS(:,:,k,bi,bj)
497 heimbach 1.129 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
498 heimbach 1.126 # ifdef NONLIN_FRSURF
499     cph-test
500 jmc 1.150 CADJ STORE phiHydC (:,:)
501 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
502 jmc 1.150 CADJ STORE phiHydF (:,:)
503 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
504 jmc 1.161 CADJ STORE gU(:,:,k,bi,bj)
505 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
506 jmc 1.161 CADJ STORE gV(:,:,k,bi,bj)
507 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
508 heimbach 1.148 # ifndef ALLOW_ADAMSBASHFORTH_3
509 jmc 1.161 CADJ STORE guNm1(:,:,k,bi,bj)
510 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
511 jmc 1.161 CADJ STORE gvNm1(:,:,k,bi,bj)
512 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
513 heimbach 1.148 # else
514 jmc 1.161 CADJ STORE guNm(:,:,k,bi,bj,1)
515 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
516 jmc 1.161 CADJ STORE guNm(:,:,k,bi,bj,2)
517 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
518 jmc 1.161 CADJ STORE gvNm(:,:,k,bi,bj,1)
519 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
520 jmc 1.161 CADJ STORE gvNm(:,:,k,bi,bj,2)
521 heimbach 1.148 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
522     # endif
523 heimbach 1.126 # ifdef ALLOW_CD_CODE
524 jmc 1.161 CADJ STORE uNM1(:,:,k,bi,bj)
525 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
526 jmc 1.161 CADJ STORE vNM1(:,:,k,bi,bj)
527 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
528 jmc 1.150 CADJ STORE uVelD(:,:,k,bi,bj)
529 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
530 jmc 1.150 CADJ STORE vVelD(:,:,k,bi,bj)
531 heimbach 1.126 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
532     # endif
533 jmc 1.168 # endif /* NONLIN_FRSURF */
534 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
535    
536 jmc 1.165 C-- Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0
537 jmc 1.128 IF ( implicitIntGravWave ) THEN
538     CALL CALC_PHI_HYD(
539     I bi,bj,iMin,iMax,jMin,jMax,k,
540     I gT, gS,
541     U phiHydF,
542     O phiHydC, dPhiHydX, dPhiHydY,
543     I myTime, myIter, myThid )
544     ELSE
545     CALL CALC_PHI_HYD(
546 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
547     I theta, salt,
548 jmc 1.94 U phiHydF,
549     O phiHydC, dPhiHydX, dPhiHydY,
550 jmc 1.92 I myTime, myIter, myThid )
551 jmc 1.128 ENDIF
552 jmc 1.145 #ifdef ALLOW_DIAGNOSTICS
553     IF ( dPhiHydDiagIsOn ) THEN
554     tmpFac = -1. _d 0
555     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydX, tmpFac, 1,
556     & 'Um_dPHdx', k, 1, 2, bi, bj, myThid )
557     CALL DIAGNOSTICS_SCALE_FILL( dPhiHydY, tmpFac, 1,
558     & 'Vm_dPHdy', k, 1, 2, bi, bj, myThid )
559     ENDIF
560     #endif /* ALLOW_DIAGNOSTICS */
561 mlosch 1.89
562 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
563 jmc 1.96 C and step forward storing the result in gU, gV, etc...
564 adcroft 1.58 IF ( momStepping ) THEN
565 jmc 1.168 #ifdef ALLOW_AUTODIFF
566     DO j=1-OLy,sNy+OLy
567     DO i=1-OLx,sNx+OLx
568     guDissip(i,j) = 0. _d 0
569     gvDissip(i,j) = 0. _d 0
570     ENDDO
571     ENDDO
572     #endif /* ALLOW_AUTODIFF */
573 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
574 jmc 1.168 # if (defined NONLIN_FRSURF) && (defined ALLOW_MOM_FLUXFORM)
575     # ifndef DISABLE_RSTAR_CODE
576 jmc 1.150 CADJ STORE dWtransC(:,:,bi,bj)
577 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
578 jmc 1.150 CADJ STORE dWtransU(:,:,bi,bj)
579 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
580 jmc 1.150 CADJ STORE dWtransV(:,:,bi,bj)
581 heimbach 1.138 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
582     # endif
583 jmc 1.168 # endif /* NONLIN_FRSURF and ALLOW_MOM_FLUXFORM */
584     # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
585     CADJ STORE fVerU(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
586     CADJ STORE fVerV(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
587     # endif
588 jmc 1.162 #endif /* ALLOW_AUTODIFF_TAMC */
589 heimbach 1.132 IF (.NOT. vectorInvariantMomentum) THEN
590 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
591 heimbach 1.132 CALL MOM_FLUXFORM(
592 jmc 1.162 I bi,bj,k,iMin,iMax,jMin,jMax,
593 jmc 1.121 I KappaRU, KappaRV,
594 jmc 1.162 U fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
595     O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
596 jmc 1.121 O guDissip, gvDissip,
597 adcroft 1.80 I myTime, myIter, myThid)
598 adcroft 1.79 #endif
599 heimbach 1.132 ELSE
600 edhill 1.105 #ifdef ALLOW_MOM_VECINV
601 heimbach 1.126 CALL MOM_VECINV(
602 jmc 1.162 I bi,bj,k,iMin,iMax,jMin,jMax,
603 jmc 1.121 I KappaRU, KappaRV,
604 jmc 1.162 I fVerU(1-OLx,1-OLy,kUp), fVerV(1-OLx,1-OLy,kUp),
605     O fVerU(1-OLx,1-OLy,kDown), fVerV(1-OLx,1-OLy,kDown),
606 jmc 1.110 O guDissip, gvDissip,
607 adcroft 1.80 I myTime, myIter, myThid)
608 heimbach 1.132 #endif
609 heimbach 1.126 ENDIF
610 jmc 1.165
611 jmc 1.167 #ifdef ALLOW_SMAG_3D
612     IF ( useSmag3D ) THEN
613     CALL MOM_CALC_SMAG_3D(
614     I str11, str22, str33, str12, str13, str23,
615     O viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
616     I smag3D_hLsC, smag3D_hLsW, smag3D_hLsS, smag3D_hLsZ,
617     I k, bi, bj, myThid )
618     CALL MOM_UV_SMAG_3D(
619     I str11, str22, str12, str13, str23,
620     I viscAh3d_00, viscAh3d_12, viscAh3d_13, viscAh3d_23,
621     O addDissU, addDissV,
622     I iMin,iMax,jMin,jMax, k, bi, bj, myThid )
623     DO j= jMin,jMax
624     DO i= iMin,iMax
625     guDissip(i,j) = guDissip(i,j) + addDissU(i,j)
626     gvDissip(i,j) = gvDissip(i,j) + addDissV(i,j)
627     ENDDO
628     ENDDO
629     ENDIF
630     #endif /* ALLOW_SMAG_3D */
631    
632 adcroft 1.58 CALL TIMESTEP(
633 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
634 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
635 jmc 1.110 I guDissip, gvDissip,
636 jmc 1.96 I myTime, myIter, myThid)
637 adcroft 1.58
638     ENDIF
639    
640     C-- end of dynamics k loop (1:Nr)
641     ENDDO
642    
643 jmc 1.106 C-- Implicit Vertical advection & viscosity
644 gforget 1.147 #if (defined (INCLUDE_IMPLVERTADV_CODE) && \
645     defined (ALLOW_MOM_COMMON) && !(defined ALLOW_AUTODIFF_TAMC))
646 jmc 1.106 IF ( momImplVertAdv ) THEN
647 jmc 1.136 CALL MOM_U_IMPLICIT_R( kappaRU,
648 jmc 1.106 I bi, bj, myTime, myIter, myThid )
649     CALL MOM_V_IMPLICIT_R( kappaRV,
650     I bi, bj, myTime, myIter, myThid )
651     ELSEIF ( implicitViscosity ) THEN
652     #else /* INCLUDE_IMPLVERTADV_CODE */
653     IF ( implicitViscosity ) THEN
654     #endif /* INCLUDE_IMPLVERTADV_CODE */
655 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
656 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
657 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
658 adcroft 1.42 CALL IMPLDIFF(
659     I bi, bj, iMin, iMax, jMin, jMax,
660 jmc 1.160 I -1, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
661 jmc 1.96 U gU,
662 adcroft 1.42 I myThid )
663 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
664 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
665 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
666 adcroft 1.42 CALL IMPLDIFF(
667     I bi, bj, iMin, iMax, jMin, jMax,
668 jmc 1.160 I -2, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
669 jmc 1.96 U gV,
670 adcroft 1.42 I myThid )
671 jmc 1.106 ENDIF
672 heimbach 1.49
673 mlosch 1.159 #ifdef ALLOW_OBCS
674 adcroft 1.58 C-- Apply open boundary conditions
675 jmc 1.151 IF ( useOBCS ) THEN
676 jmc 1.160 C-- but first save intermediate velocities to be used in the
677 mlosch 1.159 C next time step for the Stevens boundary conditions
678 jmc 1.160 CALL OBCS_SAVE_UV_N(
679     I bi, bj, iMin, iMax, jMin, jMax, 0,
680 mlosch 1.159 I gU, gV, myThid )
681 jmc 1.151 CALL OBCS_APPLY_UV( bi, bj, 0, gU, gV, myThid )
682 jmc 1.106 ENDIF
683 mlosch 1.159 #endif /* ALLOW_OBCS */
684 heimbach 1.49
685 edhill 1.102 #ifdef ALLOW_CD_CODE
686 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
687 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
688 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
689 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
690 adcroft 1.42 CALL IMPLDIFF(
691     I bi, bj, iMin, iMax, jMin, jMax,
692 jmc 1.160 I 0, KappaRU, recip_hFacW(1-OLx,1-OLy,1,bi,bj),
693 adcroft 1.42 U vVelD,
694     I myThid )
695 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
696 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
697 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
698 adcroft 1.42 CALL IMPLDIFF(
699     I bi, bj, iMin, iMax, jMin, jMax,
700 jmc 1.160 I 0, KappaRV, recip_hFacS(1-OLx,1-OLy,1,bi,bj),
701 adcroft 1.42 U uVelD,
702     I myThid )
703 jmc 1.106 ENDIF
704 edhill 1.102 #endif /* ALLOW_CD_CODE */
705 jmc 1.106 C-- End implicit Vertical advection & viscosity
706 heimbach 1.109
707 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
708    
709 jmc 1.122 #ifdef ALLOW_NONHYDROSTATIC
710     C-- Step forward W field in N-H algorithm
711 jmc 1.136 IF ( nonHydrostatic ) THEN
712 jmc 1.122 #ifdef ALLOW_DEBUG
713 jmc 1.153 IF (debugMode) CALL DEBUG_CALL('CALC_GW', myThid )
714 jmc 1.122 #endif
715     CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
716 baylor 1.135 CALL CALC_GW(
717 jmc 1.136 I bi,bj, KappaRU, KappaRV,
718 jmc 1.167 I str13, str23, str33,
719     I viscAh3d_00, viscAh3d_13, viscAh3d_23,
720 jmc 1.136 I myTime, myIter, myThid )
721     ENDIF
722     IF ( nonHydrostatic.OR.implicitIntGravWave )
723     & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
724     IF ( nonHydrostatic )
725 jmc 1.128 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
726 jmc 1.122 #endif
727    
728     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
729    
730 jmc 1.136 C- end of bi,bj loops
731     ENDDO
732     ENDDO
733    
734     #ifdef ALLOW_OBCS
735 jmc 1.152 IF (useOBCS) THEN
736 jmc 1.155 CALL OBCS_EXCHANGES( myThid )
737 jmc 1.152 ENDIF
738 jmc 1.136 #endif
739    
740 mlosch 1.90 Cml(
741     C In order to compare the variance of phiHydLow of a p/z-coordinate
742     C run with etaH of a z/p-coordinate run the drift of phiHydLow
743     C has to be removed by something like the following subroutine:
744 jmc 1.144 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
745     C & 'phiHydLow', myTime, myThid )
746 mlosch 1.90 Cml)
747 adcroft 1.69
748 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
749 jmc 1.130 IF ( useDiagnostics ) THEN
750 jmc 1.113
751     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
752 jmc 1.120 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
753 molod 1.116
754 jmc 1.120 tmpFac = 1. _d 0
755     CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
756     & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
757 molod 1.116
758 jmc 1.120 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
759     & 'PHIBOTSQ',0, 1,0,1,1,myThid)
760 jmc 1.113
761     ENDIF
762     #endif /* ALLOW_DIAGNOSTICS */
763 jmc 1.136
764 edhill 1.104 #ifdef ALLOW_DEBUG
765 jmc 1.158 IF ( debugLevel .GE. debLevD ) THEN
766 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
767 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
768 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
769     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
770     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
771     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
772 jmc 1.115 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
773     CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
774     CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
775     CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
776     #ifndef ALLOW_ADAMSBASHFORTH_3
777     CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
778     CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
779     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
780     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
781     #endif
782 adcroft 1.70 ENDIF
783 adcroft 1.69 #endif
784 cnh 1.1
785 jmc 1.125 #ifdef DYNAMICS_GUGV_EXCH_CHECK
786 jmc 1.144 C- jmc: For safety checking only: This Exchange here should not change
787     C the solution. If solution changes, it means something is wrong,
788 jmc 1.125 C but it does not mean that it is less wrong with this exchange.
789 jmc 1.158 IF ( debugLevel .GE. debLevE ) THEN
790 jmc 1.125 CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
791     ENDIF
792     #endif
793    
794 jmc 1.123 #ifdef ALLOW_DEBUG
795 jmc 1.153 IF (debugMode) CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
796 jmc 1.123 #endif
797    
798 cnh 1.1 RETURN
799     END

  ViewVC Help
Powered by ViewVC 1.1.22