/[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.50 - (hide annotations) (download)
Tue Oct 7 04:31:30 2003 UTC (20 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint51i_pre
Changes since 1.49: +7 -2 lines
fix Problem with bulk_force & therm_seaice.

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

  ViewVC Help
Powered by ViewVC 1.1.22