/[MITgcm]/MITgcm/model/src/thermodynamics.F
ViewVC logotype

Contents of /MITgcm/model/src/thermodynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.128 - (show annotations) (download)
Tue Aug 10 17:58:30 2010 UTC (13 years, 10 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint62j
Changes since 1.127: +13 -1 lines
Adjoint related modifications -- allowing the
use of implicit vertical advection in adjoint model.

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

  ViewVC Help
Powered by ViewVC 1.1.22