/[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.4 - (hide annotations) (download)
Mon Sep 10 01:22:48 2001 UTC (22 years, 9 months ago) by adcroft
Branch: MAIN
Changes since 1.3: +21 -1 lines
Added multi-dimensional form of advection
 o available only for single step schemes (ie. can't be used with ABII)
 o stable for max(cfl_u,cfl_v,cfl_w)<=1  (without cfl_u+cfl_v+cfl_w <=1)
 o selected using multiDimAdvection=.T.  (default)
 o had to hack some existing routines to work on local arrays

1 adcroft 1.4 C $Header: /u/gcmpack/models/MITgcmUV/model/src/thermodynamics.F,v 1.3 2001/08/15 15:51:46 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     ENDDO
242     ENDDO
243     ENDDO
244    
245     iMin = 1-OLx+1
246     iMax = sNx+OLx
247     jMin = 1-OLy+1
248     jMax = sNy+OLy
249    
250    
251     #ifdef ALLOW_AUTODIFF_TAMC
252     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
253     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
254     #endif /* ALLOW_AUTODIFF_TAMC */
255    
256     C-- Start of diagnostic loop
257     DO k=Nr,1,-1
258    
259     #ifdef ALLOW_AUTODIFF_TAMC
260     C? Patrick, is this formula correct now that we change the loop range?
261     C? Do we still need this?
262     cph kkey formula corrected.
263     cph Needed for rhok, rhokm1, in the case useGMREDI.
264     kkey = (ikey-1)*Nr + k
265     CADJ STORE rhokm1(:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
266     CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
267     #endif /* ALLOW_AUTODIFF_TAMC */
268    
269     C-- Integrate continuity vertically for vertical velocity
270     CALL INTEGRATE_FOR_W(
271     I bi, bj, k, uVel, vVel,
272     O wVel,
273     I myThid )
274    
275     #ifdef ALLOW_OBCS
276     #ifdef ALLOW_NONHYDROSTATIC
277     C-- Apply OBC to W if in N-H mode
278     IF (useOBCS.AND.nonHydrostatic) THEN
279     CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
280     ENDIF
281     #endif /* ALLOW_NONHYDROSTATIC */
282     #endif /* ALLOW_OBCS */
283    
284     C-- Calculate gradients of potential density for isoneutral
285     C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
286     c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
287     IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
288     #ifdef ALLOW_AUTODIFF_TAMC
289     CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
290     CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
291     #endif /* ALLOW_AUTODIFF_TAMC */
292     CALL FIND_RHO(
293     I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
294     I theta, salt,
295     O rhoK,
296     I myThid )
297     IF (k.GT.1) THEN
298     #ifdef ALLOW_AUTODIFF_TAMC
299     CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
300     CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
301     #endif /* ALLOW_AUTODIFF_TAMC */
302     CALL FIND_RHO(
303     I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
304     I theta, salt,
305     O rhoKm1,
306     I myThid )
307     ENDIF
308     CALL GRAD_SIGMA(
309     I bi, bj, iMin, iMax, jMin, jMax, k,
310     I rhoK, rhoKm1, rhoK,
311     O sigmaX, sigmaY, sigmaR,
312     I myThid )
313     ENDIF
314    
315     C-- Implicit Vertical Diffusion for Convection
316     c ==> should use sigmaR !!!
317     IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
318     CALL CALC_IVDC(
319     I bi, bj, iMin, iMax, jMin, jMax, k,
320     I rhoKm1, rhoK,
321     U ConvectCount, KappaRT, KappaRS,
322     I myTime, myIter, myThid)
323     ENDIF
324    
325     C-- end of diagnostic k loop (Nr:1)
326     ENDDO
327    
328     #ifdef ALLOW_AUTODIFF_TAMC
329     cph avoids recomputation of integrate_for_w
330     CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
331     #endif /* ALLOW_AUTODIFF_TAMC */
332    
333     #ifdef ALLOW_OBCS
334     C-- Calculate future values on open boundaries
335     IF (useOBCS) THEN
336     CALL OBCS_CALC( bi, bj, myTime+deltaT,
337     I uVel, vVel, wVel, theta, salt,
338     I myThid )
339     ENDIF
340     #endif /* ALLOW_OBCS */
341    
342     C-- Determines forcing terms based on external fields
343     C relaxation terms, etc.
344     CALL EXTERNAL_FORCING_SURF(
345     I bi, bj, iMin, iMax, jMin, jMax,
346     I myThid )
347     #ifdef ALLOW_AUTODIFF_TAMC
348     cph needed for KPP
349     CADJ STORE surfacetendencyU(:,:,bi,bj)
350     CADJ & = comlev1_bibj, key=ikey, byte=isbyte
351     CADJ STORE surfacetendencyV(:,:,bi,bj)
352     CADJ & = comlev1_bibj, key=ikey, byte=isbyte
353     CADJ STORE surfacetendencyS(:,:,bi,bj)
354     CADJ & = comlev1_bibj, key=ikey, byte=isbyte
355     CADJ STORE surfacetendencyT(:,:,bi,bj)
356     CADJ & = comlev1_bibj, key=ikey, byte=isbyte
357     #endif /* ALLOW_AUTODIFF_TAMC */
358    
359     #ifdef ALLOW_GMREDI
360    
361     #ifdef ALLOW_AUTODIFF_TAMC
362     CADJ STORE sigmaX(:,:,:) = comlev1, key=ikey, byte=isbyte
363     CADJ STORE sigmaY(:,:,:) = comlev1, key=ikey, byte=isbyte
364     CADJ STORE sigmaR(:,:,:) = comlev1, key=ikey, byte=isbyte
365     #endif /* ALLOW_AUTODIFF_TAMC */
366     C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
367     IF (useGMRedi) THEN
368     DO k=1,Nr
369     CALL GMREDI_CALC_TENSOR(
370     I bi, bj, iMin, iMax, jMin, jMax, k,
371     I sigmaX, sigmaY, sigmaR,
372     I myThid )
373     ENDDO
374     #ifdef ALLOW_AUTODIFF_TAMC
375     ELSE
376     DO k=1, Nr
377     CALL GMREDI_CALC_TENSOR_DUMMY(
378     I bi, bj, iMin, iMax, jMin, jMax, k,
379     I sigmaX, sigmaY, sigmaR,
380     I myThid )
381     ENDDO
382     #endif /* ALLOW_AUTODIFF_TAMC */
383     ENDIF
384    
385     #ifdef ALLOW_AUTODIFF_TAMC
386     CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
387     CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
388     CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
389     #endif /* ALLOW_AUTODIFF_TAMC */
390    
391     #endif /* ALLOW_GMREDI */
392    
393     #ifdef ALLOW_KPP
394     C-- Compute KPP mixing coefficients
395     IF (useKPP) THEN
396     CALL KPP_CALC(
397     I bi, bj, myTime, myThid )
398     #ifdef ALLOW_AUTODIFF_TAMC
399     ELSE
400     CALL KPP_CALC_DUMMY(
401     I bi, bj, myTime, myThid )
402     #endif /* ALLOW_AUTODIFF_TAMC */
403     ENDIF
404    
405     #ifdef ALLOW_AUTODIFF_TAMC
406     CADJ STORE KPPghat (:,:,:,bi,bj)
407     CADJ & , KPPviscAz (:,:,:,bi,bj)
408     CADJ & , KPPdiffKzT(:,:,:,bi,bj)
409     CADJ & , KPPdiffKzS(:,:,:,bi,bj)
410     CADJ & , KPPfrac (:,: ,bi,bj)
411     CADJ & = comlev1_bibj, key=ikey, byte=isbyte
412     #endif /* ALLOW_AUTODIFF_TAMC */
413    
414     #endif /* ALLOW_KPP */
415    
416     #ifdef ALLOW_AUTODIFF_TAMC
417     CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
418     CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
419     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
420     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
421     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
422     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
423     #ifdef ALLOW_PASSIVE_TRACER
424     CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
425     #endif
426     #endif /* ALLOW_AUTODIFF_TAMC */
427    
428     #ifdef ALLOW_AIM
429     C AIM - atmospheric intermediate model, physics package code.
430     C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
431     IF ( useAIM ) THEN
432     CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
433     CALL AIM_DO_ATMOS_PHYSICS( phiHyd, bi, bj, myTime, myThid )
434     CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
435     ENDIF
436     #endif /* ALLOW_AIM */
437 adcroft 1.4
438     C-- Some advection schemes are better calculated using a multi-dimensional
439     C method in the absence of any other terms and, if used, is done here.
440     IF (multiDimAdvection) THEN
441     IF (tempStepping .AND.
442     & tempAdvScheme.NE.ENUM_CENTERED_2ND .AND.
443     & tempAdvScheme.NE.ENUM_UPWIND_3RD .AND.
444     & tempAdvScheme.NE.ENUM_CENTERED_4TH )
445     & CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,theta,
446     U gT,
447     I myTime,myIter,myThid)
448     IF (saltStepping .AND.
449     & saltAdvScheme.NE.ENUM_CENTERED_2ND .AND.
450     & saltAdvScheme.NE.ENUM_UPWIND_3RD .AND.
451     & saltAdvScheme.NE.ENUM_CENTERED_4TH )
452     & CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,salt,
453     U gS,
454     I myTime,myIter,myThid)
455     ENDIF
456 adcroft 1.1
457    
458     C-- Start of thermodynamics loop
459     DO k=Nr,1,-1
460     #ifdef ALLOW_AUTODIFF_TAMC
461     C? Patrick Is this formula correct?
462     cph Yes, but I rewrote it.
463     cph Also, the KappaR? need the index and subscript k!
464     kkey = (ikey-1)*Nr + k
465     #endif /* ALLOW_AUTODIFF_TAMC */
466    
467     C-- km1 Points to level above k (=k-1)
468     C-- kup Cycles through 1,2 to point to layer above
469     C-- kDown Cycles through 2,1 to point to current layer
470    
471     km1 = MAX(1,k-1)
472     kup = 1+MOD(k+1,2)
473     kDown= 1+MOD(k,2)
474    
475     iMin = 1-OLx
476     iMax = sNx+OLx
477     jMin = 1-OLy
478     jMax = sNy+OLy
479    
480     C-- Get temporary terms used by tendency routines
481     CALL CALC_COMMON_FACTORS (
482     I bi,bj,iMin,iMax,jMin,jMax,k,
483     O xA,yA,uTrans,vTrans,rTrans,maskUp,
484     I myThid)
485    
486     #ifdef ALLOW_AUTODIFF_TAMC
487     CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
488     CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
489     #endif /* ALLOW_AUTODIFF_TAMC */
490    
491     #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
492     C-- Calculate the total vertical diffusivity
493     CALL CALC_DIFFUSIVITY(
494     I bi,bj,iMin,iMax,jMin,jMax,k,
495     I maskUp,
496 heimbach 1.2 O KappaRT,KappaRS,
497 adcroft 1.1 I myThid)
498     #endif
499    
500     iMin = 1-OLx+2
501     iMax = sNx+OLx-1
502     jMin = 1-OLy+2
503     jMax = sNy+OLy-1
504    
505     C-- Calculate active tracer tendencies (gT,gS,...)
506     C and step forward storing result in gTnm1, gSnm1, etc.
507     IF ( tempStepping ) THEN
508     CALL CALC_GT(
509     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
510     I xA,yA,uTrans,vTrans,rTrans,maskUp,
511     I KappaRT,
512     U fVerT,
513     I myTime, myThid)
514     CALL TIMESTEP_TRACER(
515 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
516 adcroft 1.1 I theta, gT,
517     U gTnm1,
518     I myIter, myThid)
519     ENDIF
520     IF ( saltStepping ) THEN
521     CALL CALC_GS(
522     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
523     I xA,yA,uTrans,vTrans,rTrans,maskUp,
524     I KappaRS,
525     U fVerS,
526     I myTime, myThid)
527     CALL TIMESTEP_TRACER(
528 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
529 adcroft 1.1 I salt, gS,
530     U gSnm1,
531     I myIter, myThid)
532     ENDIF
533     #ifdef ALLOW_PASSIVE_TRACER
534     IF ( tr1Stepping ) THEN
535     CALL CALC_GTR1(
536     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
537     I xA,yA,uTrans,vTrans,rTrans,maskUp,
538     I KappaRT,
539     U fVerTr1,
540     I myTime, myThid)
541     CALL TIMESTEP_TRACER(
542 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,tracerAdvScheme,
543 adcroft 1.1 I Tr1, gTr1,
544     U gTr1NM1,
545     I myIter, myThid)
546     ENDIF
547     #endif
548    
549     #ifdef ALLOW_OBCS
550     C-- Apply open boundary conditions
551     IF (useOBCS) THEN
552     CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
553     END IF
554     #endif /* ALLOW_OBCS */
555    
556     C-- Freeze water
557     IF (allowFreezing) THEN
558     #ifdef ALLOW_AUTODIFF_TAMC
559     CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
560     CADJ & , key = kkey, byte = isbyte
561     #endif /* ALLOW_AUTODIFF_TAMC */
562     CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
563     END IF
564    
565     C-- end of thermodynamic k loop (Nr:1)
566     ENDDO
567    
568    
569     #ifdef ALLOW_AUTODIFF_TAMC
570     C? Patrick? What about this one?
571     cph Keys iikey and idkey don't seem to be needed
572     cph since storing occurs on different tape for each
573     cph impldiff call anyways.
574     cph Thus, common block comlev1_impl isn't needed either.
575     cph Storing below needed in the case useGMREDI.
576     iikey = (ikey-1)*maximpl
577     #endif /* ALLOW_AUTODIFF_TAMC */
578    
579     C-- Implicit diffusion
580     IF (implicitDiffusion) THEN
581    
582     IF (tempStepping) THEN
583     #ifdef ALLOW_AUTODIFF_TAMC
584     idkey = iikey + 1
585     CADJ STORE gTNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
586     #endif /* ALLOW_AUTODIFF_TAMC */
587     CALL IMPLDIFF(
588     I bi, bj, iMin, iMax, jMin, jMax,
589     I deltaTtracer, KappaRT, recip_HFacC,
590     U gTNm1,
591     I myThid )
592     ENDIF
593    
594     IF (saltStepping) THEN
595     #ifdef ALLOW_AUTODIFF_TAMC
596     idkey = iikey + 2
597     CADJ STORE gSNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
598     #endif /* ALLOW_AUTODIFF_TAMC */
599     CALL IMPLDIFF(
600     I bi, bj, iMin, iMax, jMin, jMax,
601     I deltaTtracer, KappaRS, recip_HFacC,
602     U gSNm1,
603     I myThid )
604     ENDIF
605    
606     #ifdef ALLOW_PASSIVE_TRACER
607     IF (tr1Stepping) THEN
608     #ifdef ALLOW_AUTODIFF_TAMC
609     CADJ STORE gTr1Nm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
610     #endif /* ALLOW_AUTODIFF_TAMC */
611     CALL IMPLDIFF(
612     I bi, bj, iMin, iMax, jMin, jMax,
613     I deltaTtracer, KappaRT, recip_HFacC,
614     U gTr1Nm1,
615     I myThid )
616     ENDIF
617     #endif
618    
619     #ifdef ALLOW_OBCS
620     C-- Apply open boundary conditions
621     IF (useOBCS) THEN
622     DO K=1,Nr
623     CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
624     ENDDO
625     END IF
626     #endif /* ALLOW_OBCS */
627    
628     C-- End If implicitDiffusion
629     ENDIF
630    
631     Ccs-
632     ENDDO
633     ENDDO
634    
635     #ifdef ALLOW_AIM
636     IF ( useAIM ) THEN
637     CALL AIM_AIM2DYN_EXCHANGES( myTime, myThid )
638     ENDIF
639     _EXCH_XYZ_R8(gTnm1,myThid)
640     _EXCH_XYZ_R8(gSnm1,myThid)
641     #else
642     IF (staggerTimeStep.AND.useCubedSphereExchange) THEN
643     _EXCH_XYZ_R8(gTnm1,myThid)
644     _EXCH_XYZ_R8(gSnm1,myThid)
645     ENDIF
646     #endif /* ALLOW_AIM */
647    
648     RETURN
649     END

  ViewVC Help
Powered by ViewVC 1.1.22