/[MITgcm]/MITgcm/model/src/thermodynamics.F
ViewVC logotype

Annotation of /MITgcm/model/src/thermodynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.3 - (hide annotations) (download)
Wed Aug 15 15:51:46 2001 UTC (22 years, 9 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre8
Changes since 1.2: +4 -9 lines
Added run-time control of advection schemes.
 - advection scheme determines method of forward integration.
 - unfortunately, we have to use integers in "data" since ENUM_CENTERED_2ND
   doesn't mean anything to fortran
 - defaults are centered second
 - output differs due to these mods! This is due to the g77 optimization.
   I have tested that using -ffloat-store, these mods do not affect
   the output so am confident about changes.

                T           S           U           V
C D M    c        m  s        m  s        m  s        m  s
n p a R  g  m  m  e  .  m  m  e  .  m  m  e  .  m  m  e  .
f n k u  2  i  a  a  d  i  a  a  d  i  a  a  d  i  a  a  d
g d e n  d  n  x  n  .  n  x  n  .  n  x  n  .  n  x  n  .

Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 22 16 16 pass  adjustment.128x64x1
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  adjustment.cs-32x32x1
Y Y N N -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- N/O   aim.5l_Equatorial_Channel
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  aim.5l_LatLon
Y Y N N -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- N/O   aim.5l_zon-ave
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 22 16 pass  exp0
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 22 16 pass  exp1
Y Y Y Y 13 16 16 16 16 16 16 16 13 16 13 13 13 13 13 13 16 pass  exp2
Y Y Y Y 12 16 16 13 16 16 16 16 16 13 16 16 16 16 13 13 16 FAIL  exp4
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 22 16 16 16 22 16 pass  exp5
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  front_relax
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 22 16 16 16 11 16 pass  hs94.128x64x5
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 22 16 pass  hs94.1x64x5
Y Y Y Y 14 16 16 16 16 16 16 16 16 16 13 16 16 16 16 16 13 pass  hs94.cs-32x32x5
Y Y Y Y 16 16 16 22 16 16 16 16 16 16 16 22 16 16 16 16 16 pass  internal_wave
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  natl_box
Y Y Y Y 16 16 16 16 16 13 16 13 16 16 16 16 16 16 16 16 16 pass  solid-body.cs-32x32x1

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

  ViewVC Help
Powered by ViewVC 1.1.22