1 |
C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.104 2006/06/07 01:55:13 heimbach 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 |
#ifdef ALLOW_GENERIC_ADVDIFF |
81 |
#include "GAD.h" |
82 |
#endif |
83 |
#ifdef ALLOW_PTRACERS |
84 |
#include "PTRACERS_SIZE.h" |
85 |
#include "PTRACERS.h" |
86 |
#endif |
87 |
#ifdef ALLOW_TIMEAVE |
88 |
#include "TIMEAVE_STATV.h" |
89 |
#endif |
90 |
|
91 |
#ifdef ALLOW_AUTODIFF_TAMC |
92 |
# include "tamc.h" |
93 |
# include "tamc_keys.h" |
94 |
# include "FFIELDS.h" |
95 |
# include "EOS.h" |
96 |
# ifdef ALLOW_KPP |
97 |
# include "KPP.h" |
98 |
# endif |
99 |
# ifdef ALLOW_GMREDI |
100 |
# include "GMREDI.h" |
101 |
# endif |
102 |
# ifdef ALLOW_EBM |
103 |
# include "EBM.h" |
104 |
# endif |
105 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
106 |
|
107 |
C !INPUT/OUTPUT PARAMETERS: |
108 |
C == Routine arguments == |
109 |
C myTime - Current time in simulation |
110 |
C myIter - Current iteration number in simulation |
111 |
C myThid - Thread number for this instance of the routine. |
112 |
_RL myTime |
113 |
INTEGER myIter |
114 |
INTEGER myThid |
115 |
|
116 |
#ifdef ALLOW_GENERIC_ADVDIFF |
117 |
C !LOCAL VARIABLES: |
118 |
C == Local variables |
119 |
C xA, yA - Per block temporaries holding face areas |
120 |
C uFld, vFld, wFld - Local copy of velocity field (3 components) |
121 |
C uTrans, vTrans, rTrans - Per block temporaries holding flow transport |
122 |
C o uTrans: Zonal transport |
123 |
C o vTrans: Meridional transport |
124 |
C o rTrans: Vertical transport |
125 |
C rTransKp1 o vertical volume transp. at interface k+1 |
126 |
C maskUp o maskUp: land/water mask for W points |
127 |
C fVer[STUV] o fVer: Vertical flux term - note fVer |
128 |
C is "pipelined" in the vertical |
129 |
C so we need an fVer for each |
130 |
C variable. |
131 |
C kappaRT, - Total diffusion in vertical at level k, for T and S |
132 |
C kappaRS (background + spatially varying, isopycnal term). |
133 |
C kappaRTr - Total diffusion in vertical at level k, |
134 |
C for each passive Tracer |
135 |
C kappaRk - Total diffusion in vertical, all levels, 1 tracer |
136 |
C useVariableK = T when vertical diffusion is not constant |
137 |
C iMin, iMax - Ranges and sub-block indices on which calculations |
138 |
C jMin, jMax are applied. |
139 |
C bi, bj |
140 |
C k, kup, - Index for layer above and below. kup and kDown |
141 |
C kDown, km1 are switched with layer to be the appropriate |
142 |
C index into fVerTerm. |
143 |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
144 |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
145 |
_RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
146 |
_RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
147 |
_RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
148 |
_RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
149 |
_RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
150 |
_RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
151 |
_RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
152 |
_RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
153 |
_RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
154 |
_RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
155 |
_RL kappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
156 |
_RL kappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
157 |
#ifdef ALLOW_PTRACERS |
158 |
_RL fVerP (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num) |
159 |
_RL kappaRTr(1-Olx:sNx+Olx,1-Oly:sNy+Oly,PTRACERS_num) |
160 |
#endif |
161 |
_RL kappaRk (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
162 |
INTEGER iMin, iMax |
163 |
INTEGER jMin, jMax |
164 |
INTEGER bi, bj |
165 |
INTEGER i, j |
166 |
INTEGER k, km1, kup, kDown |
167 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
168 |
INTEGER iterNb, m1, m2 |
169 |
#endif |
170 |
#ifdef ALLOW_TIMEAVE |
171 |
LOGICAL useVariableK |
172 |
#endif |
173 |
#ifdef ALLOW_PTRACERS |
174 |
INTEGER iTracer, ip |
175 |
#endif |
176 |
|
177 |
CEOP |
178 |
|
179 |
#ifdef ALLOW_DEBUG |
180 |
IF ( debugLevel .GE. debLevB ) |
181 |
& CALL DEBUG_ENTER('THERMODYNAMICS',myThid) |
182 |
#endif |
183 |
|
184 |
#ifdef ALLOW_AUTODIFF_TAMC |
185 |
C-- dummy statement to end declaration part |
186 |
ikey = 1 |
187 |
itdkey = 1 |
188 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
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 |
CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS |
200 |
CHPF$& ,utrans,vtrans,xA,yA |
201 |
CHPF$& ,kappaRT,kappaRS |
202 |
CHPF$& ) |
203 |
# ifdef ALLOW_PTRACERS |
204 |
CHPF$ INDEPENDENT, NEW (fVerP,kappaRTr) |
205 |
# endif |
206 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
207 |
|
208 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
209 |
|
210 |
#ifdef ALLOW_AUTODIFF_TAMC |
211 |
act1 = bi - myBxLo(myThid) |
212 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
213 |
act2 = bj - myByLo(myThid) |
214 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
215 |
act3 = myThid - 1 |
216 |
max3 = nTx*nTy |
217 |
act4 = ikey_dynamics - 1 |
218 |
itdkey = (act1 + 1) + act2*max1 |
219 |
& + act3*max1*max2 |
220 |
& + act4*max1*max2*max3 |
221 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
222 |
|
223 |
C-- Set up work arrays with valid (i.e. not NaN) values |
224 |
C These inital values do not alter the numerical results. They |
225 |
C just ensure that all memory references are to valid floating |
226 |
C point numbers. This prevents spurious hardware signals due to |
227 |
C uninitialised but inert locations. |
228 |
|
229 |
DO j=1-OLy,sNy+OLy |
230 |
DO i=1-OLx,sNx+OLx |
231 |
xA(i,j) = 0. _d 0 |
232 |
yA(i,j) = 0. _d 0 |
233 |
uTrans(i,j) = 0. _d 0 |
234 |
vTrans(i,j) = 0. _d 0 |
235 |
rTrans (i,j) = 0. _d 0 |
236 |
rTransKp1(i,j) = 0. _d 0 |
237 |
fVerT (i,j,1) = 0. _d 0 |
238 |
fVerT (i,j,2) = 0. _d 0 |
239 |
fVerS (i,j,1) = 0. _d 0 |
240 |
fVerS (i,j,2) = 0. _d 0 |
241 |
kappaRT(i,j) = 0. _d 0 |
242 |
kappaRS(i,j) = 0. _d 0 |
243 |
ENDDO |
244 |
ENDDO |
245 |
|
246 |
DO k=1,Nr |
247 |
DO j=1-OLy,sNy+OLy |
248 |
DO i=1-OLx,sNx+OLx |
249 |
C This is currently also used by IVDC and Diagnostics |
250 |
kappaRk(i,j,k) = 0. _d 0 |
251 |
C- tracer tendency needs to be set to zero (moved here from gad_calc_rhs): |
252 |
gT(i,j,k,bi,bj) = 0. _d 0 |
253 |
gS(i,j,k,bi,bj) = 0. _d 0 |
254 |
ENDDO |
255 |
ENDDO |
256 |
ENDDO |
257 |
|
258 |
#ifdef ALLOW_PTRACERS |
259 |
IF ( usePTRACERS ) THEN |
260 |
DO ip=1,PTRACERS_num |
261 |
DO j=1-OLy,sNy+OLy |
262 |
DO i=1-OLx,sNx+OLx |
263 |
fVerP (i,j,1,ip) = 0. _d 0 |
264 |
fVerP (i,j,2,ip) = 0. _d 0 |
265 |
kappaRTr(i,j,ip) = 0. _d 0 |
266 |
ENDDO |
267 |
ENDDO |
268 |
ENDDO |
269 |
C- set tracer tendency to zero: |
270 |
DO iTracer=1,PTRACERS_num |
271 |
DO k=1,Nr |
272 |
DO j=1-OLy,sNy+OLy |
273 |
DO i=1-OLx,sNx+OLx |
274 |
gPTr(i,j,k,bi,bj,itracer) = 0. _d 0 |
275 |
ENDDO |
276 |
ENDDO |
277 |
ENDDO |
278 |
ENDDO |
279 |
ENDIF |
280 |
#endif |
281 |
|
282 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
283 |
C- Apply AB on T,S : |
284 |
iterNb = myIter |
285 |
IF (staggerTimeStep) iterNb = myIter - 1 |
286 |
m1 = 1 + MOD(iterNb+1,2) |
287 |
m2 = 1 + MOD( iterNb ,2) |
288 |
C compute T^n+1/2 (stored in gtNm) extrapolating T forward in time |
289 |
IF ( AdamsBashforth_T ) CALL ADAMS_BASHFORTH3( |
290 |
I bi, bj, 0, |
291 |
U theta, gtNm, |
292 |
I tempStartAB, iterNb, myThid ) |
293 |
C compute S^n+1/2 (stored in gsNm) extrapolating S forward in time |
294 |
IF ( AdamsBashforth_S ) CALL ADAMS_BASHFORTH3( |
295 |
I bi, bj, 0, |
296 |
U salt, gsNm, |
297 |
I saltStartAB, iterNb, myThid ) |
298 |
#endif /* ALLOW_ADAMSBASHFORTH_3 */ |
299 |
|
300 |
c iMin = 1-OLx |
301 |
c iMax = sNx+OLx |
302 |
c jMin = 1-OLy |
303 |
c jMax = sNy+OLy |
304 |
|
305 |
#ifdef ALLOW_AUTODIFF_TAMC |
306 |
cph avoids recomputation of integrate_for_w |
307 |
CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
308 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
309 |
|
310 |
C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h |
311 |
C-- MOST of THERMODYNAMICS will be disabled |
312 |
#ifndef SINGLE_LAYER_MODE |
313 |
|
314 |
#ifdef ALLOW_AUTODIFF_TAMC |
315 |
CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
316 |
CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
317 |
CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
318 |
CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
319 |
cph-test |
320 |
# ifdef ALLOW_DEPTH_CONTROL |
321 |
CADJ STORE gtnm1(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
322 |
CADJ STORE gsnm1(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
323 |
# endif |
324 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
325 |
|
326 |
#ifndef DISABLE_MULTIDIM_ADVECTION |
327 |
C-- Some advection schemes are better calculated using a multi-dimensional |
328 |
C method in the absence of any other terms and, if used, is done here. |
329 |
C |
330 |
C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h |
331 |
C The default is to use multi-dimensinal advection for non-linear advection |
332 |
C schemes. However, for the sake of efficiency of the adjoint it is necessary |
333 |
C to be able to exclude this scheme to avoid excessive storage and |
334 |
C recomputation. It *is* differentiable, if you need it. |
335 |
C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to |
336 |
C disable this section of code. |
337 |
IF (tempMultiDimAdvec) THEN |
338 |
#ifdef ALLOW_DEBUG |
339 |
IF ( debugLevel .GE. debLevB ) |
340 |
& CALL DEBUG_CALL('GAD_ADVECTION',myThid) |
341 |
#endif |
342 |
CALL GAD_ADVECTION( |
343 |
I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme, |
344 |
I GAD_TEMPERATURE, |
345 |
I uVel, vVel, wVel, theta, |
346 |
O gT, |
347 |
I bi,bj,myTime,myIter,myThid) |
348 |
ENDIF |
349 |
IF (saltMultiDimAdvec) THEN |
350 |
#ifdef ALLOW_DEBUG |
351 |
IF ( debugLevel .GE. debLevB ) |
352 |
& CALL DEBUG_CALL('GAD_ADVECTION',myThid) |
353 |
#endif |
354 |
CALL GAD_ADVECTION( |
355 |
I saltImplVertAdv, saltAdvScheme, saltVertAdvScheme, |
356 |
I GAD_SALINITY, |
357 |
I uVel, vVel, wVel, salt, |
358 |
O gS, |
359 |
I bi,bj,myTime,myIter,myThid) |
360 |
ENDIF |
361 |
|
362 |
C Since passive tracers are configurable separately from T,S we |
363 |
C call the multi-dimensional method for PTRACERS regardless |
364 |
C of whether multiDimAdvection is set or not. |
365 |
#ifdef ALLOW_PTRACERS |
366 |
IF ( usePTRACERS ) THEN |
367 |
#ifdef ALLOW_DEBUG |
368 |
IF ( debugLevel .GE. debLevB ) |
369 |
& CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid) |
370 |
#endif |
371 |
CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid ) |
372 |
ENDIF |
373 |
#endif /* ALLOW_PTRACERS */ |
374 |
#endif /* DISABLE_MULTIDIM_ADVECTION */ |
375 |
|
376 |
#ifdef ALLOW_DEBUG |
377 |
IF ( debugLevel .GE. debLevB ) |
378 |
& CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid) |
379 |
#endif |
380 |
|
381 |
C-- Start of thermodynamics loop |
382 |
DO k=Nr,1,-1 |
383 |
#ifdef ALLOW_AUTODIFF_TAMC |
384 |
C? Patrick Is this formula correct? |
385 |
cph Yes, but I rewrote it. |
386 |
cph Also, the kappaR? need the index and subscript k! |
387 |
kkey = (itdkey-1)*Nr + k |
388 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
389 |
|
390 |
C-- km1 Points to level above k (=k-1) |
391 |
C-- kup Cycles through 1,2 to point to layer above |
392 |
C-- kDown Cycles through 2,1 to point to current layer |
393 |
|
394 |
km1 = MAX(1,k-1) |
395 |
kup = 1+MOD(k+1,2) |
396 |
kDown= 1+MOD(k,2) |
397 |
|
398 |
iMin = 1-OLx |
399 |
iMax = sNx+OLx |
400 |
jMin = 1-OLy |
401 |
jMax = sNy+OLy |
402 |
|
403 |
IF (k.EQ.Nr) THEN |
404 |
DO j=1-Oly,sNy+Oly |
405 |
DO i=1-Olx,sNx+Olx |
406 |
rTransKp1(i,j) = 0. _d 0 |
407 |
ENDDO |
408 |
ENDDO |
409 |
ELSE |
410 |
DO j=1-Oly,sNy+Oly |
411 |
DO i=1-Olx,sNx+Olx |
412 |
rTransKp1(i,j) = rTrans(i,j) |
413 |
ENDDO |
414 |
ENDDO |
415 |
ENDIF |
416 |
#ifdef ALLOW_AUTODIFF_TAMC |
417 |
CADJ STORE rTransKp1(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte |
418 |
#endif |
419 |
|
420 |
C-- Get temporary terms used by tendency routines : |
421 |
C- Calculate horizontal "volume transport" through tracer cell face |
422 |
CALL CALC_COMMON_FACTORS ( |
423 |
I uVel, vVel, |
424 |
O uFld, vFld, uTrans, vTrans, xA, yA, |
425 |
I k,bi,bj, myThid ) |
426 |
|
427 |
C- Calculate vertical "volume transport" through tracer cell face |
428 |
IF (k.EQ.1) THEN |
429 |
C- Surface interface : |
430 |
DO j=1-Oly,sNy+Oly |
431 |
DO i=1-Olx,sNx+Olx |
432 |
wFld(i,j) = 0. _d 0 |
433 |
maskUp(i,j) = 0. _d 0 |
434 |
rTrans(i,j) = 0. _d 0 |
435 |
ENDDO |
436 |
ENDDO |
437 |
ELSE |
438 |
C- Interior interface : |
439 |
DO j=1-Oly,sNy+Oly |
440 |
DO i=1-Olx,sNx+Olx |
441 |
wFld(i,j) = wVel(i,j,k,bi,bj) |
442 |
maskUp(i,j) = maskC(i,j,k-1,bi,bj)*maskC(i,j,k,bi,bj) |
443 |
rTrans(i,j) = wFld(i,j)*rA(i,j,bi,bj)*maskUp(i,j) |
444 |
ENDDO |
445 |
ENDDO |
446 |
ENDIF |
447 |
|
448 |
#ifdef ALLOW_GMREDI |
449 |
C-- Residual transp = Bolus transp + Eulerian transp |
450 |
IF (useGMRedi) THEN |
451 |
CALL GMREDI_CALC_UVFLOW( |
452 |
U uFld, vFld, uTrans, vTrans, |
453 |
I k, bi, bj, myThid ) |
454 |
IF (K.GE.2) CALL GMREDI_CALC_WFLOW( |
455 |
U wFld, rTrans, |
456 |
I k, bi, bj, myThid ) |
457 |
ENDIF |
458 |
# ifdef ALLOW_AUTODIFF_TAMC |
459 |
CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte |
460 |
# ifdef GM_BOLUS_ADVEC |
461 |
CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte |
462 |
CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte |
463 |
# endif |
464 |
# endif /* ALLOW_AUTODIFF_TAMC */ |
465 |
#endif /* ALLOW_GMREDI */ |
466 |
|
467 |
#ifdef INCLUDE_CALC_DIFFUSIVITY_CALL |
468 |
C-- Calculate the total vertical diffusivity |
469 |
IF ( .NOT.implicitDiffusion ) THEN |
470 |
CALL CALC_DIFFUSIVITY( |
471 |
I bi,bj,iMin,iMax,jMin,jMax,k, |
472 |
I maskUp, |
473 |
O kappaRT,kappaRS, |
474 |
I myThid) |
475 |
ENDIF |
476 |
# ifdef ALLOW_AUTODIFF_TAMC |
477 |
CADJ STORE kappaRT(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte |
478 |
CADJ STORE kappaRS(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte |
479 |
# endif /* ALLOW_AUTODIFF_TAMC */ |
480 |
#endif |
481 |
|
482 |
iMin = 1-OLx+2 |
483 |
iMax = sNx+OLx-1 |
484 |
jMin = 1-OLy+2 |
485 |
jMax = sNy+OLy-1 |
486 |
|
487 |
C-- Calculate active tracer tendencies (gT,gS,...) |
488 |
C and step forward storing result in gT, gS, etc. |
489 |
C-- |
490 |
# ifdef ALLOW_AUTODIFF_TAMC |
491 |
# if ((defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)) && (defined ALLOW_GMREDI) |
492 |
CADJ STORE kux(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
493 |
CADJ STORE kvy(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
494 |
# ifdef GM_EXTRA_DIAGONAL |
495 |
CADJ STORE kuz(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
496 |
CADJ STORE kvz(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
497 |
# endif |
498 |
# endif |
499 |
# endif /* ALLOW_AUTODIFF_TAMC */ |
500 |
IF ( tempStepping ) THEN |
501 |
#ifdef ALLOW_AUTODIFF_TAMC |
502 |
CADJ STORE gTnm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
503 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
504 |
CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
505 |
# endif |
506 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
507 |
CALL CALC_GT( |
508 |
I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown, |
509 |
I xA, yA, maskUp, uFld, vFld, wFld, |
510 |
I uTrans, vTrans, rTrans, rTransKp1, |
511 |
I kappaRT, |
512 |
U fVerT, |
513 |
I myTime,myIter,myThid) |
514 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
515 |
IF ( AdamsBashforth_T ) THEN |
516 |
CALL TIMESTEP_TRACER( |
517 |
I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme, |
518 |
I gtNm(1-Olx,1-Oly,1,1,1,m2), |
519 |
U gT, |
520 |
I myIter, myThid) |
521 |
ELSE |
522 |
#endif |
523 |
CALL TIMESTEP_TRACER( |
524 |
I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme, |
525 |
I theta, |
526 |
U gT, |
527 |
I myIter, myThid) |
528 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
529 |
ENDIF |
530 |
#endif |
531 |
ENDIF |
532 |
|
533 |
IF ( saltStepping ) THEN |
534 |
#ifdef ALLOW_AUTODIFF_TAMC |
535 |
CADJ STORE gSnm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
536 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
537 |
CADJ STORE gS(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
538 |
# endif |
539 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
540 |
CALL CALC_GS( |
541 |
I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown, |
542 |
I xA, yA, maskUp, uFld, vFld, wFld, |
543 |
I uTrans, vTrans, rTrans, rTransKp1, |
544 |
I kappaRS, |
545 |
U fVerS, |
546 |
I myTime,myIter,myThid) |
547 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
548 |
IF ( AdamsBashforth_S ) THEN |
549 |
CALL TIMESTEP_TRACER( |
550 |
I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme, |
551 |
I gsNm(1-Olx,1-Oly,1,1,1,m2), |
552 |
U gS, |
553 |
I myIter, myThid) |
554 |
ELSE |
555 |
#endif |
556 |
CALL TIMESTEP_TRACER( |
557 |
I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme, |
558 |
I salt, |
559 |
U gS, |
560 |
I myIter, myThid) |
561 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
562 |
ENDIF |
563 |
#endif |
564 |
ENDIF |
565 |
|
566 |
#ifdef ALLOW_PTRACERS |
567 |
IF ( usePTRACERS ) THEN |
568 |
IF ( .NOT.implicitDiffusion ) THEN |
569 |
CALL PTRACERS_CALC_DIFF( |
570 |
I bi,bj,iMin,iMax,jMin,jMax,k, |
571 |
I maskUp, |
572 |
O kappaRTr, |
573 |
I myThid) |
574 |
ENDIF |
575 |
# ifdef ALLOW_AUTODIFF_TAMC |
576 |
CADJ STORE kappaRTr(:,:,:) = comlev1_bibj_k, key=kkey, byte=isbyte |
577 |
# endif /* ALLOW_AUTODIFF_TAMC */ |
578 |
CALL PTRACERS_INTEGRATE( |
579 |
I bi,bj,k, |
580 |
I xA, yA, maskUp, uFld, vFld, wFld, |
581 |
I uTrans, vTrans, rTrans, rTransKp1, |
582 |
I kappaRTr, |
583 |
U fVerP, |
584 |
I myTime,myIter,myThid) |
585 |
ENDIF |
586 |
#endif /* ALLOW_PTRACERS */ |
587 |
|
588 |
#ifdef ALLOW_OBCS |
589 |
C-- Apply open boundary conditions |
590 |
IF (useOBCS) THEN |
591 |
CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid ) |
592 |
END IF |
593 |
#endif /* ALLOW_OBCS */ |
594 |
|
595 |
C-- Freeze water |
596 |
C this bit of code is left here for backward compatibility. |
597 |
C freezing at surface level has been moved to FORWARD_STEP |
598 |
IF ( useOldFreezing .AND. .NOT. useSEAICE |
599 |
& .AND. .NOT.(useThSIce.AND.k.EQ.1) ) THEN |
600 |
#ifdef ALLOW_AUTODIFF_TAMC |
601 |
CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k |
602 |
CADJ & , key = kkey, byte = isbyte |
603 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
604 |
CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid ) |
605 |
ENDIF |
606 |
|
607 |
C-- end of thermodynamic k loop (Nr:1) |
608 |
ENDDO |
609 |
|
610 |
C All explicit advection/diffusion/sources should now be |
611 |
C done. The updated tracer field is in gPtr. Accumalate |
612 |
C explicit tendency and also reset gPtr to initial tracer |
613 |
C field for implicit matrix calculation |
614 |
|
615 |
#ifdef ALLOW_MATRIX |
616 |
IF (useMATRIX) |
617 |
& CALL MATRIX_STORE_TENDENCY_EXP(bi,bj, myTime,myIter,myThid) |
618 |
#endif |
619 |
|
620 |
iMin = 1 |
621 |
iMax = sNx |
622 |
jMin = 1 |
623 |
jMax = sNy |
624 |
|
625 |
C-- Implicit vertical advection & diffusion |
626 |
IF ( tempStepping .AND. implicitDiffusion ) THEN |
627 |
CALL CALC_3D_DIFFUSIVITY( |
628 |
I bi,bj,iMin,iMax,jMin,jMax, |
629 |
I GAD_TEMPERATURE, useGMredi, useKPP, |
630 |
O kappaRk, |
631 |
I myThid) |
632 |
ENDIF |
633 |
#ifdef INCLUDE_IMPLVERTADV_CODE |
634 |
IF ( tempImplVertAdv ) THEN |
635 |
CALL GAD_IMPLICIT_R( |
636 |
I tempImplVertAdv, tempAdvScheme, GAD_TEMPERATURE, |
637 |
I kappaRk, wVel, theta, |
638 |
U gT, |
639 |
I bi, bj, myTime, myIter, myThid ) |
640 |
ELSEIF ( tempStepping .AND. implicitDiffusion ) THEN |
641 |
#else /* INCLUDE_IMPLVERTADV_CODE */ |
642 |
IF ( tempStepping .AND. implicitDiffusion ) THEN |
643 |
#endif /* INCLUDE_IMPLVERTADV_CODE */ |
644 |
#ifdef ALLOW_AUTODIFF_TAMC |
645 |
CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
646 |
CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte |
647 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
648 |
CALL IMPLDIFF( |
649 |
I bi, bj, iMin, iMax, jMin, jMax, |
650 |
I GAD_TEMPERATURE, kappaRk, recip_hFacC, |
651 |
U gT, |
652 |
I myThid ) |
653 |
ENDIF |
654 |
|
655 |
#ifdef ALLOW_TIMEAVE |
656 |
useVariableK = useKPP .OR. usePP81 .OR. useMY82 .OR. useGGL90 |
657 |
& .OR. useGMredi .OR. ivdc_kappa.NE.0. |
658 |
IF (taveFreq.GT.0. .AND. useVariableK ) THEN |
659 |
IF (implicitDiffusion) THEN |
660 |
CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRk, |
661 |
I Nr, 3, deltaTclock, bi, bj, myThid) |
662 |
c ELSE |
663 |
c CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT, |
664 |
c I Nr, 3, deltaTclock, bi, bj, myThid) |
665 |
ENDIF |
666 |
ENDIF |
667 |
#endif /* ALLOW_TIMEAVE */ |
668 |
|
669 |
IF ( saltStepping .AND. implicitDiffusion ) THEN |
670 |
CALL CALC_3D_DIFFUSIVITY( |
671 |
I bi,bj,iMin,iMax,jMin,jMax, |
672 |
I GAD_SALINITY, useGMredi, useKPP, |
673 |
O kappaRk, |
674 |
I myThid) |
675 |
ENDIF |
676 |
|
677 |
#ifdef INCLUDE_IMPLVERTADV_CODE |
678 |
IF ( saltImplVertAdv ) THEN |
679 |
CALL GAD_IMPLICIT_R( |
680 |
I saltImplVertAdv, saltAdvScheme, GAD_SALINITY, |
681 |
I kappaRk, wVel, salt, |
682 |
U gS, |
683 |
I bi, bj, myTime, myIter, myThid ) |
684 |
ELSEIF ( saltStepping .AND. implicitDiffusion ) THEN |
685 |
#else /* INCLUDE_IMPLVERTADV_CODE */ |
686 |
IF ( saltStepping .AND. implicitDiffusion ) THEN |
687 |
#endif /* INCLUDE_IMPLVERTADV_CODE */ |
688 |
#ifdef ALLOW_AUTODIFF_TAMC |
689 |
CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
690 |
CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte |
691 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
692 |
CALL IMPLDIFF( |
693 |
I bi, bj, iMin, iMax, jMin, jMax, |
694 |
I GAD_SALINITY, kappaRk, recip_hFacC, |
695 |
U gS, |
696 |
I myThid ) |
697 |
ENDIF |
698 |
|
699 |
#ifdef ALLOW_PTRACERS |
700 |
IF ( usePTRACERS ) THEN |
701 |
C-- Vertical advection/diffusion (implicit) for passive tracers |
702 |
CALL PTRACERS_IMPLICIT( |
703 |
U kappaRk, |
704 |
I bi, bj, myTime, myIter, myThid ) |
705 |
ENDIF |
706 |
#endif /* ALLOW_PTRACERS */ |
707 |
|
708 |
#ifdef ALLOW_OBCS |
709 |
C-- Apply open boundary conditions |
710 |
IF ( ( implicitDiffusion |
711 |
& .OR. tempImplVertAdv |
712 |
& .OR. saltImplVertAdv |
713 |
& ) .AND. useOBCS ) THEN |
714 |
DO K=1,Nr |
715 |
CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid ) |
716 |
ENDDO |
717 |
ENDIF |
718 |
#endif /* ALLOW_OBCS */ |
719 |
|
720 |
#endif /* SINGLE_LAYER_MODE */ |
721 |
|
722 |
C-- end bi,bj loops. |
723 |
ENDDO |
724 |
ENDDO |
725 |
|
726 |
#ifdef ALLOW_DEBUG |
727 |
If (debugMode) THEN |
728 |
CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid) |
729 |
CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid) |
730 |
CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid) |
731 |
CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid) |
732 |
CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid) |
733 |
CALL DEBUG_STATS_RL(Nr,gT,'Gt (THERMODYNAMICS)',myThid) |
734 |
CALL DEBUG_STATS_RL(Nr,gS,'Gs (THERMODYNAMICS)',myThid) |
735 |
#ifndef ALLOW_ADAMSBASHFORTH_3 |
736 |
CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (THERMODYNAMICS)',myThid) |
737 |
CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (THERMODYNAMICS)',myThid) |
738 |
#endif |
739 |
#ifdef ALLOW_PTRACERS |
740 |
IF ( usePTRACERS ) THEN |
741 |
CALL PTRACERS_DEBUG(myThid) |
742 |
ENDIF |
743 |
#endif /* ALLOW_PTRACERS */ |
744 |
ENDIF |
745 |
#endif /* ALLOW_DEBUG */ |
746 |
|
747 |
#ifdef ALLOW_DEBUG |
748 |
IF ( debugLevel .GE. debLevB ) |
749 |
& CALL DEBUG_LEAVE('THERMODYNAMICS',myThid) |
750 |
#endif |
751 |
|
752 |
#endif /* ALLOW_GENERIC_ADVDIFF */ |
753 |
|
754 |
RETURN |
755 |
END |