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