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