/[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.66 - (hide annotations) (download)
Sun Mar 25 22:33:52 2001 UTC (23 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint38, c37_adj
Changes since 1.65: +58 -21 lines
Modifications and additions to enable automatic differentiation.
Detailed info's in doc/notes_c37_adj.txt

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

  ViewVC Help
Powered by ViewVC 1.1.22