/[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.5 - (hide annotations) (download)
Mon Sep 10 16:35:27 2001 UTC (22 years, 9 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint40pre9
Changes since 1.4: +24 -17 lines
Include initialisation of gTracer fields to break
dependency for TAF.

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

  ViewVC Help
Powered by ViewVC 1.1.22