1 |
C $Header: /u/gcmpack/MITgcm/model/src/temp_integrate.F,v 1.3 2013/11/27 23:58:25 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 |
|
10 |
CBOP |
11 |
C !ROUTINE: TEMP_INTEGRATE |
12 |
C !INTERFACE: |
13 |
SUBROUTINE TEMP_INTEGRATE( |
14 |
I bi, bj, |
15 |
I uFld, vFld, wFld, |
16 |
U KappaRk, |
17 |
I myTime, myIter, myThid ) |
18 |
C !DESCRIPTION: \bv |
19 |
C *==========================================================* |
20 |
C | SUBROUTINE TEMP_INTEGRATE |
21 |
C | o Calculate tendency for temperature |
22 |
C | and integrates forward in time. |
23 |
C *==========================================================* |
24 |
C | A procedure called EXTERNAL_FORCING_T is called from |
25 |
C | here. These procedures can be used to add per problem |
26 |
C | heat flux source terms. |
27 |
C | Note: Although it is slightly counter-intuitive the |
28 |
C | EXTERNAL_FORCING routine is not the place to put |
29 |
C | file I/O. Instead files that are required to |
30 |
C | calculate the external source terms are generally |
31 |
C | read during the model main loop. This makes the |
32 |
C | logistics of multi-processing simpler and also |
33 |
C | makes the adjoint generation simpler. It also |
34 |
C | allows for I/O to overlap computation where that |
35 |
C | is supported by hardware. |
36 |
C | Aside from the problem specific term the code here |
37 |
C | forms the tendency terms due to advection and mixing |
38 |
C | The baseline implementation here uses a centered |
39 |
C | difference form for the advection term and a tensorial |
40 |
C | divergence of a flux form for the diffusive term. The |
41 |
C | diffusive term is formulated so that isopycnal mixing |
42 |
C | and GM-style subgrid-scale terms can be incorporated by |
43 |
C | simply setting the diffusion tensor terms appropriately. |
44 |
C *==========================================================* |
45 |
C \ev |
46 |
|
47 |
C !USES: |
48 |
IMPLICIT NONE |
49 |
C == GLobal variables == |
50 |
#include "SIZE.h" |
51 |
#include "DYNVARS.h" |
52 |
#include "EEPARAMS.h" |
53 |
#include "PARAMS.h" |
54 |
#include "RESTART.h" |
55 |
#ifdef ALLOW_GENERIC_ADVDIFF |
56 |
# include "GAD.h" |
57 |
# include "GAD_SOM_VARS.h" |
58 |
#endif |
59 |
#ifdef ALLOW_AUTODIFF_TAMC |
60 |
# include "tamc.h" |
61 |
# include "tamc_keys.h" |
62 |
#endif |
63 |
|
64 |
C !INPUT/OUTPUT PARAMETERS: |
65 |
C == Routine arguments == |
66 |
C bi, bj, :: tile indices |
67 |
C uFld,vFld :: Local copy of horizontal velocity field |
68 |
C wFld :: Local copy of vertical velocity field |
69 |
C KappaRk :: Vertical diffusion for Tempertature |
70 |
C myTime :: current time |
71 |
C myIter :: current iteration number |
72 |
C myThid :: my Thread Id. number |
73 |
INTEGER bi, bj |
74 |
_RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
75 |
_RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
76 |
_RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
77 |
_RL KappaRk(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
78 |
_RL myTime |
79 |
INTEGER myIter |
80 |
INTEGER myThid |
81 |
CEOP |
82 |
|
83 |
#ifdef ALLOW_GENERIC_ADVDIFF |
84 |
C !LOCAL VARIABLES: |
85 |
C iMin,iMax :: Range of points for which calculation |
86 |
C jMin,jMax :: results will be set. |
87 |
C k :: vertical index |
88 |
C kM1 :: =k-1 for k>1, =1 for k=1 |
89 |
C kUp :: index into 2 1/2D array, toggles between 1|2 |
90 |
C kDown :: index into 2 1/2D array, toggles between 2|1 |
91 |
C xA :: Tracer cell face area normal to X |
92 |
C yA :: Tracer cell face area normal to X |
93 |
C maskUp :: Land/water mask for Wvel points (interface k) |
94 |
C uTrans :: Zonal volume transport through cell face |
95 |
C vTrans :: Meridional volume transport through cell face |
96 |
C rTrans :: Vertical volume transport at interface k |
97 |
C rTransKp :: Vertical volume transport at inteface k+1 |
98 |
C fVerT :: Flux of temperature (T) in the vertical direction |
99 |
C at the upper(U) and lower(D) faces of a cell. |
100 |
INTEGER iMin, iMax, jMin, jMax |
101 |
INTEGER i, j, k |
102 |
INTEGER kUp, kDown, kM1 |
103 |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
104 |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
105 |
_RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
106 |
_RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
107 |
_RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
108 |
_RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
109 |
_RL rTransKp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
110 |
_RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
111 |
_RL gt_AB (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
112 |
LOGICAL calcAdvection |
113 |
INTEGER iterNb |
114 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
115 |
INTEGER m1, m2 |
116 |
#endif |
117 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
118 |
|
119 |
C- Loop ranges for daughter routines |
120 |
iMin = 1-OLx+2 |
121 |
iMax = sNx+OLx-1 |
122 |
jMin = 1-OLy+2 |
123 |
jMax = sNy+OLy-1 |
124 |
|
125 |
iterNb = myIter |
126 |
IF (staggerTimeStep) iterNb = myIter -1 |
127 |
|
128 |
#ifdef ALLOW_AUTODIFF_TAMC |
129 |
act1 = bi - myBxLo(myThid) |
130 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
131 |
act2 = bj - myByLo(myThid) |
132 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
133 |
act3 = myThid - 1 |
134 |
max3 = nTx*nTy |
135 |
act4 = ikey_dynamics - 1 |
136 |
itdkey = (act1 + 1) + act2*max1 |
137 |
& + act3*max1*max2 |
138 |
& + act4*max1*max2*max3 |
139 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
140 |
|
141 |
C- Tracer tendency needs to be set to zero (moved here from gad_calc_rhs): |
142 |
DO k=1,Nr |
143 |
DO j=1-OLy,sNy+OLy |
144 |
DO i=1-OLx,sNx+OLx |
145 |
gT(i,j,k,bi,bj) = 0. _d 0 |
146 |
ENDDO |
147 |
ENDDO |
148 |
ENDDO |
149 |
DO j=1-OLy,sNy+OLy |
150 |
DO i=1-OLx,sNx+OLx |
151 |
fVerT(i,j,1) = 0. _d 0 |
152 |
fVerT(i,j,2) = 0. _d 0 |
153 |
ENDDO |
154 |
ENDDO |
155 |
#ifdef ALLOW_AUTODIFF |
156 |
DO k=1,Nr |
157 |
DO j=1-OLy,sNy+OLy |
158 |
DO i=1-OLx,sNx+OLx |
159 |
kappaRk(i,j,k) = 0. _d 0 |
160 |
ENDDO |
161 |
ENDDO |
162 |
ENDDO |
163 |
#endif /* ALLOW_AUTODIFF */ |
164 |
|
165 |
#ifdef INCLUDE_CALC_DIFFUSIVITY_CALL |
166 |
IF ( .NOT. implicitDiffusion ) THEN |
167 |
CALL CALC_3D_DIFFUSIVITY( |
168 |
I bi, bj, iMin, iMax, jMin, jMax, |
169 |
I GAD_TEMPERATURE, useGMredi, useKPP, |
170 |
O kappaRk, |
171 |
I myThid ) |
172 |
ENDIF |
173 |
#endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */ |
174 |
|
175 |
#ifndef DISABLE_MULTIDIM_ADVECTION |
176 |
C-- Some advection schemes are better calculated using a multi-dimensional |
177 |
C method in the absence of any other terms and, if used, is done here. |
178 |
C |
179 |
C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h |
180 |
C The default is to use multi-dimensinal advection for non-linear advection |
181 |
C schemes. However, for the sake of efficiency of the adjoint it is necessary |
182 |
C to be able to exclude this scheme to avoid excessive storage and |
183 |
C recomputation. It *is* differentiable, if you need it. |
184 |
C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to |
185 |
C disable this section of code. |
186 |
#ifdef GAD_ALLOW_TS_SOM_ADV |
187 |
# ifdef ALLOW_AUTODIFF_TAMC |
188 |
CADJ STORE som_T = comlev1_bibj, key=itdkey, byte=isbyte |
189 |
# endif |
190 |
IF ( tempSOM_Advection ) THEN |
191 |
# ifdef ALLOW_DEBUG |
192 |
IF (debugMode) CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid) |
193 |
# endif |
194 |
CALL GAD_SOM_ADVECT( |
195 |
I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme, |
196 |
I GAD_TEMPERATURE, dTtracerLev, |
197 |
I uFld, vFld, wFld, theta, |
198 |
U som_T, |
199 |
O gT, |
200 |
I bi, bj, myTime, myIter, myThid ) |
201 |
ELSEIF (tempMultiDimAdvec) THEN |
202 |
#else /* GAD_ALLOW_TS_SOM_ADV */ |
203 |
IF (tempMultiDimAdvec) THEN |
204 |
#endif /* GAD_ALLOW_TS_SOM_ADV */ |
205 |
# ifdef ALLOW_DEBUG |
206 |
IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid) |
207 |
# endif |
208 |
CALL GAD_ADVECTION( |
209 |
I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme, |
210 |
I GAD_TEMPERATURE, dTtracerLev, |
211 |
I uFld, vFld, wFld, theta, |
212 |
O gT, |
213 |
I bi, bj, myTime, myIter, myThid ) |
214 |
ENDIF |
215 |
#endif /* DISABLE_MULTIDIM_ADVECTION */ |
216 |
|
217 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
218 |
|
219 |
C- Start vertical index (k) loop (Nr:1) |
220 |
calcAdvection = tempAdvection .AND. .NOT.tempMultiDimAdvec |
221 |
DO k=Nr,1,-1 |
222 |
#ifdef ALLOW_AUTODIFF_TAMC |
223 |
kkey = (itdkey-1)*Nr + k |
224 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
225 |
kM1 = MAX(1,k-1) |
226 |
kUp = 1+MOD(k+1,2) |
227 |
kDown= 1+MOD(k,2) |
228 |
|
229 |
#ifdef ALLOW_AUTODIFF_TAMC |
230 |
CADJ STORE fVerT(:,:,:) = comlev1_bibj_k, key=kkey, |
231 |
CADJ & byte=isbyte, kind = isbyte |
232 |
CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, |
233 |
CADJ & byte=isbyte, kind = isbyte |
234 |
# ifdef ALLOW_ADAMSBASHFORTH_3 |
235 |
CADJ STORE gtNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey, |
236 |
CADJ & byte=isbyte, kind = isbyte |
237 |
CADJ STORE gtNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey, |
238 |
CADJ & byte=isbyte, kind = isbyte |
239 |
# else |
240 |
CADJ STORE gtNm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, |
241 |
CADJ & byte=isbyte, kind = isbyte |
242 |
# endif |
243 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
244 |
CALL CALC_ADV_FLOW( |
245 |
I uFld, vFld, wFld, |
246 |
U rTrans, |
247 |
O uTrans, vTrans, rTransKp, |
248 |
O maskUp, xA, yA, |
249 |
I k, bi, bj, myThid ) |
250 |
|
251 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
252 |
m1 = 1 + MOD(iterNb+1,2) |
253 |
m2 = 1 + MOD( iterNb ,2) |
254 |
CALL GAD_CALC_RHS( |
255 |
I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown, |
256 |
I xA, yA, maskUp, uFld(1-OLx,1-OLy,k), |
257 |
I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k), |
258 |
I uTrans, vTrans, rTrans, rTransKp, |
259 |
I diffKhT, diffK4T, KappaRk(1-OLx,1-OLy,k), diffKr4T, |
260 |
I gtNm(1-OLx,1-OLy,1,1,1,m2), theta, dTtracerLev, |
261 |
I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme, |
262 |
I calcAdvection, tempImplVertAdv, AdamsBashforth_T, |
263 |
I tempVertDiff4, useGMRedi, useKPP, |
264 |
U fVerT, gT, |
265 |
I myTime, myIter, myThid ) |
266 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
267 |
CALL GAD_CALC_RHS( |
268 |
I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown, |
269 |
I xA, yA, maskUp, uFld(1-OLx,1-OLy,k), |
270 |
I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k), |
271 |
I uTrans, vTrans, rTrans, rTransKp, |
272 |
I diffKhT, diffK4T, KappaRk(1-OLx,1-OLy,k), diffKr4T, |
273 |
I gtNm1, theta, dTtracerLev, |
274 |
I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme, |
275 |
I calcAdvection, tempImplVertAdv, AdamsBashforth_T, |
276 |
I tempVertDiff4, useGMRedi, useKPP, |
277 |
U fVerT, gT, |
278 |
I myTime, myIter, myThid ) |
279 |
#endif |
280 |
|
281 |
C-- External thermal forcing term(s) inside Adams-Bashforth: |
282 |
IF ( tempForcing .AND. tracForcingOutAB.NE.1 ) |
283 |
& CALL EXTERNAL_FORCING_T( |
284 |
I iMin, iMax, jMin, jMax, bi, bj, k, |
285 |
I myTime, myThid ) |
286 |
|
287 |
IF ( AdamsBashforthGt ) THEN |
288 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
289 |
CALL ADAMS_BASHFORTH3( |
290 |
I bi, bj, k, Nr, |
291 |
U gT, gtNm, gt_AB, |
292 |
I tempStartAB, iterNb, myThid ) |
293 |
#else |
294 |
CALL ADAMS_BASHFORTH2( |
295 |
I bi, bj, k, Nr, |
296 |
U gT, gtNm1, gt_AB, |
297 |
I tempStartAB, iterNb, myThid ) |
298 |
#endif |
299 |
#ifdef ALLOW_DIAGNOSTICS |
300 |
IF ( useDiagnostics ) THEN |
301 |
CALL DIAGNOSTICS_FILL(gt_AB,'AB_gT ',k,1,2,bi,bj,myThid) |
302 |
ENDIF |
303 |
#endif /* ALLOW_DIAGNOSTICS */ |
304 |
ENDIF |
305 |
|
306 |
C-- External thermal forcing term(s) outside Adams-Bashforth: |
307 |
IF ( tempForcing .AND. tracForcingOutAB.EQ.1 ) |
308 |
& CALL EXTERNAL_FORCING_T( |
309 |
I iMin, iMax, jMin, jMax, bi, bj, k, |
310 |
I myTime, myThid ) |
311 |
|
312 |
#ifdef NONLIN_FRSURF |
313 |
IF (nonlinFreeSurf.GT.0) THEN |
314 |
CALL FREESURF_RESCALE_G( |
315 |
I bi, bj, k, |
316 |
U gT, |
317 |
I myThid ) |
318 |
IF ( AdamsBashforthGt ) THEN |
319 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
320 |
# ifdef ALLOW_AUTODIFF_TAMC |
321 |
CADJ STORE gtNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey, |
322 |
CADJ & byte=isbyte, kind = isbyte |
323 |
CADJ STORE gtNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey, |
324 |
CADJ & byte=isbyte, kind = isbyte |
325 |
# endif |
326 |
CALL FREESURF_RESCALE_G( |
327 |
I bi, bj, k, |
328 |
U gtNm(1-OLx,1-OLy,1,1,1,1), |
329 |
I myThid ) |
330 |
CALL FREESURF_RESCALE_G( |
331 |
I bi, bj, k, |
332 |
U gtNm(1-OLx,1-OLy,1,1,1,2), |
333 |
I myThid ) |
334 |
#else |
335 |
CALL FREESURF_RESCALE_G( |
336 |
I bi, bj, k, |
337 |
U gtNm1, |
338 |
I myThid ) |
339 |
#endif |
340 |
ENDIF |
341 |
ENDIF |
342 |
#endif /* NONLIN_FRSURF */ |
343 |
|
344 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
345 |
IF ( AdamsBashforth_T ) THEN |
346 |
CALL TIMESTEP_TRACER( |
347 |
I bi, bj, k, dTtracerLev(k), |
348 |
I gtNm(1-OLx,1-OLy,1,1,1,m2), |
349 |
U gT, |
350 |
I myIter, myThid ) |
351 |
ELSE |
352 |
#endif |
353 |
CALL TIMESTEP_TRACER( |
354 |
I bi, bj, k, dTtracerLev(k), |
355 |
I theta, |
356 |
U gT, |
357 |
I myIter, myThid ) |
358 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
359 |
ENDIF |
360 |
#endif |
361 |
|
362 |
C- end of vertical index (k) loop (Nr:1) |
363 |
ENDDO |
364 |
|
365 |
#endif /* ALLOW_GENERIC_ADVDIFF */ |
366 |
|
367 |
RETURN |
368 |
END |