/[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.61 - (hide annotations) (download)
Wed Feb 7 21:48:02 2001 UTC (23 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint35
Changes since 1.60: +8 -15 lines
remove unused array "rVel"

1 jmc 1.61 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.60 2001/02/07 16:28:54 adcroft Exp $
2 adcroft 1.60 C $Name: $
3 cnh 1.1
4 adcroft 1.24 #include "CPP_OPTIONS.h"
5 cnh 1.1
6 cnh 1.8 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
7 cnh 1.1 C /==========================================================\
8     C | SUBROUTINE DYNAMICS |
9     C | o Controlling routine for the explicit part of the model |
10     C | dynamics. |
11     C |==========================================================|
12     C | This routine evaluates the "dynamics" terms for each |
13     C | block of ocean in turn. Because the blocks of ocean have |
14     C | overlap regions they are independent of one another. |
15     C | If terms involving lateral integrals are needed in this |
16     C | routine care will be needed. Similarly finite-difference |
17     C | operations with stencils wider than the overlap region |
18     C | require special consideration. |
19     C | Notes |
20     C | ===== |
21     C | C*P* comments indicating place holders for which code is |
22     C | presently being developed. |
23     C \==========================================================/
24 adcroft 1.40 IMPLICIT NONE
25 cnh 1.1
26     C == Global variables ===
27     #include "SIZE.h"
28     #include "EEPARAMS.h"
29     #include "CG2D.h"
30 adcroft 1.6 #include "PARAMS.h"
31 adcroft 1.3 #include "DYNVARS.h"
32 adcroft 1.42 #include "GRID.h"
33 heimbach 1.49
34     #ifdef ALLOW_AUTODIFF_TAMC
35 heimbach 1.53 # include "tamc.h"
36     # include "tamc_keys.h"
37     #endif /* ALLOW_AUTODIFF_TAMC */
38 heimbach 1.49
39 adcroft 1.58 #ifdef ALLOW_KPP
40     # include "KPP.h"
41     #endif
42    
43 cnh 1.1 C == Routine arguments ==
44 cnh 1.8 C myTime - Current time in simulation
45     C myIter - Current iteration number in simulation
46 cnh 1.1 C myThid - Thread number for this instance of the routine.
47 cnh 1.8 _RL myTime
48     INTEGER myIter
49 adcroft 1.47 INTEGER myThid
50 cnh 1.1
51     C == Local variables
52     C xA, yA - Per block temporaries holding face areas
53 cnh 1.38 C uTrans, vTrans, rTrans - Per block temporaries holding flow
54     C transport
55 jmc 1.61 C o uTrans: Zonal transport
56 cnh 1.1 C o vTrans: Meridional transport
57 cnh 1.30 C o rTrans: Vertical transport
58 cnh 1.1 C maskC,maskUp o maskC: land/water mask for tracer cells
59     C o maskUp: land/water mask for W points
60 adcroft 1.58 C fVer[STUV] o fVer: Vertical flux term - note fVer
61 cnh 1.1 C is "pipelined" in the vertical
62     C so we need an fVer for each
63     C variable.
64 adcroft 1.58 C rhoK, rhoKM1 - Density at current level, and level above
65 cnh 1.31 C phiHyd - Hydrostatic part of the potential phiHydi.
66 cnh 1.38 C In z coords phiHydiHyd is the hydrostatic
67     C pressure anomaly
68     C In p coords phiHydiHyd is the geopotential
69     C surface height
70 cnh 1.30 C anomaly.
71     C etaSurfX, - Holds surface elevation gradient in X and Y.
72     C etaSurfY
73     C KappaRT, - Total diffusion in vertical for T and S.
74 cnh 1.38 C KappaRS (background + spatially varying, isopycnal term).
75 cnh 1.30 C iMin, iMax - Ranges and sub-block indices on which calculations
76     C jMin, jMax are applied.
77 cnh 1.1 C bi, bj
78 heimbach 1.53 C k, kup, - Index for layer above and below. kup and kDown
79     C kDown, km1 are switched with layer to be the appropriate
80 cnh 1.38 C index into fVerTerm.
81 cnh 1.30 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86     _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
89     _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
90     _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
91     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
92 cnh 1.31 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
93 cnh 1.30 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94     _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 cnh 1.31 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
96     _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
97 adcroft 1.42 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
98     _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
99 adcroft 1.50 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
100     _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
101     _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
102 adcroft 1.12
103 adcroft 1.52 C This is currently also used by IVDC and Diagnostics
104     C #ifdef INCLUDE_CONVECT_CALL
105 adcroft 1.45 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
106 adcroft 1.52 C #endif
107 adcroft 1.45
108 cnh 1.1 INTEGER iMin, iMax
109     INTEGER jMin, jMax
110     INTEGER bi, bj
111     INTEGER i, j
112 heimbach 1.53 INTEGER k, km1, kup, kDown
113 cnh 1.1
114 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
115     INTEGER isbyte
116     PARAMETER( isbyte = 4 )
117    
118     INTEGER act1, act2, act3, act4
119     INTEGER max1, max2, max3
120 heimbach 1.51 INTEGER iikey, kkey
121 heimbach 1.49 INTEGER maximpl
122 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
123 heimbach 1.49
124 adcroft 1.11 C--- The algorithm...
125     C
126     C "Correction Step"
127     C =================
128     C Here we update the horizontal velocities with the surface
129     C pressure such that the resulting flow is either consistent
130     C with the free-surface evolution or the rigid-lid:
131     C U[n] = U* + dt x d/dx P
132     C V[n] = V* + dt x d/dy P
133     C
134     C "Calculation of Gs"
135     C ===================
136     C This is where all the accelerations and tendencies (ie.
137 heimbach 1.53 C physics, parameterizations etc...) are calculated
138 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
139 cnh 1.27 C b = b(rho, theta)
140 adcroft 1.11 C K31 = K31 ( rho )
141 jmc 1.61 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
142     C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
143     C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
144     C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
145 adcroft 1.11 C
146 adcroft 1.12 C "Time-stepping" or "Prediction"
147 adcroft 1.11 C ================================
148     C The models variables are stepped forward with the appropriate
149     C time-stepping scheme (currently we use Adams-Bashforth II)
150     C - For momentum, the result is always *only* a "prediction"
151     C in that the flow may be divergent and will be "corrected"
152     C later with a surface pressure gradient.
153     C - Normally for tracers the result is the new field at time
154     C level [n+1} *BUT* in the case of implicit diffusion the result
155     C is also *only* a prediction.
156     C - We denote "predictors" with an asterisk (*).
157     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
158     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
159     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
160     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
161 adcroft 1.12 C With implicit diffusion:
162 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
163     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
164 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
165     C (1 + dt * K * d_zz) salt[n] = salt*
166 adcroft 1.11 C---
167    
168 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
169     C-- dummy statement to end declaration part
170     ikey = 1
171 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
172 heimbach 1.49
173 cnh 1.1 C-- Set up work arrays with valid (i.e. not NaN) values
174     C These inital values do not alter the numerical results. They
175     C just ensure that all memory references are to valid floating
176     C point numbers. This prevents spurious hardware signals due to
177     C uninitialised but inert locations.
178     DO j=1-OLy,sNy+OLy
179     DO i=1-OLx,sNx+OLx
180 adcroft 1.5 xA(i,j) = 0. _d 0
181     yA(i,j) = 0. _d 0
182     uTrans(i,j) = 0. _d 0
183     vTrans(i,j) = 0. _d 0
184 heimbach 1.53 DO k=1,Nr
185 adcroft 1.58 phiHyd(i,j,k) = 0. _d 0
186 adcroft 1.45 KappaRU(i,j,k) = 0. _d 0
187     KappaRV(i,j,k) = 0. _d 0
188 adcroft 1.50 sigmaX(i,j,k) = 0. _d 0
189     sigmaY(i,j,k) = 0. _d 0
190     sigmaR(i,j,k) = 0. _d 0
191 cnh 1.1 ENDDO
192 cnh 1.30 rhoKM1 (i,j) = 0. _d 0
193     rhok (i,j) = 0. _d 0
194     maskC (i,j) = 0. _d 0
195 cnh 1.1 ENDDO
196     ENDDO
197    
198 cnh 1.35
199 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
200     C-- HPF directive to help TAMC
201 heimbach 1.53 CHPF$ INDEPENDENT
202     #endif /* ALLOW_AUTODIFF_TAMC */
203 heimbach 1.49
204 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
205 heimbach 1.49
206     #ifdef ALLOW_AUTODIFF_TAMC
207     C-- HPF directive to help TAMC
208 jmc 1.61 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS,fVerU,fVerV
209 heimbach 1.53 CHPF$& ,phiHyd,utrans,vtrans,maskc,xA,yA
210     CHPF$& ,KappaRT,KappaRS,KappaRU,KappaRV
211     CHPF$& )
212     #endif /* ALLOW_AUTODIFF_TAMC */
213 heimbach 1.49
214 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
215    
216 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
217     act1 = bi - myBxLo(myThid)
218     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
219    
220     act2 = bj - myByLo(myThid)
221     max2 = myByHi(myThid) - myByLo(myThid) + 1
222    
223     act3 = myThid - 1
224     max3 = nTx*nTy
225    
226     act4 = ikey_dynamics - 1
227    
228     ikey = (act1 + 1) + act2*max1
229     & + act3*max1*max2
230     & + act4*max1*max2*max3
231 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
232 heimbach 1.49
233 cnh 1.7 C-- Set up work arrays that need valid initial values
234     DO j=1-OLy,sNy+OLy
235     DO i=1-OLx,sNx+OLx
236 cnh 1.27 rTrans(i,j) = 0. _d 0
237 cnh 1.30 fVerT (i,j,1) = 0. _d 0
238     fVerT (i,j,2) = 0. _d 0
239     fVerS (i,j,1) = 0. _d 0
240     fVerS (i,j,2) = 0. _d 0
241     fVerU (i,j,1) = 0. _d 0
242     fVerU (i,j,2) = 0. _d 0
243     fVerV (i,j,1) = 0. _d 0
244     fVerV (i,j,2) = 0. _d 0
245 cnh 1.7 ENDDO
246     ENDDO
247    
248 adcroft 1.45 DO k=1,Nr
249     DO j=1-OLy,sNy+OLy
250     DO i=1-OLx,sNx+OLx
251     #ifdef INCLUDE_CONVECT_CALL
252     ConvectCount(i,j,k) = 0.
253     #endif
254     KappaRT(i,j,k) = 0. _d 0
255     KappaRS(i,j,k) = 0. _d 0
256     ENDDO
257     ENDDO
258     ENDDO
259    
260 cnh 1.1 iMin = 1-OLx+1
261     iMax = sNx+OLx
262     jMin = 1-OLy+1
263     jMax = sNy+OLy
264 cnh 1.35
265 adcroft 1.5
266 adcroft 1.58 C-- Start of diagnostic loop
267     DO k=Nr,1,-1
268 heimbach 1.49
269     #ifdef ALLOW_AUTODIFF_TAMC
270 adcroft 1.58 C? Patrick, is this formula correct now that we change the loop range?
271     C? Do we still need this?
272 heimbach 1.51 kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1
273 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
274 heimbach 1.49
275 adcroft 1.58 C-- Integrate continuity vertically for vertical velocity
276     CALL INTEGRATE_FOR_W(
277     I bi, bj, k, uVel, vVel,
278     O wVel,
279     I myThid )
280    
281     #ifdef ALLOW_OBCS
282     #ifdef ALLOW_NONHYDROSTATIC
283 adcroft 1.60 C-- Apply OBC to W if in N-H mode
284 adcroft 1.58 IF (useOBCS.AND.nonHydrostatic) THEN
285     CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
286     ENDIF
287     #endif /* ALLOW_NONHYDROSTATIC */
288     #endif /* ALLOW_OBCS */
289    
290     C-- Calculate gradients of potential density for isoneutral
291     C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
292     c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
293     IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
294     CALL FIND_RHO(
295     I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
296     I theta, salt,
297     O rhoK,
298 cnh 1.30 I myThid )
299 adcroft 1.58 IF (k.GT.1) CALL FIND_RHO(
300 heimbach 1.53 I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
301 adcroft 1.58 I theta, salt,
302     O rhoKm1,
303 cnh 1.30 I myThid )
304 adcroft 1.58 CALL GRAD_SIGMA(
305 heimbach 1.53 I bi, bj, iMin, iMax, jMin, jMax, k,
306 adcroft 1.58 I rhoK, rhoKm1, rhoK,
307 adcroft 1.50 O sigmaX, sigmaY, sigmaR,
308     I myThid )
309 adcroft 1.58 ENDIF
310 heimbach 1.49
311 adcroft 1.58 C-- Implicit Vertical Diffusion for Convection
312     c ==> should use sigmaR !!!
313     IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
314     CALL CALC_IVDC(
315     I bi, bj, iMin, iMax, jMin, jMax, k,
316     I rhoKm1, rhoK,
317     U ConvectCount, KappaRT, KappaRS,
318     I myTime, myIter, myThid)
319     END IF
320 heimbach 1.53
321 adcroft 1.58 C-- end of diagnostic k loop (Nr:1)
322 heimbach 1.49 ENDDO
323    
324 adcroft 1.58 #ifdef ALLOW_OBCS
325     C-- Calculate future values on open boundaries
326     IF (useOBCS) THEN
327     CALL OBCS_CALC( bi, bj, myTime+deltaT,
328     I uVel, vVel, wVel, theta, salt,
329     I myThid )
330     ENDIF
331     #endif /* ALLOW_OBCS */
332    
333     C-- Determines forcing terms based on external fields
334     C relaxation terms, etc.
335     CALL EXTERNAL_FORCING_SURF(
336 heimbach 1.54 I bi, bj, iMin, iMax, jMin, jMax,
337     I myThid )
338    
339 adcroft 1.58 #ifdef ALLOW_GMREDI
340     C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
341 heimbach 1.53 IF (useGMRedi) THEN
342 adcroft 1.58 DO k=1,Nr
343 heimbach 1.53 CALL GMREDI_CALC_TENSOR(
344     I bi, bj, iMin, iMax, jMin, jMax, k,
345 adcroft 1.50 I sigmaX, sigmaY, sigmaR,
346     I myThid )
347 heimbach 1.53 ENDDO
348 heimbach 1.55 #ifdef ALLOW_AUTODIFF_TAMC
349     ELSE
350     DO k=1, Nr
351     CALL GMREDI_CALC_TENSOR_DUMMY(
352     I bi, bj, iMin, iMax, jMin, jMax, k,
353     I sigmaX, sigmaY, sigmaR,
354     I myThid )
355     ENDDO
356     #endif /* ALLOW_AUTODIFF_TAMC */
357 heimbach 1.53 ENDIF
358 adcroft 1.58 #endif /* ALLOW_GMREDI */
359 heimbach 1.53
360 adcroft 1.58 #ifdef ALLOW_KPP
361     C-- Compute KPP mixing coefficients
362 heimbach 1.53 IF (useKPP) THEN
363     CALL KPP_CALC(
364 heimbach 1.54 I bi, bj, myTime, myThid )
365 adcroft 1.58 ENDIF
366     #endif /* ALLOW_KPP */
367 heimbach 1.53
368     #ifdef ALLOW_AUTODIFF_TAMC
369 adcroft 1.58 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
370     CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
371     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
372     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
373     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
374     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
375     #endif /* ALLOW_AUTODIFF_TAMC */
376    
377     #ifdef ALLOW_AIM
378     C AIM - atmospheric intermediate model, physics package code.
379     C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
380     IF ( useAIM ) THEN
381     CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
382     CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid )
383     CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
384 heimbach 1.53 ENDIF
385 adcroft 1.58 #endif /* ALLOW_AIM */
386    
387 heimbach 1.53
388 adcroft 1.58 C-- Start of thermodynamics loop
389     DO k=Nr,1,-1
390    
391     C-- km1 Points to level above k (=k-1)
392     C-- kup Cycles through 1,2 to point to layer above
393     C-- kDown Cycles through 2,1 to point to current layer
394    
395     km1 = MAX(1,k-1)
396     kup = 1+MOD(k+1,2)
397     kDown= 1+MOD(k,2)
398    
399     iMin = 1-OLx+2
400     iMax = sNx+OLx-1
401     jMin = 1-OLy+2
402     jMax = sNy+OLy-1
403 cnh 1.1
404 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
405 adcroft 1.58 CPatrick Is this formula correct?
406 heimbach 1.51 kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
407 adcroft 1.58 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
408     CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
409     CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
410 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
411 heimbach 1.49
412 cnh 1.1 C-- Get temporary terms used by tendency routines
413     CALL CALC_COMMON_FACTORS (
414 heimbach 1.53 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
415 jmc 1.61 O xA,yA,uTrans,vTrans,rTrans,maskC,maskUp,
416 cnh 1.1 I myThid)
417 heimbach 1.49
418 cnh 1.38 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
419 adcroft 1.12 C-- Calculate the total vertical diffusivity
420     CALL CALC_DIFFUSIVITY(
421 heimbach 1.53 I bi,bj,iMin,iMax,jMin,jMax,k,
422 adcroft 1.58 I maskC,maskup,
423 adcroft 1.42 O KappaRT,KappaRS,KappaRU,KappaRV,
424 adcroft 1.12 I myThid)
425 cnh 1.38 #endif
426 adcroft 1.58
427     C-- Calculate active tracer tendencies (gT,gS,...)
428     C and step forward storing result in gTnm1, gSnm1, etc.
429 cnh 1.9 IF ( tempStepping ) THEN
430 adcroft 1.58 CALL CALC_GT(
431 heimbach 1.53 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
432 cnh 1.30 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
433 adcroft 1.50 I KappaRT,
434 adcroft 1.58 U fVerT,
435 cnh 1.37 I myTime, myThid)
436 adcroft 1.58 CALL TIMESTEP_TRACER(
437     I bi,bj,iMin,iMax,jMin,jMax,k,
438     I theta, gT,
439     U gTnm1,
440     I myIter, myThid)
441 cnh 1.9 ENDIF
442 adcroft 1.18 IF ( saltStepping ) THEN
443 adcroft 1.58 CALL CALC_GS(
444 heimbach 1.53 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
445 cnh 1.30 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
446 adcroft 1.50 I KappaRS,
447 adcroft 1.58 U fVerS,
448 cnh 1.37 I myTime, myThid)
449 adcroft 1.58 CALL TIMESTEP_TRACER(
450     I bi,bj,iMin,iMax,jMin,jMax,k,
451     I salt, gS,
452     U gSnm1,
453     I myIter, myThid)
454 adcroft 1.18 ENDIF
455 adcroft 1.58
456     #ifdef ALLOW_OBCS
457 adcroft 1.41 C-- Apply open boundary conditions
458 adcroft 1.58 IF (useOBCS) THEN
459     CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
460     END IF
461     #endif /* ALLOW_OBCS */
462 heimbach 1.54
463 adcroft 1.41 C-- Freeze water
464 heimbach 1.49 IF (allowFreezing) THEN
465     #ifdef ALLOW_AUTODIFF_TAMC
466 adcroft 1.58 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
467     CADJ & , key = kkey, byte = isbyte
468 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
469     CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
470 heimbach 1.49 END IF
471 adcroft 1.48
472 adcroft 1.58 C-- end of thermodynamic k loop (Nr:1)
473     ENDDO
474 adcroft 1.45
475 adcroft 1.11
476 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
477 adcroft 1.58 CPatrick? What about this one?
478 heimbach 1.49 maximpl = 6
479 heimbach 1.51 iikey = (ikey-1)*maximpl
480 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
481 heimbach 1.51
482     C-- Implicit diffusion
483     IF (implicitDiffusion) THEN
484 heimbach 1.49
485 adcroft 1.58 IF (tempStepping) THEN
486 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
487     idkey = iikey + 1
488 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
489 heimbach 1.49 CALL IMPLDIFF(
490 adcroft 1.42 I bi, bj, iMin, iMax, jMin, jMax,
491 adcroft 1.58 I deltaTtracer, KappaRT, recip_HFacC,
492 adcroft 1.42 U gTNm1,
493     I myThid )
494 adcroft 1.58 ENDIF
495 heimbach 1.49
496     IF (saltStepping) THEN
497     #ifdef ALLOW_AUTODIFF_TAMC
498     idkey = iikey + 2
499 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
500 heimbach 1.49 CALL IMPLDIFF(
501 adcroft 1.42 I bi, bj, iMin, iMax, jMin, jMax,
502 adcroft 1.58 I deltaTtracer, KappaRS, recip_HFacC,
503 adcroft 1.42 U gSNm1,
504     I myThid )
505 adcroft 1.58 ENDIF
506    
507     #ifdef ALLOW_OBCS
508     C-- Apply open boundary conditions
509     IF (useOBCS) THEN
510     DO K=1,Nr
511     CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
512     ENDDO
513 heimbach 1.49 END IF
514 adcroft 1.58 #endif /* ALLOW_OBCS */
515 heimbach 1.49
516 adcroft 1.58 C-- End If implicitDiffusion
517 heimbach 1.53 ENDIF
518 heimbach 1.49
519 adcroft 1.58
520    
521     C-- Start of dynamics loop
522     DO k=1,Nr
523    
524     C-- km1 Points to level above k (=k-1)
525     C-- kup Cycles through 1,2 to point to layer above
526     C-- kDown Cycles through 2,1 to point to current layer
527    
528     km1 = MAX(1,k-1)
529     kup = 1+MOD(k+1,2)
530     kDown= 1+MOD(k,2)
531    
532     iMin = 1-OLx+2
533     iMax = sNx+OLx-1
534     jMin = 1-OLy+2
535     jMax = sNy+OLy-1
536    
537     C-- Integrate hydrostatic balance for phiHyd with BC of
538     C phiHyd(z=0)=0
539     C distinguishe between Stagger and Non Stagger time stepping
540     IF (staggerTimeStep) THEN
541     CALL CALC_PHI_HYD(
542     I bi,bj,iMin,iMax,jMin,jMax,k,
543     I gTnm1, gSnm1,
544     U phiHyd,
545     I myThid )
546     ELSE
547     CALL CALC_PHI_HYD(
548     I bi,bj,iMin,iMax,jMin,jMax,k,
549     I theta, salt,
550     U phiHyd,
551     I myThid )
552     ENDIF
553    
554     C-- Calculate accelerations in the momentum equations (gU, gV, ...)
555     C and step forward storing the result in gUnm1, gVnm1, etc...
556     IF ( momStepping ) THEN
557     CALL CALC_MOM_RHS(
558     I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
559     I phiHyd,KappaRU,KappaRV,
560     U fVerU, fVerV,
561     I myTime, myThid)
562     CALL TIMESTEP(
563     I bi,bj,iMin,iMax,jMin,jMax,k,phiHyd,
564     I myIter, myThid)
565    
566     #ifdef ALLOW_OBCS
567     C-- Apply open boundary conditions
568     IF (useOBCS) THEN
569     CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
570     END IF
571     #endif /* ALLOW_OBCS */
572    
573     #ifdef ALLOW_AUTODIFF_TAMC
574     #ifdef INCLUDE_CD_CODE
575     ELSE
576     DO j=1-OLy,sNy+OLy
577     DO i=1-OLx,sNx+OLx
578     guCD(i,j,k,bi,bj) = 0.0
579     gvCD(i,j,k,bi,bj) = 0.0
580     END DO
581     END DO
582     #endif /* INCLUDE_CD_CODE */
583     #endif /* ALLOW_AUTODIFF_TAMC */
584     ENDIF
585    
586    
587     C-- end of dynamics k loop (1:Nr)
588     ENDDO
589    
590    
591    
592 adcroft 1.44 C-- Implicit viscosity
593 adcroft 1.58 IF (implicitViscosity.AND.momStepping) THEN
594     #ifdef ALLOW_AUTODIFF_TAMC
595     idkey = iikey + 3
596     #endif /* ALLOW_AUTODIFF_TAMC */
597 adcroft 1.42 CALL IMPLDIFF(
598     I bi, bj, iMin, iMax, jMin, jMax,
599     I deltaTmom, KappaRU,recip_HFacW,
600     U gUNm1,
601     I myThid )
602 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
603     idkey = iikey + 4
604     #endif /* ALLOW_AUTODIFF_TAMC */
605 adcroft 1.42 CALL IMPLDIFF(
606     I bi, bj, iMin, iMax, jMin, jMax,
607     I deltaTmom, KappaRV,recip_HFacS,
608     U gVNm1,
609     I myThid )
610 heimbach 1.49
611 adcroft 1.58 #ifdef ALLOW_OBCS
612     C-- Apply open boundary conditions
613     IF (useOBCS) THEN
614     DO K=1,Nr
615     CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
616     ENDDO
617     END IF
618     #endif /* ALLOW_OBCS */
619 heimbach 1.49
620 adcroft 1.58 #ifdef INCLUDE_CD_CODE
621     #ifdef ALLOW_AUTODIFF_TAMC
622     idkey = iikey + 5
623     #endif /* ALLOW_AUTODIFF_TAMC */
624 adcroft 1.42 CALL IMPLDIFF(
625     I bi, bj, iMin, iMax, jMin, jMax,
626     I deltaTmom, KappaRU,recip_HFacW,
627     U vVelD,
628     I myThid )
629 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
630     idkey = iikey + 6
631     #endif /* ALLOW_AUTODIFF_TAMC */
632 adcroft 1.42 CALL IMPLDIFF(
633     I bi, bj, iMin, iMax, jMin, jMax,
634     I deltaTmom, KappaRV,recip_HFacS,
635     U uVelD,
636     I myThid )
637 adcroft 1.58 #endif /* INCLUDE_CD_CODE */
638     C-- End If implicitViscosity.AND.momStepping
639 heimbach 1.53 ENDIF
640 cnh 1.1
641     ENDDO
642     ENDDO
643    
644     RETURN
645     END

  ViewVC Help
Powered by ViewVC 1.1.22