/[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.73 - (hide annotations) (download)
Fri Jul 20 19:16:28 2001 UTC (22 years, 10 months ago) by adcroft
Branch: MAIN
Changes since 1.72: +3 -2 lines
Missing diag call for uVel.

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

  ViewVC Help
Powered by ViewVC 1.1.22