/[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.62 - (hide annotations) (download)
Wed Feb 14 22:51:27 2001 UTC (23 years, 3 months ago) by jmc
Branch: MAIN
Changes since 1.61: +37 -9 lines
recover (after checkpoint35) time average output

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

  ViewVC Help
Powered by ViewVC 1.1.22