/[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.86 - (hide annotations) (download)
Sun Dec 5 22:23:17 2004 UTC (19 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57, checkpoint57a_post, checkpoint57a_pre
Changes since 1.85: +6 -4 lines
implement Implicit Vertical advection for pTracers

1 jmc 1.86 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.85 2004/12/04 05:59:50 dimitri Exp $
2 adcroft 1.1 C $Name: $
3    
4 edhill 1.51 #include "PACKAGES_CONFIG.h"
5 adcroft 1.1 #include "CPP_OPTIONS.h"
6 edhill 1.51
7 jmc 1.21 #ifdef ALLOW_AUTODIFF_TAMC
8     # ifdef ALLOW_GMREDI
9     # include "GMREDI_OPTIONS.h"
10     # endif
11     # ifdef ALLOW_KPP
12     # include "KPP_OPTIONS.h"
13 heimbach 1.42 # endif
14 jmc 1.21 #endif /* ALLOW_AUTODIFF_TAMC */
15 adcroft 1.1
16 cnh 1.9 CBOP
17     C !ROUTINE: THERMODYNAMICS
18     C !INTERFACE:
19 adcroft 1.1 SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
20 cnh 1.9 C !DESCRIPTION: \bv
21     C *==========================================================*
22     C | SUBROUTINE THERMODYNAMICS
23     C | o Controlling routine for the prognostic part of the
24     C | thermo-dynamics.
25     C *===========================================================
26     C | The algorithm...
27     C |
28     C | "Correction Step"
29     C | =================
30     C | Here we update the horizontal velocities with the surface
31     C | pressure such that the resulting flow is either consistent
32     C | with the free-surface evolution or the rigid-lid:
33     C | U[n] = U* + dt x d/dx P
34     C | V[n] = V* + dt x d/dy P
35     C |
36     C | "Calculation of Gs"
37     C | ===================
38     C | This is where all the accelerations and tendencies (ie.
39     C | physics, parameterizations etc...) are calculated
40     C | rho = rho ( theta[n], salt[n] )
41     C | b = b(rho, theta)
42     C | K31 = K31 ( rho )
43     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
44     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
45     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
46     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
47     C |
48     C | "Time-stepping" or "Prediction"
49     C | ================================
50     C | The models variables are stepped forward with the appropriate
51     C | time-stepping scheme (currently we use Adams-Bashforth II)
52     C | - For momentum, the result is always *only* a "prediction"
53     C | in that the flow may be divergent and will be "corrected"
54     C | later with a surface pressure gradient.
55     C | - Normally for tracers the result is the new field at time
56     C | level [n+1} *BUT* in the case of implicit diffusion the result
57     C | is also *only* a prediction.
58     C | - We denote "predictors" with an asterisk (*).
59     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
60     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
61     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
62     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
63     C | With implicit diffusion:
64     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
65     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66     C | (1 + dt * K * d_zz) theta[n] = theta*
67     C | (1 + dt * K * d_zz) salt[n] = salt*
68     C |
69     C *==========================================================*
70     C \ev
71    
72     C !USES:
73 adcroft 1.1 IMPLICIT NONE
74     C == Global variables ===
75     #include "SIZE.h"
76     #include "EEPARAMS.h"
77     #include "PARAMS.h"
78     #include "DYNVARS.h"
79     #include "GRID.h"
80 jmc 1.83 #ifdef ALLOW_GENERIC_ADVDIFF
81 adcroft 1.4 #include "GAD.h"
82 jmc 1.83 #endif
83 stephd 1.75 #ifdef ALLOW_OFFLINE
84     #include "OFFLINE.h"
85     #endif
86 jmc 1.45 #ifdef ALLOW_PTRACERS
87 jmc 1.74 #include "PTRACERS_SIZE.h"
88 jmc 1.45 #include "PTRACERS.h"
89     #endif
90 heimbach 1.42 #ifdef ALLOW_TIMEAVE
91     #include "TIMEAVE_STATV.h"
92     #endif
93    
94 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
95     # include "tamc.h"
96     # include "tamc_keys.h"
97     # include "FFIELDS.h"
98 heimbach 1.30 # include "EOS.h"
99 adcroft 1.1 # ifdef ALLOW_KPP
100     # include "KPP.h"
101 heimbach 1.67 # endif
102     # ifdef ALLOW_GMREDI
103     # include "GMREDI.h"
104 adcroft 1.1 # endif
105 heimbach 1.68 # ifdef ALLOW_EBM
106     # include "EBM.h"
107 adcroft 1.1 # endif
108     #endif /* ALLOW_AUTODIFF_TAMC */
109    
110 stephd 1.75
111 cnh 1.9 C !INPUT/OUTPUT PARAMETERS:
112 adcroft 1.1 C == Routine arguments ==
113     C myTime - Current time in simulation
114     C myIter - Current iteration number in simulation
115     C myThid - Thread number for this instance of the routine.
116     _RL myTime
117     INTEGER myIter
118     INTEGER myThid
119    
120 jmc 1.83 #ifdef ALLOW_GENERIC_ADVDIFF
121 cnh 1.9 C !LOCAL VARIABLES:
122 adcroft 1.1 C == Local variables
123     C xA, yA - Per block temporaries holding face areas
124     C uTrans, vTrans, rTrans - Per block temporaries holding flow
125     C transport
126     C o uTrans: Zonal transport
127     C o vTrans: Meridional transport
128     C o rTrans: Vertical transport
129 jmc 1.64 C rTransKp1 o vertical volume transp. at interface k+1
130 adcroft 1.1 C maskUp o maskUp: land/water mask for W points
131     C fVer[STUV] o fVer: Vertical flux term - note fVer
132     C is "pipelined" in the vertical
133     C so we need an fVer for each
134     C variable.
135 jmc 1.81 C kappaRT, - Total diffusion in vertical at level k, for T and S
136     C kappaRS (background + spatially varying, isopycnal term).
137     C kappaRTr - Total diffusion in vertical at level k,
138     C for each passive Tracer
139     C kappaRk - Total diffusion in vertical, all levels, 1 tracer
140 jmc 1.39 C useVariableK = T when vertical diffusion is not constant
141 adcroft 1.1 C iMin, iMax - Ranges and sub-block indices on which calculations
142     C jMin, jMax are applied.
143     C bi, bj
144     C k, kup, - Index for layer above and below. kup and kDown
145     C kDown, km1 are switched with layer to be the appropriate
146     C index into fVerTerm.
147     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152 jmc 1.64 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153 adcroft 1.1 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
155     _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
156 jmc 1.81 _RL kappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
157     _RL kappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
158 heimbach 1.55 #ifdef ALLOW_PTRACERS
159     _RL fVerP (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
160 jmc 1.81 _RL kappaRTr(1-Olx:sNx+Olx,1-Oly:sNy+Oly,PTRACERS_num)
161 heimbach 1.55 #endif
162 jmc 1.81 _RL kappaRk (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
163 jmc 1.64 _RL kp1Msk
164 jmc 1.39 LOGICAL useVariableK
165 adcroft 1.1 INTEGER iMin, iMax
166     INTEGER jMin, jMax
167     INTEGER bi, bj
168     INTEGER i, j
169     INTEGER k, km1, kup, kDown
170 heimbach 1.55 INTEGER iTracer, ip
171 adcroft 1.1
172 cnh 1.9 CEOP
173 adcroft 1.40
174 edhill 1.56 #ifdef ALLOW_DEBUG
175 heimbach 1.43 IF ( debugLevel .GE. debLevB )
176 jmc 1.63 & CALL DEBUG_ENTER('THERMODYNAMICS',myThid)
177 adcroft 1.40 #endif
178 adcroft 1.1
179     #ifdef ALLOW_AUTODIFF_TAMC
180     C-- dummy statement to end declaration part
181     ikey = 1
182 heimbach 1.30 itdkey = 1
183 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
184    
185     #ifdef ALLOW_AUTODIFF_TAMC
186     C-- HPF directive to help TAMC
187     CHPF$ INDEPENDENT
188     #endif /* ALLOW_AUTODIFF_TAMC */
189    
190     DO bj=myByLo(myThid),myByHi(myThid)
191    
192     #ifdef ALLOW_AUTODIFF_TAMC
193     C-- HPF directive to help TAMC
194 heimbach 1.2 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
195 jmc 1.37 CHPF$& ,utrans,vtrans,xA,yA
196 jmc 1.81 CHPF$& ,kappaRT,kappaRS
197 adcroft 1.1 CHPF$& )
198     #endif /* ALLOW_AUTODIFF_TAMC */
199    
200     DO bi=myBxLo(myThid),myBxHi(myThid)
201    
202     #ifdef ALLOW_AUTODIFF_TAMC
203     act1 = bi - myBxLo(myThid)
204     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
205     act2 = bj - myByLo(myThid)
206     max2 = myByHi(myThid) - myByLo(myThid) + 1
207     act3 = myThid - 1
208     max3 = nTx*nTy
209     act4 = ikey_dynamics - 1
210 heimbach 1.30 itdkey = (act1 + 1) + act2*max1
211 adcroft 1.1 & + act3*max1*max2
212     & + act4*max1*max2*max3
213     #endif /* ALLOW_AUTODIFF_TAMC */
214    
215 heimbach 1.41 C-- Set up work arrays with valid (i.e. not NaN) values
216     C These inital values do not alter the numerical results. They
217     C just ensure that all memory references are to valid floating
218     C point numbers. This prevents spurious hardware signals due to
219     C uninitialised but inert locations.
220    
221 adcroft 1.1 DO j=1-OLy,sNy+OLy
222     DO i=1-OLx,sNx+OLx
223 heimbach 1.41 xA(i,j) = 0. _d 0
224     yA(i,j) = 0. _d 0
225     uTrans(i,j) = 0. _d 0
226     vTrans(i,j) = 0. _d 0
227 adcroft 1.1 rTrans (i,j) = 0. _d 0
228 jmc 1.64 rTransKp1(i,j) = 0. _d 0
229 adcroft 1.1 fVerT (i,j,1) = 0. _d 0
230     fVerT (i,j,2) = 0. _d 0
231     fVerS (i,j,1) = 0. _d 0
232     fVerS (i,j,2) = 0. _d 0
233 jmc 1.81 kappaRT(i,j) = 0. _d 0
234     kappaRS(i,j) = 0. _d 0
235 heimbach 1.55 #ifdef ALLOW_PTRACERS
236     DO ip=1,PTRACERS_num
237     fVerP (i,j,1,ip) = 0. _d 0
238     fVerP (i,j,2,ip) = 0. _d 0
239 jmc 1.81 kappaRTr(i,j,ip) = 0. _d 0
240 heimbach 1.55 ENDDO
241     #endif
242 adcroft 1.1 ENDDO
243     ENDDO
244    
245     DO k=1,Nr
246     DO j=1-OLy,sNy+OLy
247     DO i=1-OLx,sNx+OLx
248     C This is currently also used by IVDC and Diagnostics
249 jmc 1.81 kappaRk(i,j,k) = 0. _d 0
250 jmc 1.45 C- tracer tendency needs to be set to zero (moved here from gad_calc_rhs):
251 heimbach 1.30 gT(i,j,k,bi,bj) = 0. _d 0
252     gS(i,j,k,bi,bj) = 0. _d 0
253 heimbach 1.42 # ifdef ALLOW_PTRACERS
254     DO iTracer=1,PTRACERS_numInUse
255     gPTr(i,j,k,bi,bj,itracer) = 0. _d 0
256     ENDDO
257     # endif
258 adcroft 1.1 ENDDO
259     ENDDO
260     ENDDO
261    
262 jmc 1.72 c iMin = 1-OLx
263     c iMax = sNx+OLx
264     c jMin = 1-OLy
265     c jMax = sNy+OLy
266 adcroft 1.1
267     #ifdef ALLOW_AUTODIFF_TAMC
268     cph avoids recomputation of integrate_for_w
269 heimbach 1.30 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
270 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
271    
272 heimbach 1.22 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
273     C-- MOST of THERMODYNAMICS will be disabled
274     #ifndef SINGLE_LAYER_MODE
275    
276 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
277 heimbach 1.30 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
278     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
279     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
280     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
281 heimbach 1.42 #ifdef ALLOW_PTRACERS
282     cph-- moved to forward_step to avoid key computation
283     cphCADJ STORE ptracer(:,:,:,bi,bj,itracer) = comlev1_bibj,
284     cphCADJ & key=itdkey, byte=isbyte
285     #endif
286 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
287    
288 adcroft 1.12 #ifndef DISABLE_MULTIDIM_ADVECTION
289 adcroft 1.4 C-- Some advection schemes are better calculated using a multi-dimensional
290     C method in the absence of any other terms and, if used, is done here.
291 adcroft 1.13 C
292     C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
293     C The default is to use multi-dimensinal advection for non-linear advection
294     C schemes. However, for the sake of efficiency of the adjoint it is necessary
295     C to be able to exclude this scheme to avoid excessive storage and
296     C recomputation. It *is* differentiable, if you need it.
297     C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
298     C disable this section of code.
299 stephd 1.75 #ifndef ALLOW_OFFLINE
300 jmc 1.24 IF (tempMultiDimAdvec) THEN
301 edhill 1.56 #ifdef ALLOW_DEBUG
302 heimbach 1.43 IF ( debugLevel .GE. debLevB )
303     & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
304 adcroft 1.40 #endif
305 jmc 1.63 CALL GAD_ADVECTION(
306 jmc 1.69 I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme,
307     I GAD_TEMPERATURE,
308 jmc 1.63 I uVel, vVel, wVel, theta,
309     O gT,
310     I bi,bj,myTime,myIter,myThid)
311 jmc 1.23 ENDIF
312 stephd 1.75 #endif
313     #ifndef ALLOW_OFFLINE
314 jmc 1.24 IF (saltMultiDimAdvec) THEN
315 edhill 1.56 #ifdef ALLOW_DEBUG
316 heimbach 1.43 IF ( debugLevel .GE. debLevB )
317     & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
318 adcroft 1.40 #endif
319 jmc 1.63 CALL GAD_ADVECTION(
320 jmc 1.69 I saltImplVertAdv, saltAdvScheme, saltVertAdvScheme,
321     I GAD_SALINITY,
322 jmc 1.63 I uVel, vVel, wVel, salt,
323     O gS,
324     I bi,bj,myTime,myIter,myThid)
325 adcroft 1.6 ENDIF
326 stephd 1.75 #endif
327 adcroft 1.17 C Since passive tracers are configurable separately from T,S we
328     C call the multi-dimensional method for PTRACERS regardless
329     C of whether multiDimAdvection is set or not.
330     #ifdef ALLOW_PTRACERS
331     IF ( usePTRACERS ) THEN
332 edhill 1.56 #ifdef ALLOW_DEBUG
333 heimbach 1.43 IF ( debugLevel .GE. debLevB )
334     & CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
335 adcroft 1.40 #endif
336 adcroft 1.17 CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
337     ENDIF
338     #endif /* ALLOW_PTRACERS */
339 adcroft 1.12 #endif /* DISABLE_MULTIDIM_ADVECTION */
340 adcroft 1.1
341 edhill 1.56 #ifdef ALLOW_DEBUG
342 jmc 1.63 IF ( debugLevel .GE. debLevB )
343 heimbach 1.43 & CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
344 adcroft 1.40 #endif
345    
346 adcroft 1.1 C-- Start of thermodynamics loop
347     DO k=Nr,1,-1
348     #ifdef ALLOW_AUTODIFF_TAMC
349     C? Patrick Is this formula correct?
350     cph Yes, but I rewrote it.
351 jmc 1.81 cph Also, the kappaR? need the index and subscript k!
352 heimbach 1.30 kkey = (itdkey-1)*Nr + k
353 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
354    
355     C-- km1 Points to level above k (=k-1)
356     C-- kup Cycles through 1,2 to point to layer above
357     C-- kDown Cycles through 2,1 to point to current layer
358    
359     km1 = MAX(1,k-1)
360     kup = 1+MOD(k+1,2)
361     kDown= 1+MOD(k,2)
362    
363     iMin = 1-OLx
364     iMax = sNx+OLx
365     jMin = 1-OLy
366     jMax = sNy+OLy
367    
368 jmc 1.64 kp1Msk=1.
369     IF (k.EQ.Nr) kp1Msk=0.
370     DO j=1-Oly,sNy+Oly
371     DO i=1-Olx,sNx+Olx
372     rTransKp1(i,j) = kp1Msk*rTrans(i,j)
373     ENDDO
374     ENDDO
375 heimbach 1.66 #ifdef ALLOW_AUTODIFF_TAMC
376     CADJ STORE rTransKp1(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
377     #endif
378 jmc 1.64
379 adcroft 1.1 C-- Get temporary terms used by tendency routines
380     CALL CALC_COMMON_FACTORS (
381     I bi,bj,iMin,iMax,jMin,jMax,k,
382     O xA,yA,uTrans,vTrans,rTrans,maskUp,
383     I myThid)
384 jmc 1.19
385 jmc 1.64 IF (k.EQ.1) THEN
386     C- Surface interface :
387     DO j=1-Oly,sNy+Oly
388     DO i=1-Olx,sNx+Olx
389     rTrans(i,j) = 0.
390     ENDDO
391     ENDDO
392     ELSE
393     C- Interior interface :
394     DO j=1-Oly,sNy+Oly
395     DO i=1-Olx,sNx+Olx
396     rTrans(i,j) = rTrans(i,j)*maskC(i,j,k-1,bi,bj)
397     ENDDO
398     ENDDO
399     ENDIF
400    
401 jmc 1.19 #ifdef ALLOW_GMREDI
402 heimbach 1.35
403 jmc 1.19 C-- Residual transp = Bolus transp + Eulerian transp
404     IF (useGMRedi) THEN
405     CALL GMREDI_CALC_UVFLOW(
406     & uTrans, vTrans, bi, bj, k, myThid)
407     IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
408     & rTrans, bi, bj, k, myThid)
409     ENDIF
410 heimbach 1.35
411 heimbach 1.66 #ifdef ALLOW_AUTODIFF_TAMC
412     CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
413 heimbach 1.35 #ifdef GM_BOLUS_ADVEC
414     CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
415     CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
416     #endif
417     #endif /* ALLOW_AUTODIFF_TAMC */
418    
419 jmc 1.19 #endif /* ALLOW_GMREDI */
420 adcroft 1.1
421     #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
422     C-- Calculate the total vertical diffusivity
423 jmc 1.81 IF ( .NOT.implicitDiffusion ) THEN
424     CALL CALC_DIFFUSIVITY(
425     I bi,bj,iMin,iMax,jMin,jMax,k,
426     I maskUp,
427     O kappaRT,kappaRS,
428     I myThid)
429     ENDIF
430 heimbach 1.52 # ifdef ALLOW_AUTODIFF_TAMC
431 jmc 1.81 CADJ STORE kappaRT(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
432     CADJ STORE kappaRS(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
433 heimbach 1.52 # endif /* ALLOW_AUTODIFF_TAMC */
434 adcroft 1.1 #endif
435    
436     iMin = 1-OLx+2
437     iMax = sNx+OLx-1
438     jMin = 1-OLy+2
439     jMax = sNy+OLy-1
440    
441     C-- Calculate active tracer tendencies (gT,gS,...)
442     C and step forward storing result in gTnm1, gSnm1, etc.
443 stephd 1.75 #ifndef ALLOW_OFFLINE
444 adcroft 1.1 IF ( tempStepping ) THEN
445     CALL CALC_GT(
446     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
447 jmc 1.64 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
448 jmc 1.81 I kappaRT,
449 adcroft 1.1 U fVerT,
450 adcroft 1.7 I myTime,myIter,myThid)
451 adcroft 1.1 CALL TIMESTEP_TRACER(
452 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
453 adcroft 1.1 I theta, gT,
454     I myIter, myThid)
455     ENDIF
456 stephd 1.75 #endif
457 jmc 1.44
458 stephd 1.75 #ifndef ALLOW_OFFLINE
459 adcroft 1.1 IF ( saltStepping ) THEN
460     CALL CALC_GS(
461     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
462 jmc 1.64 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
463 jmc 1.81 I kappaRS,
464 adcroft 1.1 U fVerS,
465 adcroft 1.7 I myTime,myIter,myThid)
466 adcroft 1.1 CALL TIMESTEP_TRACER(
467 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
468 adcroft 1.1 I salt, gS,
469     I myIter, myThid)
470     ENDIF
471 stephd 1.75 #endif
472 adcroft 1.17 #ifdef ALLOW_PTRACERS
473     IF ( usePTRACERS ) THEN
474 jmc 1.81 IF ( .NOT.implicitDiffusion ) THEN
475     CALL PTRACERS_CALC_DIFF(
476     I bi,bj,iMin,iMax,jMin,jMax,k,
477     I maskUp,
478     O kappaRTr,
479     I myThid)
480     ENDIF
481 heimbach 1.42 CALL PTRACERS_INTEGRATE(
482 adcroft 1.17 I bi,bj,k,
483 jmc 1.64 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
484 jmc 1.80 U fVerP,
485 jmc 1.81 I kappaRTr,
486 adcroft 1.17 I myIter,myTime,myThid)
487     ENDIF
488     #endif /* ALLOW_PTRACERS */
489 adcroft 1.1
490     #ifdef ALLOW_OBCS
491     C-- Apply open boundary conditions
492     IF (useOBCS) THEN
493 adcroft 1.7 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
494 adcroft 1.1 END IF
495     #endif /* ALLOW_OBCS */
496 edhill 1.54
497 jmc 1.59 C-- Freeze water
498     C this bit of code is left here for backward compatibility.
499     C freezing at surface level has been moved to FORWARD_STEP
500 stephd 1.75 #ifndef ALLOW_OFFLINE
501 jmc 1.59 IF ( useOldFreezing .AND. .NOT. useSEAICE
502 jmc 1.61 & .AND. .NOT.(useThSIce.AND.k.EQ.1) ) THEN
503 jmc 1.59 #ifdef ALLOW_AUTODIFF_TAMC
504     CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
505     CADJ & , key = kkey, byte = isbyte
506     #endif /* ALLOW_AUTODIFF_TAMC */
507     CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
508     ENDIF
509 stephd 1.75 #endif
510 adcroft 1.1
511     C-- end of thermodynamic k loop (Nr:1)
512     ENDDO
513 cheisey 1.31
514 jmc 1.81 iMin = 1
515     iMax = sNx
516     jMin = 1
517     jMax = sNy
518 adcroft 1.1
519 jmc 1.63 C-- Implicit vertical advection & diffusion
520 stephd 1.75 #ifndef ALLOW_OFFLINE
521 jmc 1.81 IF ( tempStepping .AND. implicitDiffusion ) THEN
522     CALL CALC_3D_DIFFUSIVITY(
523     I bi,bj,iMin,iMax,jMin,jMax,
524     I GAD_TEMPERATURE, useGMredi, useKPP,
525     O kappaRk,
526     I myThid)
527     ENDIF
528 jmc 1.63 #ifdef INCLUDE_IMPLVERTADV_CODE
529 jmc 1.84 c IF ( tempImplVertAdv ) THEN
530     IF ( (dTtracerLev(1).NE.dTtracerLev(Nr)
531     & .AND. tempStepping .AND. implicitDiffusion)
532     & .OR.tempImplVertAdv ) THEN
533 jmc 1.63 CALL GAD_IMPLICIT_R(
534     I tempImplVertAdv, tempAdvScheme, GAD_TEMPERATURE,
535 jmc 1.81 I kappaRk, wVel, theta,
536 jmc 1.63 U gT,
537     I bi, bj, myTime, myIter, myThid )
538     ELSEIF ( tempStepping .AND. implicitDiffusion ) THEN
539     #else /* INCLUDE_IMPLVERTADV_CODE */
540     IF ( tempStepping .AND. implicitDiffusion ) THEN
541     #endif /* INCLUDE_IMPLVERTADV_CODE */
542 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
543 jmc 1.81 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
544 heimbach 1.30 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
545 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
546 jmc 1.63 CALL IMPLDIFF(
547 adcroft 1.1 I bi, bj, iMin, iMax, jMin, jMax,
548 jmc 1.84 I dTtracerLev(1), kappaRk, recip_hFacC,
549 adcroft 1.7 U gT,
550 adcroft 1.1 I myThid )
551 jmc 1.63 ENDIF
552 jmc 1.81 #endif /* ndef ALLOW_OFFLINE */
553    
554     #ifdef ALLOW_TIMEAVE
555 dimitri 1.85 #ifndef MINIMAL_TAVE_OUTPUT
556 jmc 1.81 useVariableK = useKPP .OR. usePP81 .OR. useMY82 .OR. useGGL90
557     & .OR. useGMredi .OR. ivdc_kappa.NE.0.
558     IF (taveFreq.GT.0. .AND. useVariableK ) THEN
559     IF (implicitDiffusion) THEN
560     CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRk,
561     I Nr, 3, deltaTclock, bi, bj, myThid)
562     c ELSE
563     c CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
564     c I Nr, 3, deltaTclock, bi, bj, myThid)
565     ENDIF
566     ENDIF
567 dimitri 1.85 #endif /* ndef MINIMAL_TAVE_OUTPUT */
568 jmc 1.81 #endif /* ALLOW_TIMEAVE */
569 adcroft 1.1
570 stephd 1.75 #ifndef ALLOW_OFFLINE
571 jmc 1.81 IF ( saltStepping .AND. implicitDiffusion ) THEN
572     CALL CALC_3D_DIFFUSIVITY(
573     I bi,bj,iMin,iMax,jMin,jMax,
574     I GAD_SALINITY, useGMredi, useKPP,
575     O kappaRk,
576     I myThid)
577     ENDIF
578    
579 jmc 1.63 #ifdef INCLUDE_IMPLVERTADV_CODE
580 jmc 1.84 c IF ( saltImplVertAdv ) THEN
581     IF ( (dTtracerLev(1).NE.dTtracerLev(Nr)
582     & .AND. saltStepping .AND. implicitDiffusion)
583     & .OR.saltImplVertAdv ) THEN
584 jmc 1.63 CALL GAD_IMPLICIT_R(
585     I saltImplVertAdv, saltAdvScheme, GAD_SALINITY,
586 jmc 1.81 I kappaRk, wVel, salt,
587 jmc 1.63 U gS,
588     I bi, bj, myTime, myIter, myThid )
589     ELSEIF ( saltStepping .AND. implicitDiffusion ) THEN
590     #else /* INCLUDE_IMPLVERTADV_CODE */
591     IF ( saltStepping .AND. implicitDiffusion ) THEN
592     #endif /* INCLUDE_IMPLVERTADV_CODE */
593 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
594 jmc 1.81 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
595 heimbach 1.30 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
596 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
597 jmc 1.63 CALL IMPLDIFF(
598 adcroft 1.1 I bi, bj, iMin, iMax, jMin, jMax,
599 jmc 1.84 I dTtracerLev(1), kappaRk, recip_hFacC,
600 adcroft 1.7 U gS,
601 adcroft 1.1 I myThid )
602 jmc 1.63 ENDIF
603 stephd 1.75 #endif
604 adcroft 1.1
605 adcroft 1.17 #ifdef ALLOW_PTRACERS
606 jmc 1.86 IF ( usePTRACERS ) THEN
607     C-- Vertical advection/diffusion (implicit) for passive tracers
608     CALL PTRACERS_IMPLICIT(
609     U kappaRk,
610     I bi, bj, myTime, myIter, myThid )
611 jmc 1.63 ENDIF
612 adcroft 1.17 #endif /* ALLOW_PTRACERS */
613    
614 adcroft 1.1 #ifdef ALLOW_OBCS
615     C-- Apply open boundary conditions
616 jmc 1.63 IF ( ( implicitDiffusion
617     & .OR. tempImplVertAdv
618     & .OR. saltImplVertAdv
619     & ) .AND. useOBCS ) THEN
620 adcroft 1.1 DO K=1,Nr
621 adcroft 1.7 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
622 adcroft 1.1 ENDDO
623 jmc 1.63 ENDIF
624 adcroft 1.1 #endif /* ALLOW_OBCS */
625    
626 jmc 1.39 #ifdef ALLOW_TIMEAVE
627 jmc 1.79 IF ( taveFreq.GT. 0. _d 0 .AND. fluidIsWater ) THEN
628 jmc 1.73 CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
629     ENDIF
630 dimitri 1.85 #ifndef MINIMAL_TAVE_OUTPUT
631 jmc 1.39 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
632 jmc 1.70 CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
633 jmc 1.39 I Nr, deltaTclock, bi, bj, myThid)
634     ENDIF
635 dimitri 1.85 #endif /* ndef MINIMAL_TAVE_OUTPUT */
636 jmc 1.39 #endif /* ALLOW_TIMEAVE */
637    
638 heimbach 1.22 #endif /* SINGLE_LAYER_MODE */
639 adcroft 1.1
640 jmc 1.39 C-- end bi,bj loops.
641 adcroft 1.1 ENDDO
642     ENDDO
643 adcroft 1.17
644 edhill 1.56 #ifdef ALLOW_DEBUG
645 adcroft 1.17 If (debugMode) THEN
646     CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
647     CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
648     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
649     CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
650     CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
651     CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
652     CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
653     CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
654     CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
655 adcroft 1.18 #ifdef ALLOW_PTRACERS
656     IF ( usePTRACERS ) THEN
657     CALL PTRACERS_DEBUG(myThid)
658     ENDIF
659     #endif /* ALLOW_PTRACERS */
660 adcroft 1.17 ENDIF
661 adcroft 1.40 #endif
662    
663 edhill 1.56 #ifdef ALLOW_DEBUG
664 heimbach 1.43 IF ( debugLevel .GE. debLevB )
665 jmc 1.63 & CALL DEBUG_LEAVE('THERMODYNAMICS',myThid)
666 adcroft 1.17 #endif
667 adcroft 1.1
668 jmc 1.83 #endif /* ALLOW_GENERIC_ADVDIFF */
669    
670 adcroft 1.1 RETURN
671     END

  ViewVC Help
Powered by ViewVC 1.1.22