/[MITgcm]/MITgcm_contrib/jscott/code_rafmod/thermodynamics.F
ViewVC logotype

Annotation of /MITgcm_contrib/jscott/code_rafmod/thermodynamics.F

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


Revision 1.2 - (hide annotations) (download)
Thu Sep 3 20:40:01 2009 UTC (15 years, 10 months ago) by jscott
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +127 -41 lines
update code for crude ML horiz mixing scheme

1 jscott 1.2 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.125 2009/06/26 23:10:09 jahn Exp $
2 jscott 1.1 C $Name: $
3    
4     #include "PACKAGES_CONFIG.h"
5     #include "CPP_OPTIONS.h"
6     #ifdef ALLOW_GENERIC_ADVDIFF
7     # include "GAD_OPTIONS.h"
8     #endif
9 jscott 1.2 #ifdef ALLOW_LONGSTEP
10     #include "LONGSTEP_OPTIONS.h"
11     #endif
12 jscott 1.1
13     #ifdef ALLOW_AUTODIFF_TAMC
14     # ifdef ALLOW_GMREDI
15     # include "GMREDI_OPTIONS.h"
16     # endif
17     # ifdef ALLOW_KPP
18     # include "KPP_OPTIONS.h"
19     # endif
20     #endif /* ALLOW_AUTODIFF_TAMC */
21    
22     CBOP
23     C !ROUTINE: THERMODYNAMICS
24     C !INTERFACE:
25     SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
26     C !DESCRIPTION: \bv
27     C *==========================================================*
28     C | SUBROUTINE THERMODYNAMICS
29     C | o Controlling routine for the prognostic part of the
30     C | thermo-dynamics.
31     C *===========================================================
32     C | The algorithm...
33     C |
34     C | "Correction Step"
35     C | =================
36     C | Here we update the horizontal velocities with the surface
37     C | pressure such that the resulting flow is either consistent
38     C | with the free-surface evolution or the rigid-lid:
39     C | U[n] = U* + dt x d/dx P
40     C | V[n] = V* + dt x d/dy P
41     C |
42     C | "Calculation of Gs"
43     C | ===================
44     C | This is where all the accelerations and tendencies (ie.
45     C | physics, parameterizations etc...) are calculated
46     C | rho = rho ( theta[n], salt[n] )
47     C | b = b(rho, theta)
48     C | K31 = K31 ( rho )
49     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
50     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
51     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
52     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
53     C |
54     C | "Time-stepping" or "Prediction"
55     C | ================================
56     C | The models variables are stepped forward with the appropriate
57     C | time-stepping scheme (currently we use Adams-Bashforth II)
58     C | - For momentum, the result is always *only* a "prediction"
59     C | in that the flow may be divergent and will be "corrected"
60     C | later with a surface pressure gradient.
61     C | - Normally for tracers the result is the new field at time
62     C | level [n+1} *BUT* in the case of implicit diffusion the result
63     C | is also *only* a prediction.
64     C | - We denote "predictors" with an asterisk (*).
65     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
66     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
67     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
68     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69     C | With implicit diffusion:
70     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
71     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
72     C | (1 + dt * K * d_zz) theta[n] = theta*
73     C | (1 + dt * K * d_zz) salt[n] = salt*
74     C |
75     C *==========================================================*
76     C \ev
77    
78     C !USES:
79     IMPLICIT NONE
80     C == Global variables ===
81     #include "SIZE.h"
82     #include "EEPARAMS.h"
83     #include "PARAMS.h"
84 jscott 1.2 #include "RESTART.h"
85 jscott 1.1 #include "DYNVARS.h"
86     #include "GRID.h"
87     #ifdef ALLOW_GENERIC_ADVDIFF
88     # include "GAD.h"
89     # include "GAD_SOM_VARS.h"
90     #endif
91     #ifdef ALLOW_PTRACERS
92 jscott 1.2 #ifndef ALLOW_LONGSTEP
93 jscott 1.1 # include "PTRACERS_SIZE.h"
94 jscott 1.2 # include "PTRACERS_PARAMS.h"
95     # include "PTRACERS_FIELDS.h"
96     #endif
97 jscott 1.1 #endif
98     #ifdef ALLOW_TIMEAVE
99     # include "TIMEAVE_STATV.h"
100     #endif
101    
102     #ifdef ALLOW_AUTODIFF_TAMC
103     # include "tamc.h"
104     # include "tamc_keys.h"
105     # include "FFIELDS.h"
106 jscott 1.2 # include "SURFACE.h"
107 jscott 1.1 # include "EOS.h"
108     # ifdef ALLOW_KPP
109     # include "KPP.h"
110     # endif
111     # ifdef ALLOW_GMREDI
112     # include "GMREDI.h"
113     # endif
114     # ifdef ALLOW_EBM
115     # include "EBM.h"
116     # endif
117     #endif /* ALLOW_AUTODIFF_TAMC */
118    
119     C !INPUT/OUTPUT PARAMETERS:
120     C == Routine arguments ==
121     C myTime - Current time in simulation
122     C myIter - Current iteration number in simulation
123     C myThid - Thread number for this instance of the routine.
124     _RL myTime
125     INTEGER myIter
126     INTEGER myThid
127    
128     #ifdef ALLOW_GENERIC_ADVDIFF
129     C !LOCAL VARIABLES:
130     C == Local variables
131     C xA, yA - Per block temporaries holding face areas
132     C uFld, vFld, wFld - Local copy of velocity field (3 components)
133     C uTrans, vTrans, rTrans - Per block temporaries holding flow transport
134     C o uTrans: Zonal transport
135     C o vTrans: Meridional transport
136     C o rTrans: Vertical transport
137     C rTransKp1 o vertical volume transp. at interface k+1
138     C maskUp o maskUp: land/water mask for W points
139     C fVer[STUV] o fVer: Vertical flux term - note fVer
140     C is "pipelined" in the vertical
141     C so we need an fVer for each
142     C variable.
143     C kappaRT, - Total diffusion in vertical at level k, for T and S
144     C kappaRS (background + spatially varying, isopycnal term).
145 jscott 1.2 C kappaRTr - Total diffusion in vertical at level k,
146 jscott 1.1 C for each passive Tracer
147     C kappaRk - Total diffusion in vertical, all levels, 1 tracer
148     C useVariableK = T when vertical diffusion is not constant
149     C iMin, iMax - Ranges and sub-block indices on which calculations
150     C jMin, jMax are applied.
151     C bi, bj
152     C k, kup, - Index for layer above and below. kup and kDown
153 jscott 1.2 C kDown, km1 are switched with layer to be the appropriate
154 jscott 1.1 C index into fVerTerm.
155     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
156     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
157     _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
158     _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
159     _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
160     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
161     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
162     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
163     _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
164     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
165     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
166     _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
167     _RL kappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
168     _RL kappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
169     #ifdef ALLOW_PTRACERS
170 jscott 1.2 #ifndef ALLOW_LONGSTEP
171 jscott 1.1 _RL fVerP (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
172     _RL kappaRTr(1-Olx:sNx+Olx,1-Oly:sNy+Oly,PTRACERS_num)
173     #endif
174 jscott 1.2 #endif
175 jscott 1.1 _RL kappaRk (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
176     _RL diffKh3d_x (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
177     _RL diffKh3d_y (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
178 jscott 1.2 INTEGER iMin, iMax
179 jscott 1.1 INTEGER jMin, jMax
180     INTEGER bi, bj
181     INTEGER i, j
182     INTEGER k, km1, kup, kDown
183     #ifdef ALLOW_ADAMSBASHFORTH_3
184     INTEGER iterNb, m1, m2
185     #endif
186     #ifdef ALLOW_TIMEAVE
187     LOGICAL useVariableK
188     #endif
189     #ifdef ALLOW_PTRACERS
190 jscott 1.2 #ifndef ALLOW_LONGSTEP
191 jscott 1.1 INTEGER iTracer, ip
192     #endif
193 jscott 1.2 #endif
194 jscott 1.1
195     CEOP
196    
197     #ifdef ALLOW_DEBUG
198     IF ( debugLevel .GE. debLevB )
199     & CALL DEBUG_ENTER('THERMODYNAMICS',myThid)
200     #endif
201    
202     #ifdef ALLOW_AUTODIFF_TAMC
203     C-- dummy statement to end declaration part
204     ikey = 1
205     itdkey = 1
206     #endif /* ALLOW_AUTODIFF_TAMC */
207    
208     #ifdef ALLOW_AUTODIFF_TAMC
209     C-- HPF directive to help TAMC
210     CHPF$ INDEPENDENT
211     #endif /* ALLOW_AUTODIFF_TAMC */
212    
213     C-- Compute correction at the surface for Lin Free Surf.
214 jscott 1.2 #ifdef ALLOW_AUTODIFF_TAMC
215     TsurfCor = 0. _d 0
216     SsurfCor = 0. _d 0
217     #endif
218     IF (linFSConserveTr) THEN
219     #ifdef ALLOW_AUTODIFF_TAMC
220     CADJ STORE theta,salt,wvel = comlev1, key = ikey_dynamics, byte=isbyte
221     #endif
222     CALL CALC_WSURF_TR(theta,salt,wVel,
223     & myTime,myIter,myThid)
224     ENDIF
225 jscott 1.1
226     DO bj=myByLo(myThid),myByHi(myThid)
227    
228     #ifdef ALLOW_AUTODIFF_TAMC
229     C-- HPF directive to help TAMC
230     CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
231     CHPF$& ,utrans,vtrans,xA,yA
232     CHPF$& ,kappaRT,kappaRS
233     CHPF$& )
234     # ifdef ALLOW_PTRACERS
235 jscott 1.2 # ifndef ALLOW_LONGSTEP
236 jscott 1.1 CHPF$ INDEPENDENT, NEW (fVerP,kappaRTr)
237     # endif
238 jscott 1.2 # endif
239 jscott 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
240    
241     DO bi=myBxLo(myThid),myBxHi(myThid)
242    
243     #ifdef ALLOW_AUTODIFF_TAMC
244     act1 = bi - myBxLo(myThid)
245     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
246     act2 = bj - myByLo(myThid)
247     max2 = myByHi(myThid) - myByLo(myThid) + 1
248     act3 = myThid - 1
249     max3 = nTx*nTy
250     act4 = ikey_dynamics - 1
251     itdkey = (act1 + 1) + act2*max1
252     & + act3*max1*max2
253     & + act4*max1*max2*max3
254     #endif /* ALLOW_AUTODIFF_TAMC */
255    
256     C-- Set up work arrays with valid (i.e. not NaN) values
257     C These inital values do not alter the numerical results. They
258     C just ensure that all memory references are to valid floating
259     C point numbers. This prevents spurious hardware signals due to
260     C uninitialised but inert locations.
261    
262     DO j=1-OLy,sNy+OLy
263     DO i=1-OLx,sNx+OLx
264     xA(i,j) = 0. _d 0
265     yA(i,j) = 0. _d 0
266     uTrans(i,j) = 0. _d 0
267     vTrans(i,j) = 0. _d 0
268     rTrans (i,j) = 0. _d 0
269     rTransKp1(i,j) = 0. _d 0
270     fVerT (i,j,1) = 0. _d 0
271     fVerT (i,j,2) = 0. _d 0
272     fVerS (i,j,1) = 0. _d 0
273     fVerS (i,j,2) = 0. _d 0
274     kappaRT(i,j) = 0. _d 0
275     kappaRS(i,j) = 0. _d 0
276     ENDDO
277     ENDDO
278    
279     DO k=1,Nr
280     DO j=1-OLy,sNy+OLy
281     DO i=1-OLx,sNx+OLx
282     C This is currently also used by IVDC and Diagnostics
283     kappaRk(i,j,k) = 0. _d 0
284     C- tracer tendency needs to be set to zero (moved here from gad_calc_rhs):
285     gT(i,j,k,bi,bj) = 0. _d 0
286     gS(i,j,k,bi,bj) = 0. _d 0
287     ENDDO
288     ENDDO
289     ENDDO
290    
291     #ifdef ALLOW_PTRACERS
292 jscott 1.2 #ifndef ALLOW_LONGSTEP
293 jscott 1.1 IF ( usePTRACERS ) THEN
294     DO ip=1,PTRACERS_num
295     DO j=1-OLy,sNy+OLy
296     DO i=1-OLx,sNx+OLx
297     fVerP (i,j,1,ip) = 0. _d 0
298     fVerP (i,j,2,ip) = 0. _d 0
299     kappaRTr(i,j,ip) = 0. _d 0
300     ENDDO
301     ENDDO
302     ENDDO
303     C- set tracer tendency to zero:
304     DO iTracer=1,PTRACERS_num
305     DO k=1,Nr
306     DO j=1-OLy,sNy+OLy
307     DO i=1-OLx,sNx+OLx
308     gPTr(i,j,k,bi,bj,itracer) = 0. _d 0
309     ENDDO
310     ENDDO
311     ENDDO
312     ENDDO
313     ENDIF
314     #endif
315 jscott 1.2 #endif
316 jscott 1.1
317     #ifdef ALLOW_ADAMSBASHFORTH_3
318     C- Apply AB on T,S :
319     iterNb = myIter
320     IF (staggerTimeStep) iterNb = myIter - 1
321     m1 = 1 + MOD(iterNb+1,2)
322     m2 = 1 + MOD( iterNb ,2)
323     C compute T^n+1/2 (stored in gtNm) extrapolating T forward in time
324     IF ( AdamsBashforth_T ) CALL ADAMS_BASHFORTH3(
325     I bi, bj, 0,
326     U theta, gtNm,
327     I tempStartAB, iterNb, myThid )
328     C compute S^n+1/2 (stored in gsNm) extrapolating S forward in time
329     IF ( AdamsBashforth_S ) CALL ADAMS_BASHFORTH3(
330     I bi, bj, 0,
331     U salt, gsNm,
332     I saltStartAB, iterNb, myThid )
333     #endif /* ALLOW_ADAMSBASHFORTH_3 */
334    
335     c iMin = 1-OLx
336     c iMax = sNx+OLx
337     c jMin = 1-OLy
338     c jMax = sNy+OLy
339    
340     #ifdef ALLOW_AUTODIFF_TAMC
341     cph avoids recomputation of integrate_for_w
342     CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
343     #endif /* ALLOW_AUTODIFF_TAMC */
344    
345 jscott 1.2 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
346 jscott 1.1 C-- MOST of THERMODYNAMICS will be disabled
347     #ifndef SINGLE_LAYER_MODE
348    
349     #ifdef ALLOW_AUTODIFF_TAMC
350     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
351     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
352     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
353     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
354     # if ((defined ALLOW_DEPTH_CONTROL) || (defined NONLIN_FRSURF))
355     CADJ STORE gtnm1(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
356     CADJ STORE gsnm1(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
357     # endif
358     #endif /* ALLOW_AUTODIFF_TAMC */
359    
360     #ifndef DISABLE_MULTIDIM_ADVECTION
361     C-- Some advection schemes are better calculated using a multi-dimensional
362     C method in the absence of any other terms and, if used, is done here.
363     C
364     C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
365     C The default is to use multi-dimensinal advection for non-linear advection
366     C schemes. However, for the sake of efficiency of the adjoint it is necessary
367     C to be able to exclude this scheme to avoid excessive storage and
368     C recomputation. It *is* differentiable, if you need it.
369     C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
370     C disable this section of code.
371 jscott 1.2 #ifdef GAD_ALLOW_TS_SOM_ADV
372 jscott 1.1 IF ( tempSOM_Advection ) THEN
373     #ifdef ALLOW_DEBUG
374     IF ( debugLevel .GE. debLevB )
375     & CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid)
376     #endif
377     CALL GAD_SOM_ADVECT(
378     I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme,
379 jscott 1.2 I GAD_TEMPERATURE, dTtracerLev,
380 jscott 1.1 I uVel, vVel, wVel, theta,
381     U som_T,
382     O gT,
383     I bi,bj,myTime,myIter,myThid)
384     ELSEIF (tempMultiDimAdvec) THEN
385 jscott 1.2 #else /* GAD_ALLOW_TS_SOM_ADV */
386 jscott 1.1 IF (tempMultiDimAdvec) THEN
387 jscott 1.2 #endif /* GAD_ALLOW_TS_SOM_ADV */
388 jscott 1.1 #ifdef ALLOW_DEBUG
389     IF ( debugLevel .GE. debLevB )
390     & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
391     #endif
392     CALL GAD_ADVECTION(
393     I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme,
394 jscott 1.2 I GAD_TEMPERATURE, dTtracerLev,
395 jscott 1.1 I uVel, vVel, wVel, theta,
396     O gT,
397     I bi,bj,myTime,myIter,myThid)
398     ENDIF
399 jscott 1.2 #ifdef GAD_ALLOW_TS_SOM_ADV
400 jscott 1.1 IF ( saltSOM_Advection ) THEN
401     #ifdef ALLOW_DEBUG
402     IF ( debugLevel .GE. debLevB )
403     & CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid)
404     #endif
405     CALL GAD_SOM_ADVECT(
406     I saltImplVertAdv, saltAdvScheme, saltVertAdvScheme,
407 jscott 1.2 I GAD_SALINITY, dTtracerLev,
408 jscott 1.1 I uVel, vVel, wVel, salt,
409     U som_S,
410     O gS,
411     I bi,bj,myTime,myIter,myThid)
412     ELSEIF (saltMultiDimAdvec) THEN
413 jscott 1.2 #else /* GAD_ALLOW_TS_SOM_ADV */
414 jscott 1.1 IF (saltMultiDimAdvec) THEN
415 jscott 1.2 #endif /* GAD_ALLOW_TS_SOM_ADV */
416 jscott 1.1 #ifdef ALLOW_DEBUG
417     IF ( debugLevel .GE. debLevB )
418     & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
419     #endif
420     CALL GAD_ADVECTION(
421     I saltImplVertAdv, saltAdvScheme, saltVertAdvScheme,
422 jscott 1.2 I GAD_SALINITY, dTtracerLev,
423 jscott 1.1 I uVel, vVel, wVel, salt,
424     O gS,
425     I bi,bj,myTime,myIter,myThid)
426     ENDIF
427    
428     C Since passive tracers are configurable separately from T,S we
429     C call the multi-dimensional method for PTRACERS regardless
430     C of whether multiDimAdvection is set or not.
431     #ifdef ALLOW_PTRACERS
432 jscott 1.2 #ifndef ALLOW_LONGSTEP
433 jscott 1.1 IF ( usePTRACERS ) THEN
434     #ifdef ALLOW_DEBUG
435     IF ( debugLevel .GE. debLevB )
436     & CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
437     #endif
438     CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
439     ENDIF
440 jscott 1.2 #endif /* ALLOW_LONGSTEP */
441 jscott 1.1 #endif /* ALLOW_PTRACERS */
442     #endif /* DISABLE_MULTIDIM_ADVECTION */
443    
444     #ifdef ALLOW_DEBUG
445     IF ( debugLevel .GE. debLevB )
446     & CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
447     #endif
448    
449     Cjrs
450 jscott 1.2 IF (diffKhT.NE.0) THEN
451 jscott 1.1 CALL CALC_MLD(bi,bj,diffkh3d_x,diffkh3d_y,myThid)
452     ENDIF
453    
454 jscott 1.2
455 jscott 1.1 C-- Start of thermodynamics loop
456     DO k=Nr,1,-1
457 jscott 1.2 #ifdef ALLOW_AUTODIFF_TAMC
458 jscott 1.1 C? Patrick Is this formula correct?
459     cph Yes, but I rewrote it.
460     cph Also, the kappaR? need the index and subscript k!
461     kkey = (itdkey-1)*Nr + k
462     #endif /* ALLOW_AUTODIFF_TAMC */
463    
464     C-- km1 Points to level above k (=k-1)
465     C-- kup Cycles through 1,2 to point to layer above
466     C-- kDown Cycles through 2,1 to point to current layer
467    
468     km1 = MAX(1,k-1)
469     kup = 1+MOD(k+1,2)
470     kDown= 1+MOD(k,2)
471    
472     iMin = 1-OLx
473     iMax = sNx+OLx
474     jMin = 1-OLy
475     jMax = sNy+OLy
476    
477     IF (k.EQ.Nr) THEN
478     DO j=1-Oly,sNy+Oly
479     DO i=1-Olx,sNx+Olx
480     rTransKp1(i,j) = 0. _d 0
481     ENDDO
482     ENDDO
483     ELSE
484     DO j=1-Oly,sNy+Oly
485     DO i=1-Olx,sNx+Olx
486     rTransKp1(i,j) = rTrans(i,j)
487     ENDDO
488     ENDDO
489     ENDIF
490     #ifdef ALLOW_AUTODIFF_TAMC
491     CADJ STORE rTransKp1(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
492     #endif
493    
494     C-- Get temporary terms used by tendency routines :
495     C- Calculate horizontal "volume transport" through tracer cell face
496     C anelastic: uTrans,vTrans are scaled by rhoFacC (~ mass transport)
497     CALL CALC_COMMON_FACTORS (
498     I uVel, vVel,
499     O uFld, vFld, uTrans, vTrans, xA, yA,
500     I k,bi,bj, myThid )
501    
502     C- Calculate vertical "volume transport" through tracer cell face
503     IF (k.EQ.1) THEN
504     C- Surface interface :
505     DO j=1-Oly,sNy+Oly
506     DO i=1-Olx,sNx+Olx
507     wFld(i,j) = 0. _d 0
508     maskUp(i,j) = 0. _d 0
509     rTrans(i,j) = 0. _d 0
510     ENDDO
511     ENDDO
512     ELSE
513     C- Interior interface :
514     C anelastic: rTrans is scaled by rhoFacF (~ mass transport)
515     DO j=1-Oly,sNy+Oly
516     DO i=1-Olx,sNx+Olx
517     wFld(i,j) = wVel(i,j,k,bi,bj)
518     maskUp(i,j) = maskC(i,j,k-1,bi,bj)*maskC(i,j,k,bi,bj)
519     rTrans(i,j) = wFld(i,j)*rA(i,j,bi,bj)*maskUp(i,j)
520     & *deepFac2F(k)*rhoFacF(k)
521     ENDDO
522     ENDDO
523     ENDIF
524    
525     #ifdef ALLOW_GMREDI
526     C-- Residual transp = Bolus transp + Eulerian transp
527     IF (useGMRedi) THEN
528     CALL GMREDI_CALC_UVFLOW(
529     U uFld, vFld, uTrans, vTrans,
530     I k, bi, bj, myThid )
531     IF (K.GE.2) THEN
532     CALL GMREDI_CALC_WFLOW(
533     U wFld, rTrans,
534     I k, bi, bj, myThid )
535     ENDIF
536     ENDIF
537     # ifdef ALLOW_AUTODIFF_TAMC
538 jscott 1.2 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
539     CADJ STORE wfld(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
540 jscott 1.1 # ifdef GM_BOLUS_ADVEC
541     CADJ STORE ufld(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
542     CADJ STORE vfld(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
543     CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
544     CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
545     # endif
546     # endif /* ALLOW_AUTODIFF_TAMC */
547     #endif /* ALLOW_GMREDI */
548    
549     #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
550     C-- Calculate the total vertical diffusivity
551     IF ( .NOT.implicitDiffusion ) THEN
552     CALL CALC_DIFFUSIVITY(
553     I bi,bj,iMin,iMax,jMin,jMax,k,
554     I maskUp,
555     O kappaRT,kappaRS,
556     I myThid)
557     ENDIF
558 jscott 1.2 # ifdef ALLOW_AUTODIFF_TAMC
559 jscott 1.1 CADJ STORE kappaRT(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
560     CADJ STORE kappaRS(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
561     # endif /* ALLOW_AUTODIFF_TAMC */
562     #endif
563    
564     iMin = 1-OLx+2
565     iMax = sNx+OLx-1
566     jMin = 1-OLy+2
567     jMax = sNy+OLy-1
568    
569     C-- Calculate active tracer tendencies (gT,gS,...)
570     C and step forward storing result in gT, gS, etc.
571     C--
572 jscott 1.2 # ifdef ALLOW_AUTODIFF_TAMC
573 jscott 1.1 # if ((defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)) && (defined ALLOW_GMREDI)
574 jscott 1.2 # ifdef GM_NON_UNITY_DIAGONAL
575 jscott 1.1 CADJ STORE kux(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
576     CADJ STORE kvy(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
577 jscott 1.2 # endif
578 jscott 1.1 # ifdef GM_EXTRA_DIAGONAL
579     CADJ STORE kuz(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
580     CADJ STORE kvz(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
581     # endif
582     # endif
583     # endif /* ALLOW_AUTODIFF_TAMC */
584     C
585 jscott 1.2 #ifdef ALLOW_AUTODIFF_TAMC
586 jscott 1.1 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
587     cph-test
588     CADJ STORE uFld(:,:), vFld(:,:), wFld(:,:)
589     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
590     CADJ STORE uTrans(:,:), vTrans(:,:)
591     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
592     CADJ STORE xA(:,:), yA(:,:)
593     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
594     # endif
595     #endif /* ALLOW_AUTODIFF_TAMC */
596     C
597     IF ( tempStepping ) THEN
598 jscott 1.2 #ifdef ALLOW_AUTODIFF_TAMC
599 jscott 1.1 CADJ STORE gTnm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
600     # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
601     CADJ STORE gt(:,:,:,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
602     CADJ STORE fvert(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
603     # endif
604     #endif /* ALLOW_AUTODIFF_TAMC */
605     CALL CALC_GT(
606     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
607     I xA, yA, maskUp, uFld, vFld, wFld,
608     I uTrans, vTrans, rTrans, rTransKp1,
609     I kappaRT,diffKh3d_x,diffKh3d_y,
610     U fVerT,
611     I myTime,myIter,myThid)
612     #ifdef ALLOW_ADAMSBASHFORTH_3
613     IF ( AdamsBashforth_T ) THEN
614     CALL TIMESTEP_TRACER(
615 jscott 1.2 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,dTtracerLev(k),
616 jscott 1.1 I gtNm(1-Olx,1-Oly,1,1,1,m2),
617     U gT,
618     I myIter, myThid)
619     ELSE
620     #endif
621     CALL TIMESTEP_TRACER(
622 jscott 1.2 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,dTtracerLev(k),
623 jscott 1.1 I theta,
624     U gT,
625     I myIter, myThid)
626     #ifdef ALLOW_ADAMSBASHFORTH_3
627     ENDIF
628     #endif
629     ENDIF
630    
631     IF ( saltStepping ) THEN
632     #ifdef ALLOW_AUTODIFF_TAMC
633     CADJ STORE gSnm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
634     # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
635     CADJ STORE gs(:,:,:,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
636     CADJ STORE fvers(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
637     # endif
638     #endif /* ALLOW_AUTODIFF_TAMC */
639     CALL CALC_GS(
640     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
641     I xA, yA, maskUp, uFld, vFld, wFld,
642     I uTrans, vTrans, rTrans, rTransKp1,
643     I kappaRS,diffKh3d_x,diffKh3d_y,
644     U fVerS,
645     I myTime,myIter,myThid)
646     #ifdef ALLOW_ADAMSBASHFORTH_3
647     IF ( AdamsBashforth_S ) THEN
648     CALL TIMESTEP_TRACER(
649 jscott 1.2 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,dTtracerLev(k),
650 jscott 1.1 I gsNm(1-Olx,1-Oly,1,1,1,m2),
651     U gS,
652     I myIter, myThid)
653     ELSE
654     #endif
655     CALL TIMESTEP_TRACER(
656 jscott 1.2 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,dTtracerLev(k),
657 jscott 1.1 I salt,
658     U gS,
659     I myIter, myThid)
660     #ifdef ALLOW_ADAMSBASHFORTH_3
661     ENDIF
662     #endif
663     ENDIF
664    
665     #ifdef ALLOW_PTRACERS
666 jscott 1.2 #ifndef ALLOW_LONGSTEP
667 jscott 1.1 IF ( usePTRACERS ) THEN
668     IF ( .NOT.implicitDiffusion ) THEN
669     CALL PTRACERS_CALC_DIFF(
670     I bi,bj,iMin,iMax,jMin,jMax,k,
671     I maskUp,
672     O kappaRTr,
673     I myThid)
674     ENDIF
675 jscott 1.2 # ifdef ALLOW_AUTODIFF_TAMC
676 jscott 1.1 CADJ STORE kappaRTr(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
677     # endif /* ALLOW_AUTODIFF_TAMC */
678     CALL PTRACERS_INTEGRATE(
679     I bi,bj,k,
680     I xA, yA, maskUp, uFld, vFld, wFld,
681     I uTrans, vTrans, rTrans, rTransKp1,
682 jscott 1.2 I kappaRTr,diffKh3d_x,diffKh3d_y,
683 jscott 1.1 U fVerP,
684     I myTime,myIter,myThid)
685     ENDIF
686 jscott 1.2 #endif /*ALLOW_LONGSTEP */
687 jscott 1.1 #endif /* ALLOW_PTRACERS */
688    
689     #ifdef ALLOW_OBCS
690     C-- Apply open boundary conditions
691     IF (useOBCS) THEN
692     CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
693     END IF
694     #endif /* ALLOW_OBCS */
695    
696     C-- Freeze water
697     C this bit of code is left here for backward compatibility.
698 jscott 1.2 C freezing at surface level has been moved to FORWARD_STEP
699 jscott 1.1 IF ( useOldFreezing .AND. .NOT. useSEAICE
700     & .AND. .NOT.(useThSIce.AND.k.EQ.1) ) THEN
701     #ifdef ALLOW_AUTODIFF_TAMC
702     CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
703     CADJ & , key = kkey, byte = isbyte
704     #endif /* ALLOW_AUTODIFF_TAMC */
705     CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
706     ENDIF
707    
708     C-- end of thermodynamic k loop (Nr:1)
709     ENDDO
710    
711 jscott 1.2 #ifdef ALLOW_DOWN_SLOPE
712     IF ( tempStepping .AND. useDOWN_SLOPE ) THEN
713     IF ( usingPCoords ) THEN
714     CALL DWNSLP_APPLY(
715     I GAD_TEMPERATURE, bi, bj, kSurfC,
716     I recip_drF, recip_hFacC, recip_rA,
717     I dTtracerLev,
718     I theta,
719     U gT,
720     I myTime, myIter, myThid )
721     ELSE
722     CALL DWNSLP_APPLY(
723     I GAD_TEMPERATURE, bi, bj, kLowC,
724     I recip_drF, recip_hFacC, recip_rA,
725     I dTtracerLev,
726     I theta,
727     U gT,
728     I myTime, myIter, myThid )
729     ENDIF
730     ENDIF
731     IF ( saltStepping .AND. useDOWN_SLOPE ) THEN
732     IF ( usingPCoords ) THEN
733     CALL DWNSLP_APPLY(
734     I GAD_SALINITY, bi, bj, kSurfC,
735     I recip_drF, recip_hFacC, recip_rA,
736     I dTtracerLev,
737     I salt,
738     U gS,
739     I myTime, myIter, myThid )
740     ELSE
741     CALL DWNSLP_APPLY(
742     I GAD_SALINITY, bi, bj, kLowC,
743     I recip_drF, recip_hFacC, recip_rA,
744     I dTtracerLev,
745     I salt,
746     U gS,
747     I myTime, myIter, myThid )
748     ENDIF
749     ENDIF
750     #ifdef ALLOW_PTRACERS
751     #ifndef ALLOW_LONGSTEP
752     IF ( usePTRACERS .AND. useDOWN_SLOPE ) THEN
753     CALL PTRACERS_DWNSLP_APPLY(
754     I bi, bj, myTime, myIter, myThid )
755     ENDIF
756     #endif /*ALLOW_LONGSTEP */
757     #endif /* ALLOW_PTRACERS */
758     #endif /* ALLOW_DOWN_SLOPE */
759    
760     C All explicit advection/diffusion/sources should now be
761     C done. The updated tracer field is in gPtr. Accumalate
762     C explicit tendency and also reset gPtr to initial tracer
763 jscott 1.1 C field for implicit matrix calculation
764    
765     #ifdef ALLOW_MATRIX
766 jscott 1.2 IF (useMATRIX)
767 jscott 1.1 & CALL MATRIX_STORE_TENDENCY_EXP(bi,bj, myTime,myIter,myThid)
768 jscott 1.2 #endif
769 jscott 1.1
770     iMin = 1
771     iMax = sNx
772     jMin = 1
773     jMax = sNy
774    
775 jscott 1.2 C-- Implicit vertical advection & diffusion
776 jscott 1.1 IF ( tempStepping .AND. implicitDiffusion ) THEN
777     CALL CALC_3D_DIFFUSIVITY(
778     I bi,bj,iMin,iMax,jMin,jMax,
779     I GAD_TEMPERATURE, useGMredi, useKPP,
780     O kappaRk,
781     I myThid)
782     ENDIF
783     #ifdef INCLUDE_IMPLVERTADV_CODE
784     IF ( tempImplVertAdv ) THEN
785     CALL GAD_IMPLICIT_R(
786     I tempImplVertAdv, tempAdvScheme, GAD_TEMPERATURE,
787 jscott 1.2 I dTtracerLev,
788 jscott 1.1 I kappaRk, wVel, theta,
789     U gT,
790     I bi, bj, myTime, myIter, myThid )
791     ELSEIF ( tempStepping .AND. implicitDiffusion ) THEN
792     #else /* INCLUDE_IMPLVERTADV_CODE */
793     IF ( tempStepping .AND. implicitDiffusion ) THEN
794     #endif /* INCLUDE_IMPLVERTADV_CODE */
795     #ifdef ALLOW_AUTODIFF_TAMC
796     CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
797     CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
798     #endif /* ALLOW_AUTODIFF_TAMC */
799     CALL IMPLDIFF(
800     I bi, bj, iMin, iMax, jMin, jMax,
801     I GAD_TEMPERATURE, kappaRk, recip_hFacC,
802     U gT,
803     I myThid )
804     ENDIF
805    
806     #ifdef ALLOW_TIMEAVE
807     useVariableK = useKPP .OR. usePP81 .OR. useMY82 .OR. useGGL90
808     & .OR. useGMredi .OR. ivdc_kappa.NE.0.
809     IF (taveFreq.GT.0. .AND. useVariableK ) THEN
810     IF (implicitDiffusion) THEN
811     CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRk,
812     I Nr, 3, deltaTclock, bi, bj, myThid)
813     c ELSE
814     c CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
815     c I Nr, 3, deltaTclock, bi, bj, myThid)
816     ENDIF
817     ENDIF
818 jscott 1.2 #endif /* ALLOW_TIMEAVE */
819 jscott 1.1
820     IF ( saltStepping .AND. implicitDiffusion ) THEN
821     CALL CALC_3D_DIFFUSIVITY(
822     I bi,bj,iMin,iMax,jMin,jMax,
823     I GAD_SALINITY, useGMredi, useKPP,
824     O kappaRk,
825     I myThid)
826     ENDIF
827    
828     #ifdef INCLUDE_IMPLVERTADV_CODE
829     IF ( saltImplVertAdv ) THEN
830     CALL GAD_IMPLICIT_R(
831     I saltImplVertAdv, saltAdvScheme, GAD_SALINITY,
832 jscott 1.2 I dTtracerLev,
833 jscott 1.1 I kappaRk, wVel, salt,
834     U gS,
835     I bi, bj, myTime, myIter, myThid )
836     ELSEIF ( saltStepping .AND. implicitDiffusion ) THEN
837     #else /* INCLUDE_IMPLVERTADV_CODE */
838     IF ( saltStepping .AND. implicitDiffusion ) THEN
839     #endif /* INCLUDE_IMPLVERTADV_CODE */
840     #ifdef ALLOW_AUTODIFF_TAMC
841     CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
842     CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
843     #endif /* ALLOW_AUTODIFF_TAMC */
844     CALL IMPLDIFF(
845     I bi, bj, iMin, iMax, jMin, jMax,
846     I GAD_SALINITY, kappaRk, recip_hFacC,
847     U gS,
848     I myThid )
849     ENDIF
850    
851     #ifdef ALLOW_PTRACERS
852 jscott 1.2 #ifndef ALLOW_LONGSTEP
853 jscott 1.1 IF ( usePTRACERS ) THEN
854     C-- Vertical advection/diffusion (implicit) for passive tracers
855     CALL PTRACERS_IMPLICIT(
856     U kappaRk,
857     I bi, bj, myTime, myIter, myThid )
858     ENDIF
859 jscott 1.2 #endif /* ALLOW_LONGSTEP */
860 jscott 1.1 #endif /* ALLOW_PTRACERS */
861    
862     #ifdef ALLOW_OBCS
863     C-- Apply open boundary conditions
864     IF ( ( implicitDiffusion
865     & .OR. tempImplVertAdv
866     & .OR. saltImplVertAdv
867     & ) .AND. useOBCS ) THEN
868     DO K=1,Nr
869     CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
870     ENDDO
871     ENDIF
872     #endif /* ALLOW_OBCS */
873    
874     #endif /* SINGLE_LAYER_MODE */
875    
876     C-- end bi,bj loops.
877     ENDDO
878     ENDDO
879    
880     #ifdef ALLOW_DEBUG
881     If (debugMode) THEN
882     CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
883     CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
884     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
885     CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
886     CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
887     CALL DEBUG_STATS_RL(Nr,gT,'Gt (THERMODYNAMICS)',myThid)
888     CALL DEBUG_STATS_RL(Nr,gS,'Gs (THERMODYNAMICS)',myThid)
889     #ifndef ALLOW_ADAMSBASHFORTH_3
890     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
891     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
892     #endif
893     #ifdef ALLOW_PTRACERS
894 jscott 1.2 #ifndef ALLOW_LONGSTEP
895 jscott 1.1 IF ( usePTRACERS ) THEN
896     CALL PTRACERS_DEBUG(myThid)
897     ENDIF
898 jscott 1.2 #endif /* ALLOW_LONGSTEP */
899 jscott 1.1 #endif /* ALLOW_PTRACERS */
900     ENDIF
901     #endif /* ALLOW_DEBUG */
902    
903     #ifdef ALLOW_DEBUG
904     IF ( debugLevel .GE. debLevB )
905     & CALL DEBUG_LEAVE('THERMODYNAMICS',myThid)
906     #endif
907    
908     #endif /* ALLOW_GENERIC_ADVDIFF */
909    
910     RETURN
911     END

  ViewVC Help
Powered by ViewVC 1.1.22