/[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.39 - (hide annotations) (download)
Thu May 1 22:30:33 2003 UTC (21 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint50c_post, checkpoint50d_post, checkpoint50d_pre
Changes since 1.38: +21 -25 lines
fix bug in convective adj. diagnostic (ConvectCount is a local array !)
and add time-average diagnostic of U*V and diff_Kr*theta

1 jmc 1.39 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.38 2003/02/28 02:20:52 heimbach Exp $
2 adcroft 1.1 C $Name: $
3    
4     #include "CPP_OPTIONS.h"
5 jmc 1.21 #ifdef ALLOW_AUTODIFF_TAMC
6     # ifdef ALLOW_GMREDI
7     # include "GMREDI_OPTIONS.h"
8     # endif
9     # ifdef ALLOW_KPP
10     # include "KPP_OPTIONS.h"
11     # endif
12 cheisey 1.31 cswdice --- add ----
13 cheisey 1.32 #ifdef ALLOW_THERM_SEAICE
14 cheisey 1.31 #include "ICE.h"
15     #endif
16     cswdice ------
17 jmc 1.21 #endif /* ALLOW_AUTODIFF_TAMC */
18 adcroft 1.1
19 cnh 1.9 CBOP
20     C !ROUTINE: THERMODYNAMICS
21     C !INTERFACE:
22 adcroft 1.1 SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
23 cnh 1.9 C !DESCRIPTION: \bv
24     C *==========================================================*
25     C | SUBROUTINE THERMODYNAMICS
26     C | o Controlling routine for the prognostic part of the
27     C | thermo-dynamics.
28     C *===========================================================
29     C | The algorithm...
30     C |
31     C | "Correction Step"
32     C | =================
33     C | Here we update the horizontal velocities with the surface
34     C | pressure such that the resulting flow is either consistent
35     C | with the free-surface evolution or the rigid-lid:
36     C | U[n] = U* + dt x d/dx P
37     C | V[n] = V* + dt x d/dy P
38     C |
39     C | "Calculation of Gs"
40     C | ===================
41     C | This is where all the accelerations and tendencies (ie.
42     C | physics, parameterizations etc...) are calculated
43     C | rho = rho ( theta[n], salt[n] )
44     C | b = b(rho, theta)
45     C | K31 = K31 ( rho )
46     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
47     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
48     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
49     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
50     C |
51     C | "Time-stepping" or "Prediction"
52     C | ================================
53     C | The models variables are stepped forward with the appropriate
54     C | time-stepping scheme (currently we use Adams-Bashforth II)
55     C | - For momentum, the result is always *only* a "prediction"
56     C | in that the flow may be divergent and will be "corrected"
57     C | later with a surface pressure gradient.
58     C | - Normally for tracers the result is the new field at time
59     C | level [n+1} *BUT* in the case of implicit diffusion the result
60     C | is also *only* a prediction.
61     C | - We denote "predictors" with an asterisk (*).
62     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
63     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
64     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
65     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66     C | With implicit diffusion:
67     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
68     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69     C | (1 + dt * K * d_zz) theta[n] = theta*
70     C | (1 + dt * K * d_zz) salt[n] = salt*
71     C |
72     C *==========================================================*
73     C \ev
74    
75     C !USES:
76 adcroft 1.1 IMPLICIT NONE
77     C == Global variables ===
78     #include "SIZE.h"
79     #include "EEPARAMS.h"
80     #include "PARAMS.h"
81     #include "DYNVARS.h"
82     #include "GRID.h"
83 adcroft 1.4 #include "GAD.h"
84 adcroft 1.1 #ifdef ALLOW_PASSIVE_TRACER
85     #include "TR1.h"
86     #endif
87     #ifdef ALLOW_AUTODIFF_TAMC
88     # include "tamc.h"
89     # include "tamc_keys.h"
90     # include "FFIELDS.h"
91 heimbach 1.30 # include "EOS.h"
92 adcroft 1.1 # ifdef ALLOW_KPP
93     # include "KPP.h"
94     # endif
95     # ifdef ALLOW_GMREDI
96     # include "GMREDI.h"
97     # endif
98     #endif /* ALLOW_AUTODIFF_TAMC */
99     #ifdef ALLOW_TIMEAVE
100     #include "TIMEAVE_STATV.h"
101     #endif
102    
103 cnh 1.9 C !INPUT/OUTPUT PARAMETERS:
104 adcroft 1.1 C == Routine arguments ==
105     C myTime - Current time in simulation
106     C myIter - Current iteration number in simulation
107     C myThid - Thread number for this instance of the routine.
108     _RL myTime
109     INTEGER myIter
110     INTEGER myThid
111    
112 cnh 1.9 C !LOCAL VARIABLES:
113 adcroft 1.1 C == Local variables
114     C xA, yA - Per block temporaries holding face areas
115     C uTrans, vTrans, rTrans - Per block temporaries holding flow
116     C transport
117     C o uTrans: Zonal transport
118     C o vTrans: Meridional transport
119     C o rTrans: Vertical transport
120     C maskUp o maskUp: land/water mask for W points
121     C fVer[STUV] o fVer: Vertical flux term - note fVer
122     C is "pipelined" in the vertical
123     C so we need an fVer for each
124     C variable.
125     C rhoK, rhoKM1 - Density at current level, and level above
126     C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
127     C phiSurfY or geopotentiel (atmos) in X and Y direction
128     C KappaRT, - Total diffusion in vertical for T and S.
129     C KappaRS (background + spatially varying, isopycnal term).
130 jmc 1.39 C useVariableK = T when vertical diffusion is not constant
131 adcroft 1.1 C iMin, iMax - Ranges and sub-block indices on which calculations
132     C jMin, jMax are applied.
133     C bi, bj
134     C k, kup, - Index for layer above and below. kup and kDown
135     C kDown, km1 are switched with layer to be the appropriate
136     C index into fVerTerm.
137     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
138     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
139     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
140     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
141     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
142     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
143     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
144     _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
145     _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
146     _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147     _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148     _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150     _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
151     _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
152     _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
153     _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
154     _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
155 cnh 1.9 C This is currently used by IVDC and Diagnostics
156 adcroft 1.1 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
157 jmc 1.39 LOGICAL useVariableK
158 adcroft 1.1 INTEGER iMin, iMax
159     INTEGER jMin, jMax
160     INTEGER bi, bj
161     INTEGER i, j
162     INTEGER k, km1, kup, kDown
163    
164 cnh 1.9 CEOP
165 adcroft 1.1
166     #ifdef ALLOW_AUTODIFF_TAMC
167     C-- dummy statement to end declaration part
168     ikey = 1
169 heimbach 1.30 itdkey = 1
170 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
171    
172     C-- Set up work arrays with valid (i.e. not NaN) values
173     C These inital values do not alter the numerical results. They
174     C just ensure that all memory references are to valid floating
175     C point numbers. This prevents spurious hardware signals due to
176     C uninitialised but inert locations.
177     DO j=1-OLy,sNy+OLy
178     DO i=1-OLx,sNx+OLx
179     xA(i,j) = 0. _d 0
180     yA(i,j) = 0. _d 0
181     uTrans(i,j) = 0. _d 0
182     vTrans(i,j) = 0. _d 0
183     rhok (i,j) = 0. _d 0
184     phiSurfX(i,j) = 0. _d 0
185     phiSurfY(i,j) = 0. _d 0
186     ENDDO
187     ENDDO
188    
189    
190     #ifdef ALLOW_AUTODIFF_TAMC
191     C-- HPF directive to help TAMC
192     CHPF$ INDEPENDENT
193     #endif /* ALLOW_AUTODIFF_TAMC */
194    
195     DO bj=myByLo(myThid),myByHi(myThid)
196    
197     #ifdef ALLOW_AUTODIFF_TAMC
198     C-- HPF directive to help TAMC
199 heimbach 1.2 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
200 jmc 1.37 CHPF$& ,utrans,vtrans,xA,yA
201 heimbach 1.2 CHPF$& ,KappaRT,KappaRS
202 adcroft 1.1 CHPF$& )
203     #endif /* ALLOW_AUTODIFF_TAMC */
204    
205     DO bi=myBxLo(myThid),myBxHi(myThid)
206    
207     #ifdef ALLOW_AUTODIFF_TAMC
208     act1 = bi - myBxLo(myThid)
209     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
210     act2 = bj - myByLo(myThid)
211     max2 = myByHi(myThid) - myByLo(myThid) + 1
212     act3 = myThid - 1
213     max3 = nTx*nTy
214     act4 = ikey_dynamics - 1
215 heimbach 1.30 itdkey = (act1 + 1) + act2*max1
216 adcroft 1.1 & + act3*max1*max2
217     & + act4*max1*max2*max3
218     #endif /* ALLOW_AUTODIFF_TAMC */
219    
220     C-- Set up work arrays that need valid initial values
221     DO j=1-OLy,sNy+OLy
222     DO i=1-OLx,sNx+OLx
223     rTrans (i,j) = 0. _d 0
224     fVerT (i,j,1) = 0. _d 0
225     fVerT (i,j,2) = 0. _d 0
226     fVerS (i,j,1) = 0. _d 0
227     fVerS (i,j,2) = 0. _d 0
228     fVerTr1(i,j,1) = 0. _d 0
229     fVerTr1(i,j,2) = 0. _d 0
230 heimbach 1.20 rhoKM1 (i,j) = 0. _d 0
231 adcroft 1.1 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 heimbach 1.20 sigmaX(i,j,k) = 0. _d 0
239     sigmaY(i,j,k) = 0. _d 0
240     sigmaR(i,j,k) = 0. _d 0
241 adcroft 1.1 ConvectCount(i,j,k) = 0.
242 heimbach 1.30 KappaRT(i,j,k) = 0. _d 0
243     KappaRS(i,j,k) = 0. _d 0
244 heimbach 1.5 #ifdef ALLOW_AUTODIFF_TAMC
245 heimbach 1.30 cph all the following init. are necessary for TAF
246     cph although some of these are re-initialised later.
247     gT(i,j,k,bi,bj) = 0. _d 0
248     gS(i,j,k,bi,bj) = 0. _d 0
249     # ifdef ALLOW_PASSIVE_TRACER
250 heimbach 1.5 gTr1(i,j,k,bi,bj) = 0. _d 0
251 heimbach 1.30 # endif
252     # ifdef ALLOW_GMREDI
253     Kwx(i,j,k,bi,bj) = 0. _d 0
254     Kwy(i,j,k,bi,bj) = 0. _d 0
255     Kwz(i,j,k,bi,bj) = 0. _d 0
256     # ifdef GM_NON_UNITY_DIAGONAL
257     Kux(i,j,k,bi,bj) = 0. _d 0
258     Kvy(i,j,k,bi,bj) = 0. _d 0
259     # endif
260     # ifdef GM_EXTRA_DIAGONAL
261     Kuz(i,j,k,bi,bj) = 0. _d 0
262     Kvz(i,j,k,bi,bj) = 0. _d 0
263     # endif
264     # ifdef GM_BOLUS_ADVEC
265     GM_PsiX(i,j,k,bi,bj) = 0. _d 0
266     GM_PsiY(i,j,k,bi,bj) = 0. _d 0
267 heimbach 1.36 # endif
268     # ifdef GM_VISBECK_VARIABLE_K
269     VisbeckK(i,j,bi,bj) = 0. _d 0
270 heimbach 1.30 # endif
271     # endif /* ALLOW_GMREDI */
272     #endif /* ALLOW_AUTODIFF_TAMC */
273 adcroft 1.1 ENDDO
274     ENDDO
275     ENDDO
276    
277 jmc 1.21 iMin = 1-OLx
278 adcroft 1.1 iMax = sNx+OLx
279 jmc 1.21 jMin = 1-OLy
280 adcroft 1.1 jMax = sNy+OLy
281    
282    
283     #ifdef ALLOW_AUTODIFF_TAMC
284 heimbach 1.30 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
285     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
286 heimbach 1.38 CADJ STORE totphihyd
287     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
288 heimbach 1.11 #ifdef ALLOW_KPP
289 heimbach 1.30 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
290     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
291 heimbach 1.11 #endif
292 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
293    
294     C-- Start of diagnostic loop
295     DO k=Nr,1,-1
296    
297     #ifdef ALLOW_AUTODIFF_TAMC
298     C? Patrick, is this formula correct now that we change the loop range?
299     C? Do we still need this?
300     cph kkey formula corrected.
301     cph Needed for rhok, rhokm1, in the case useGMREDI.
302 heimbach 1.30 kkey = (itdkey-1)*Nr + k
303 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
304    
305     C-- Integrate continuity vertically for vertical velocity
306 jmc 1.27 c CALL INTEGRATE_FOR_W(
307     c I bi, bj, k, uVel, vVel,
308     c O wVel,
309     c I myThid )
310 adcroft 1.1
311     #ifdef ALLOW_OBCS
312     #ifdef ALLOW_NONHYDROSTATIC
313     C-- Apply OBC to W if in N-H mode
314 jmc 1.27 c IF (useOBCS.AND.nonHydrostatic) THEN
315     c CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
316     c ENDIF
317 adcroft 1.1 #endif /* ALLOW_NONHYDROSTATIC */
318     #endif /* ALLOW_OBCS */
319    
320 heimbach 1.22 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
321     C-- MOST of THERMODYNAMICS will be disabled
322     #ifndef SINGLE_LAYER_MODE
323    
324 adcroft 1.1 C-- Calculate gradients of potential density for isoneutral
325     C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
326     c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
327     IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
328     #ifdef ALLOW_AUTODIFF_TAMC
329     CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
330     CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
331     #endif /* ALLOW_AUTODIFF_TAMC */
332     CALL FIND_RHO(
333 mlosch 1.26 I bi, bj, iMin, iMax, jMin, jMax, k, k,
334 adcroft 1.1 I theta, salt,
335     O rhoK,
336     I myThid )
337 heimbach 1.30
338 adcroft 1.1 IF (k.GT.1) THEN
339     #ifdef ALLOW_AUTODIFF_TAMC
340     CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
341     CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
342     #endif /* ALLOW_AUTODIFF_TAMC */
343     CALL FIND_RHO(
344 mlosch 1.26 I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
345 adcroft 1.1 I theta, salt,
346     O rhoKm1,
347     I myThid )
348     ENDIF
349     CALL GRAD_SIGMA(
350     I bi, bj, iMin, iMax, jMin, jMax, k,
351     I rhoK, rhoKm1, rhoK,
352     O sigmaX, sigmaY, sigmaR,
353     I myThid )
354     ENDIF
355    
356 heimbach 1.30 #ifdef ALLOW_AUTODIFF_TAMC
357     CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
358     CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
359     #endif /* ALLOW_AUTODIFF_TAMC */
360 adcroft 1.1 C-- Implicit Vertical Diffusion for Convection
361     c ==> should use sigmaR !!!
362     IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
363     CALL CALC_IVDC(
364     I bi, bj, iMin, iMax, jMin, jMax, k,
365     I rhoKm1, rhoK,
366     U ConvectCount, KappaRT, KappaRS,
367     I myTime, myIter, myThid)
368     ENDIF
369    
370 heimbach 1.22 #endif /* SINGLE_LAYER_MODE */
371    
372 adcroft 1.1 C-- end of diagnostic k loop (Nr:1)
373     ENDDO
374    
375     #ifdef ALLOW_AUTODIFF_TAMC
376     cph avoids recomputation of integrate_for_w
377 heimbach 1.30 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
378 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
379    
380     #ifdef ALLOW_OBCS
381     C-- Calculate future values on open boundaries
382     IF (useOBCS) THEN
383 jmc 1.16 CALL OBCS_CALC( bi, bj, myTime+deltaT, myIter+1,
384 adcroft 1.1 I uVel, vVel, wVel, theta, salt,
385     I myThid )
386     ENDIF
387     #endif /* ALLOW_OBCS */
388    
389 cheisey 1.31
390     c********************************************
391     cswdice --- add ---
392 cheisey 1.32 #ifdef ALLOW_THERM_SEAICE
393 cheisey 1.31 C-- Determines forcing terms based on external fields
394     c-- including effects from ice
395     CALL ICE_FORCING(
396     I bi, bj, iMin, iMax, jMin, jMax,
397     I myThid )
398     #else
399    
400     cswdice --- end add ---
401    
402 adcroft 1.1 C-- Determines forcing terms based on external fields
403     C relaxation terms, etc.
404 cheisey 1.31 CALL EXTERNAL_FORCING_SURF(
405 adcroft 1.1 I bi, bj, iMin, iMax, jMin, jMax,
406     I myThid )
407 cheisey 1.31 cswdice --- add ----
408     #endif
409     cswdice --- end add ---
410     c******************************************
411    
412    
413    
414    
415 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
416     cph needed for KPP
417     CADJ STORE surfacetendencyU(:,:,bi,bj)
418 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
419 adcroft 1.1 CADJ STORE surfacetendencyV(:,:,bi,bj)
420 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
421 adcroft 1.1 CADJ STORE surfacetendencyS(:,:,bi,bj)
422 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
423 adcroft 1.1 CADJ STORE surfacetendencyT(:,:,bi,bj)
424 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
425 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
426    
427 heimbach 1.22 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
428     C-- MOST of THERMODYNAMICS will be disabled
429     #ifndef SINGLE_LAYER_MODE
430    
431 adcroft 1.1 #ifdef ALLOW_GMREDI
432    
433 heimbach 1.22 #ifdef ALLOW_AUTODIFF_TAMC
434 heimbach 1.34 cph storing here is needed only for one GMREDI_OPTIONS:
435     cph define GM_BOLUS_ADVEC
436     cph but I've avoided the #ifdef for now, in case more things change
437     CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
438     CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
439     CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
440 heimbach 1.22 #endif /* ALLOW_AUTODIFF_TAMC */
441 heimbach 1.35
442 adcroft 1.1 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
443     IF (useGMRedi) THEN
444 jmc 1.15 CALL GMREDI_CALC_TENSOR(
445     I bi, bj, iMin, iMax, jMin, jMax,
446 adcroft 1.1 I sigmaX, sigmaY, sigmaR,
447     I myThid )
448     #ifdef ALLOW_AUTODIFF_TAMC
449     ELSE
450 jmc 1.15 CALL GMREDI_CALC_TENSOR_DUMMY(
451     I bi, bj, iMin, iMax, jMin, jMax,
452 adcroft 1.1 I sigmaX, sigmaY, sigmaR,
453     I myThid )
454     #endif /* ALLOW_AUTODIFF_TAMC */
455     ENDIF
456    
457     #ifdef ALLOW_AUTODIFF_TAMC
458 heimbach 1.30 CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
459     CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
460     CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
461 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
462    
463     #endif /* ALLOW_GMREDI */
464    
465     #ifdef ALLOW_KPP
466     C-- Compute KPP mixing coefficients
467     IF (useKPP) THEN
468     CALL KPP_CALC(
469     I bi, bj, myTime, myThid )
470     #ifdef ALLOW_AUTODIFF_TAMC
471     ELSE
472     CALL KPP_CALC_DUMMY(
473     I bi, bj, myTime, myThid )
474     #endif /* ALLOW_AUTODIFF_TAMC */
475     ENDIF
476    
477     #ifdef ALLOW_AUTODIFF_TAMC
478     CADJ STORE KPPghat (:,:,:,bi,bj)
479     CADJ & , KPPdiffKzT(:,:,:,bi,bj)
480     CADJ & , KPPdiffKzS(:,:,:,bi,bj)
481     CADJ & , KPPfrac (:,: ,bi,bj)
482 heimbach 1.30 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
483 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
484    
485     #endif /* ALLOW_KPP */
486    
487     #ifdef ALLOW_AUTODIFF_TAMC
488 heimbach 1.30 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
489     CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
490     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
491     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
492     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
493     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
494 adcroft 1.1 #ifdef ALLOW_PASSIVE_TRACER
495 heimbach 1.30 CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
496 adcroft 1.1 #endif
497     #endif /* ALLOW_AUTODIFF_TAMC */
498    
499     #ifdef ALLOW_AIM
500     C AIM - atmospheric intermediate model, physics package code.
501     IF ( useAIM ) THEN
502 jmc 1.33 CALL TIMER_START('AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
503     CALL AIM_DO_PHYSICS( bi, bj, myTime, myIter, myThid )
504     CALL TIMER_STOP( 'AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
505 adcroft 1.1 ENDIF
506     #endif /* ALLOW_AIM */
507 jmc 1.14
508 adcroft 1.12 #ifndef DISABLE_MULTIDIM_ADVECTION
509 adcroft 1.4 C-- Some advection schemes are better calculated using a multi-dimensional
510     C method in the absence of any other terms and, if used, is done here.
511 adcroft 1.13 C
512     C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
513     C The default is to use multi-dimensinal advection for non-linear advection
514     C schemes. However, for the sake of efficiency of the adjoint it is necessary
515     C to be able to exclude this scheme to avoid excessive storage and
516     C recomputation. It *is* differentiable, if you need it.
517     C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
518     C disable this section of code.
519 jmc 1.24 IF (tempMultiDimAdvec) THEN
520 heimbach 1.11 CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,
521     U theta,gT,
522 adcroft 1.6 I myTime,myIter,myThid)
523 jmc 1.23 ENDIF
524 jmc 1.24 IF (saltMultiDimAdvec) THEN
525 heimbach 1.11 CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,
526     U salt,gS,
527 adcroft 1.6 I myTime,myIter,myThid)
528     ENDIF
529 adcroft 1.17 C Since passive tracers are configurable separately from T,S we
530     C call the multi-dimensional method for PTRACERS regardless
531     C of whether multiDimAdvection is set or not.
532     #ifdef ALLOW_PTRACERS
533     IF ( usePTRACERS ) THEN
534     CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
535     ENDIF
536     #endif /* ALLOW_PTRACERS */
537 adcroft 1.12 #endif /* DISABLE_MULTIDIM_ADVECTION */
538 adcroft 1.1
539     C-- Start of thermodynamics loop
540     DO k=Nr,1,-1
541     #ifdef ALLOW_AUTODIFF_TAMC
542     C? Patrick Is this formula correct?
543     cph Yes, but I rewrote it.
544     cph Also, the KappaR? need the index and subscript k!
545 heimbach 1.30 kkey = (itdkey-1)*Nr + k
546 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
547    
548     C-- km1 Points to level above k (=k-1)
549     C-- kup Cycles through 1,2 to point to layer above
550     C-- kDown Cycles through 2,1 to point to current layer
551    
552     km1 = MAX(1,k-1)
553     kup = 1+MOD(k+1,2)
554     kDown= 1+MOD(k,2)
555    
556     iMin = 1-OLx
557     iMax = sNx+OLx
558     jMin = 1-OLy
559     jMax = sNy+OLy
560    
561     C-- Get temporary terms used by tendency routines
562     CALL CALC_COMMON_FACTORS (
563     I bi,bj,iMin,iMax,jMin,jMax,k,
564     O xA,yA,uTrans,vTrans,rTrans,maskUp,
565     I myThid)
566 jmc 1.19
567     #ifdef ALLOW_GMREDI
568 heimbach 1.35
569 jmc 1.19 C-- Residual transp = Bolus transp + Eulerian transp
570     IF (useGMRedi) THEN
571     CALL GMREDI_CALC_UVFLOW(
572     & uTrans, vTrans, bi, bj, k, myThid)
573     IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
574     & rTrans, bi, bj, k, myThid)
575     ENDIF
576 heimbach 1.35
577     #ifdef ALLOW_AUTODIFF_TAMC
578     #ifdef GM_BOLUS_ADVEC
579     CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
580     CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
581     CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
582     #endif
583     #endif /* ALLOW_AUTODIFF_TAMC */
584    
585 jmc 1.19 #endif /* ALLOW_GMREDI */
586 adcroft 1.1
587     #ifdef ALLOW_AUTODIFF_TAMC
588     CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
589     CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
590     #endif /* ALLOW_AUTODIFF_TAMC */
591    
592     #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
593     C-- Calculate the total vertical diffusivity
594     CALL CALC_DIFFUSIVITY(
595     I bi,bj,iMin,iMax,jMin,jMax,k,
596     I maskUp,
597 heimbach 1.2 O KappaRT,KappaRS,
598 adcroft 1.1 I myThid)
599     #endif
600    
601     iMin = 1-OLx+2
602     iMax = sNx+OLx-1
603     jMin = 1-OLy+2
604     jMax = sNy+OLy-1
605    
606     C-- Calculate active tracer tendencies (gT,gS,...)
607     C and step forward storing result in gTnm1, gSnm1, etc.
608     IF ( tempStepping ) THEN
609     CALL CALC_GT(
610     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
611     I xA,yA,uTrans,vTrans,rTrans,maskUp,
612     I KappaRT,
613     U fVerT,
614 adcroft 1.7 I myTime,myIter,myThid)
615 adcroft 1.1 CALL TIMESTEP_TRACER(
616 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
617 adcroft 1.1 I theta, gT,
618     I myIter, myThid)
619     ENDIF
620 cheisey 1.32 cswdice ---- add ---
621     #ifdef ALLOW_THERM_SEAICE
622     if (k.eq.1) then
623     call ICE_FREEZE(bi, bj, iMin, iMax, jMin, jMax, myThid )
624     endif
625     #endif
626     cswdice -- end add ---
627 adcroft 1.1 IF ( saltStepping ) THEN
628     CALL CALC_GS(
629     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
630     I xA,yA,uTrans,vTrans,rTrans,maskUp,
631     I KappaRS,
632     U fVerS,
633 adcroft 1.7 I myTime,myIter,myThid)
634 adcroft 1.1 CALL TIMESTEP_TRACER(
635 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
636 adcroft 1.1 I salt, gS,
637     I myIter, myThid)
638     ENDIF
639     #ifdef ALLOW_PASSIVE_TRACER
640     IF ( tr1Stepping ) THEN
641     CALL CALC_GTR1(
642     I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
643     I xA,yA,uTrans,vTrans,rTrans,maskUp,
644     I KappaRT,
645     U fVerTr1,
646 heimbach 1.8 I myTime,myIter,myThid)
647 adcroft 1.1 CALL TIMESTEP_TRACER(
648 adcroft 1.3 I bi,bj,iMin,iMax,jMin,jMax,k,tracerAdvScheme,
649 adcroft 1.1 I Tr1, gTr1,
650 heimbach 1.8 I myIter,myThid)
651 adcroft 1.1 ENDIF
652     #endif
653 adcroft 1.17 #ifdef ALLOW_PTRACERS
654     IF ( usePTRACERS ) THEN
655     CALL PTRACERS_INTEGERATE(
656     I bi,bj,k,
657     I xA,yA,uTrans,vTrans,rTrans,maskUp,
658     X KappaRS,
659     I myIter,myTime,myThid)
660     ENDIF
661     #endif /* ALLOW_PTRACERS */
662 adcroft 1.1
663     #ifdef ALLOW_OBCS
664     C-- Apply open boundary conditions
665     IF (useOBCS) THEN
666 adcroft 1.7 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
667 adcroft 1.1 END IF
668     #endif /* ALLOW_OBCS */
669    
670     C-- Freeze water
671 heimbach 1.29 IF ( allowFreezing .AND. .NOT. useSEAICE ) THEN
672 adcroft 1.1 #ifdef ALLOW_AUTODIFF_TAMC
673 heimbach 1.11 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
674 adcroft 1.1 CADJ & , key = kkey, byte = isbyte
675     #endif /* ALLOW_AUTODIFF_TAMC */
676     CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
677     END IF
678    
679     C-- end of thermodynamic k loop (Nr:1)
680     ENDDO
681 cheisey 1.31
682     cswdice -- add ---
683 cheisey 1.32 #ifdef ALLOW_THERM_SEAICE
684 cheisey 1.31 c timeaveraging for ice model values
685     CALL ICE_AVE(bi,bj,iMin,iMax,jMin,jMax,myThid )
686     #endif
687     cswdice --- end add ---
688    
689    
690    
691 adcroft 1.1
692     C-- Implicit diffusion
693     IF (implicitDiffusion) THEN
694    
695     IF (tempStepping) THEN
696     #ifdef ALLOW_AUTODIFF_TAMC
697 heimbach 1.30 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
698 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
699     CALL IMPLDIFF(
700     I bi, bj, iMin, iMax, jMin, jMax,
701     I deltaTtracer, KappaRT, recip_HFacC,
702 adcroft 1.7 U gT,
703 adcroft 1.1 I myThid )
704     ENDIF
705    
706     IF (saltStepping) THEN
707     #ifdef ALLOW_AUTODIFF_TAMC
708 heimbach 1.30 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
709 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
710     CALL IMPLDIFF(
711     I bi, bj, iMin, iMax, jMin, jMax,
712     I deltaTtracer, KappaRS, recip_HFacC,
713 adcroft 1.7 U gS,
714 adcroft 1.1 I myThid )
715     ENDIF
716    
717     #ifdef ALLOW_PASSIVE_TRACER
718     IF (tr1Stepping) THEN
719     #ifdef ALLOW_AUTODIFF_TAMC
720 heimbach 1.30 CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
721 adcroft 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
722     CALL IMPLDIFF(
723     I bi, bj, iMin, iMax, jMin, jMax,
724     I deltaTtracer, KappaRT, recip_HFacC,
725 adcroft 1.7 U gTr1,
726 adcroft 1.1 I myThid )
727     ENDIF
728     #endif
729    
730 adcroft 1.17 #ifdef ALLOW_PTRACERS
731     C Vertical diffusion (implicit) for passive tracers
732     IF ( usePTRACERS ) THEN
733     CALL PTRACERS_IMPLDIFF( bi,bj,KappaRS,myThid )
734     ENDIF
735     #endif /* ALLOW_PTRACERS */
736    
737 adcroft 1.1 #ifdef ALLOW_OBCS
738     C-- Apply open boundary conditions
739     IF (useOBCS) THEN
740     DO K=1,Nr
741 adcroft 1.7 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
742 adcroft 1.1 ENDDO
743     END IF
744     #endif /* ALLOW_OBCS */
745    
746     C-- End If implicitDiffusion
747     ENDIF
748 heimbach 1.22
749 jmc 1.39 #ifdef ALLOW_TIMEAVE
750     IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
751     CALL TIMEAVE_CUMUL_1T(ConvectCountTave, ConvectCount,
752     I Nr, deltaTclock, bi, bj, myThid)
753     ENDIF
754     useVariableK = useKPP .OR. useGMredi .OR. ivdc_kappa.NE.0.
755     IF (taveFreq.GT.0. .AND. useVariableK ) THEN
756     IF (implicitDiffusion) THEN
757     CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRT,
758     I Nr, 3, deltaTclock, bi, bj, myThid)
759     ELSE
760     CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
761     I Nr, 3, deltaTclock, bi, bj, myThid)
762     ENDIF
763     ENDIF
764     #endif /* ALLOW_TIMEAVE */
765    
766 heimbach 1.22 #endif /* SINGLE_LAYER_MODE */
767 adcroft 1.1
768 jmc 1.39 C-- end bi,bj loops.
769 adcroft 1.1 ENDDO
770     ENDDO
771 adcroft 1.17
772     #ifndef DISABLE_DEBUGMODE
773     If (debugMode) THEN
774     CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
775     CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
776     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
777     CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
778     CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
779     CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
780     CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
781     CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
782     CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
783 adcroft 1.18 #ifdef ALLOW_PTRACERS
784     IF ( usePTRACERS ) THEN
785     CALL PTRACERS_DEBUG(myThid)
786     ENDIF
787     #endif /* ALLOW_PTRACERS */
788 adcroft 1.17 ENDIF
789     #endif
790 adcroft 1.1
791     RETURN
792     END

  ViewVC Help
Powered by ViewVC 1.1.22