/[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.2 - (show annotations) (download)
Mon Aug 13 18:05:26 2001 UTC (22 years, 9 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint40pre7, checkpoint40pre6
Changes since 1.1: +4 -19 lines
Modifications related to split into thermodynamics.F, dynamics.F
o missing initialisations in dynamics.F added
o some fields no longer needed in dynamics/thermodynamics deleted
o split of calc_diffusivity.F into calc_viscosity.F
  (plus split of kpp_calc_diff.F into kpp_calc_visc.F)
o Modifications of some store directives for TAF

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/thermodynamics.F,v 1.1 2001/08/03 19:06:11 adcroft Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
7 C /==========================================================\
8 C | SUBROUTINE THERMODYNAMICS |
9 C | o Controlling routine for the prognostic part of the |
10 C | thermo-dynamics. |
11 C |==========================================================|
12 C \==========================================================/
13 IMPLICIT NONE
14
15 C == Global variables ===
16 #include "SIZE.h"
17 #include "EEPARAMS.h"
18 #include "PARAMS.h"
19 #include "DYNVARS.h"
20 #include "GRID.h"
21 #ifdef ALLOW_PASSIVE_TRACER
22 #include "TR1.h"
23 #endif
24
25 #ifdef ALLOW_AUTODIFF_TAMC
26 # include "tamc.h"
27 # include "tamc_keys.h"
28 # include "FFIELDS.h"
29 # ifdef ALLOW_KPP
30 # include "KPP.h"
31 # endif
32 # ifdef ALLOW_GMREDI
33 # include "GMREDI.h"
34 # endif
35 #endif /* ALLOW_AUTODIFF_TAMC */
36
37 #ifdef ALLOW_TIMEAVE
38 #include "TIMEAVE_STATV.h"
39 #endif
40
41 C == Routine arguments ==
42 C myTime - Current time in simulation
43 C myIter - Current iteration number in simulation
44 C myThid - Thread number for this instance of the routine.
45 _RL myTime
46 INTEGER myIter
47 INTEGER myThid
48
49 C == Local variables
50 C xA, yA - Per block temporaries holding face areas
51 C uTrans, vTrans, rTrans - Per block temporaries holding flow
52 C transport
53 C o uTrans: Zonal transport
54 C o vTrans: Meridional transport
55 C o rTrans: Vertical transport
56 C maskUp o maskUp: land/water mask for W points
57 C fVer[STUV] o fVer: Vertical flux term - note fVer
58 C is "pipelined" in the vertical
59 C so we need an fVer for each
60 C variable.
61 C rhoK, rhoKM1 - Density at current level, and level above
62 C phiHyd - Hydrostatic part of the potential phiHydi.
63 C In z coords phiHydiHyd is the hydrostatic
64 C Potential (=pressure/rho0) anomaly
65 C In p coords phiHydiHyd is the geopotential
66 C surface height anomaly.
67 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
68 C phiSurfY or geopotentiel (atmos) in X and Y direction
69 C KappaRT, - Total diffusion in vertical for T and S.
70 C KappaRS (background + spatially varying, isopycnal term).
71 C iMin, iMax - Ranges and sub-block indices on which calculations
72 C jMin, jMax are applied.
73 C bi, bj
74 C k, kup, - Index for layer above and below. kup and kDown
75 C kDown, km1 are switched with layer to be the appropriate
76 C index into fVerTerm.
77 C tauAB - Adams-Bashforth timestepping weight: 0=forward ; 1/2=Adams-Bashf.
78 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
85 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
86 _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
87 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
88 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
93 _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
94 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
95 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
96 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
97 _RL tauAB
98
99 C This is currently used by IVDC and Diagnostics
100 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
101
102 INTEGER iMin, iMax
103 INTEGER jMin, jMax
104 INTEGER bi, bj
105 INTEGER i, j
106 INTEGER k, km1, kup, kDown
107
108 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
109 c CHARACTER*(MAX_LEN_MBUF) suff
110 c LOGICAL DIFFERENT_MULTIPLE
111 c EXTERNAL DIFFERENT_MULTIPLE
112 Cjmc(end)
113
114 C--- The algorithm...
115 C
116 C "Correction Step"
117 C =================
118 C Here we update the horizontal velocities with the surface
119 C pressure such that the resulting flow is either consistent
120 C with the free-surface evolution or the rigid-lid:
121 C U[n] = U* + dt x d/dx P
122 C V[n] = V* + dt x d/dy P
123 C
124 C "Calculation of Gs"
125 C ===================
126 C This is where all the accelerations and tendencies (ie.
127 C physics, parameterizations etc...) are calculated
128 C rho = rho ( theta[n], salt[n] )
129 C b = b(rho, theta)
130 C K31 = K31 ( rho )
131 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
132 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
133 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
134 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
135 C
136 C "Time-stepping" or "Prediction"
137 C ================================
138 C The models variables are stepped forward with the appropriate
139 C time-stepping scheme (currently we use Adams-Bashforth II)
140 C - For momentum, the result is always *only* a "prediction"
141 C in that the flow may be divergent and will be "corrected"
142 C later with a surface pressure gradient.
143 C - Normally for tracers the result is the new field at time
144 C level [n+1} *BUT* in the case of implicit diffusion the result
145 C is also *only* a prediction.
146 C - We denote "predictors" with an asterisk (*).
147 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
148 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
149 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
150 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
151 C With implicit diffusion:
152 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
153 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
154 C (1 + dt * K * d_zz) theta[n] = theta*
155 C (1 + dt * K * d_zz) salt[n] = salt*
156 C---
157
158 #ifdef ALLOW_AUTODIFF_TAMC
159 C-- dummy statement to end declaration part
160 ikey = 1
161 #endif /* ALLOW_AUTODIFF_TAMC */
162
163 C-- Set up work arrays with valid (i.e. not NaN) values
164 C These inital values do not alter the numerical results. They
165 C just ensure that all memory references are to valid floating
166 C point numbers. This prevents spurious hardware signals due to
167 C uninitialised but inert locations.
168 DO j=1-OLy,sNy+OLy
169 DO i=1-OLx,sNx+OLx
170 xA(i,j) = 0. _d 0
171 yA(i,j) = 0. _d 0
172 uTrans(i,j) = 0. _d 0
173 vTrans(i,j) = 0. _d 0
174 DO k=1,Nr
175 phiHyd(i,j,k) = 0. _d 0
176 sigmaX(i,j,k) = 0. _d 0
177 sigmaY(i,j,k) = 0. _d 0
178 sigmaR(i,j,k) = 0. _d 0
179 ENDDO
180 rhoKM1 (i,j) = 0. _d 0
181 rhok (i,j) = 0. _d 0
182 phiSurfX(i,j) = 0. _d 0
183 phiSurfY(i,j) = 0. _d 0
184 ENDDO
185 ENDDO
186
187
188 #ifdef ALLOW_AUTODIFF_TAMC
189 C-- HPF directive to help TAMC
190 CHPF$ INDEPENDENT
191 #endif /* ALLOW_AUTODIFF_TAMC */
192
193 DO bj=myByLo(myThid),myByHi(myThid)
194
195 #ifdef ALLOW_AUTODIFF_TAMC
196 C-- HPF directive to help TAMC
197 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
198 CHPF$& ,phiHyd,utrans,vtrans,xA,yA
199 CHPF$& ,KappaRT,KappaRS
200 CHPF$& )
201 #endif /* ALLOW_AUTODIFF_TAMC */
202
203 DO bi=myBxLo(myThid),myBxHi(myThid)
204
205 #ifdef ALLOW_AUTODIFF_TAMC
206 act1 = bi - myBxLo(myThid)
207 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
208
209 act2 = bj - myByLo(myThid)
210 max2 = myByHi(myThid) - myByLo(myThid) + 1
211
212 act3 = myThid - 1
213 max3 = nTx*nTy
214
215 act4 = ikey_dynamics - 1
216
217 ikey = (act1 + 1) + act2*max1
218 & + act3*max1*max2
219 & + act4*max1*max2*max3
220 #endif /* ALLOW_AUTODIFF_TAMC */
221
222 C-- Set up work arrays that need valid initial values
223 DO j=1-OLy,sNy+OLy
224 DO i=1-OLx,sNx+OLx
225 rTrans (i,j) = 0. _d 0
226 fVerT (i,j,1) = 0. _d 0
227 fVerT (i,j,2) = 0. _d 0
228 fVerS (i,j,1) = 0. _d 0
229 fVerS (i,j,2) = 0. _d 0
230 fVerTr1(i,j,1) = 0. _d 0
231 fVerTr1(i,j,2) = 0. _d 0
232 ENDDO
233 ENDDO
234
235 DO k=1,Nr
236 DO j=1-OLy,sNy+OLy
237 DO i=1-OLx,sNx+OLx
238 C This is currently also used by IVDC and Diagnostics
239 ConvectCount(i,j,k) = 0.
240 KappaRT(i,j,k) = 0. _d 0
241 KappaRS(i,j,k) = 0. _d 0
242 ENDDO
243 ENDDO
244 ENDDO
245
246 iMin = 1-OLx+1
247 iMax = sNx+OLx
248 jMin = 1-OLy+1
249 jMax = sNy+OLy
250
251
252 #ifdef ALLOW_AUTODIFF_TAMC
253 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
254 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
255 #endif /* ALLOW_AUTODIFF_TAMC */
256
257 C-- Start of diagnostic loop
258 DO k=Nr,1,-1
259
260 #ifdef ALLOW_AUTODIFF_TAMC
261 C? Patrick, is this formula correct now that we change the loop range?
262 C? Do we still need this?
263 cph kkey formula corrected.
264 cph Needed for rhok, rhokm1, in the case useGMREDI.
265 kkey = (ikey-1)*Nr + k
266 CADJ STORE rhokm1(:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
267 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
268 #endif /* ALLOW_AUTODIFF_TAMC */
269
270 C-- Integrate continuity vertically for vertical velocity
271 CALL INTEGRATE_FOR_W(
272 I bi, bj, k, uVel, vVel,
273 O wVel,
274 I myThid )
275
276 #ifdef ALLOW_OBCS
277 #ifdef ALLOW_NONHYDROSTATIC
278 C-- Apply OBC to W if in N-H mode
279 IF (useOBCS.AND.nonHydrostatic) THEN
280 CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
281 ENDIF
282 #endif /* ALLOW_NONHYDROSTATIC */
283 #endif /* ALLOW_OBCS */
284
285 C-- Calculate gradients of potential density for isoneutral
286 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
287 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
288 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
289 #ifdef ALLOW_AUTODIFF_TAMC
290 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
291 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
292 #endif /* ALLOW_AUTODIFF_TAMC */
293 CALL FIND_RHO(
294 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
295 I theta, salt,
296 O rhoK,
297 I myThid )
298 IF (k.GT.1) THEN
299 #ifdef ALLOW_AUTODIFF_TAMC
300 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
301 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
302 #endif /* ALLOW_AUTODIFF_TAMC */
303 CALL FIND_RHO(
304 I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
305 I theta, salt,
306 O rhoKm1,
307 I myThid )
308 ENDIF
309 CALL GRAD_SIGMA(
310 I bi, bj, iMin, iMax, jMin, jMax, k,
311 I rhoK, rhoKm1, rhoK,
312 O sigmaX, sigmaY, sigmaR,
313 I myThid )
314 ENDIF
315
316 C-- Implicit Vertical Diffusion for Convection
317 c ==> should use sigmaR !!!
318 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
319 CALL CALC_IVDC(
320 I bi, bj, iMin, iMax, jMin, jMax, k,
321 I rhoKm1, rhoK,
322 U ConvectCount, KappaRT, KappaRS,
323 I myTime, myIter, myThid)
324 ENDIF
325
326 C-- end of diagnostic k loop (Nr:1)
327 ENDDO
328
329 #ifdef ALLOW_AUTODIFF_TAMC
330 cph avoids recomputation of integrate_for_w
331 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
332 #endif /* ALLOW_AUTODIFF_TAMC */
333
334 #ifdef ALLOW_OBCS
335 C-- Calculate future values on open boundaries
336 IF (useOBCS) THEN
337 CALL OBCS_CALC( bi, bj, myTime+deltaT,
338 I uVel, vVel, wVel, theta, salt,
339 I myThid )
340 ENDIF
341 #endif /* ALLOW_OBCS */
342
343 C-- Determines forcing terms based on external fields
344 C relaxation terms, etc.
345 CALL EXTERNAL_FORCING_SURF(
346 I bi, bj, iMin, iMax, jMin, jMax,
347 I myThid )
348 #ifdef ALLOW_AUTODIFF_TAMC
349 cph needed for KPP
350 CADJ STORE surfacetendencyU(:,:,bi,bj)
351 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
352 CADJ STORE surfacetendencyV(:,:,bi,bj)
353 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
354 CADJ STORE surfacetendencyS(:,:,bi,bj)
355 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
356 CADJ STORE surfacetendencyT(:,:,bi,bj)
357 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
358 #endif /* ALLOW_AUTODIFF_TAMC */
359
360 #ifdef ALLOW_GMREDI
361
362 #ifdef ALLOW_AUTODIFF_TAMC
363 CADJ STORE sigmaX(:,:,:) = comlev1, key=ikey, byte=isbyte
364 CADJ STORE sigmaY(:,:,:) = comlev1, key=ikey, byte=isbyte
365 CADJ STORE sigmaR(:,:,:) = comlev1, key=ikey, byte=isbyte
366 #endif /* ALLOW_AUTODIFF_TAMC */
367 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
368 IF (useGMRedi) THEN
369 DO k=1,Nr
370 CALL GMREDI_CALC_TENSOR(
371 I bi, bj, iMin, iMax, jMin, jMax, k,
372 I sigmaX, sigmaY, sigmaR,
373 I myThid )
374 ENDDO
375 #ifdef ALLOW_AUTODIFF_TAMC
376 ELSE
377 DO k=1, Nr
378 CALL GMREDI_CALC_TENSOR_DUMMY(
379 I bi, bj, iMin, iMax, jMin, jMax, k,
380 I sigmaX, sigmaY, sigmaR,
381 I myThid )
382 ENDDO
383 #endif /* ALLOW_AUTODIFF_TAMC */
384 ENDIF
385
386 #ifdef ALLOW_AUTODIFF_TAMC
387 CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
388 CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
389 CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
390 #endif /* ALLOW_AUTODIFF_TAMC */
391
392 #endif /* ALLOW_GMREDI */
393
394 #ifdef ALLOW_KPP
395 C-- Compute KPP mixing coefficients
396 IF (useKPP) THEN
397 CALL KPP_CALC(
398 I bi, bj, myTime, myThid )
399 #ifdef ALLOW_AUTODIFF_TAMC
400 ELSE
401 CALL KPP_CALC_DUMMY(
402 I bi, bj, myTime, myThid )
403 #endif /* ALLOW_AUTODIFF_TAMC */
404 ENDIF
405
406 #ifdef ALLOW_AUTODIFF_TAMC
407 CADJ STORE KPPghat (:,:,:,bi,bj)
408 CADJ & , KPPviscAz (:,:,:,bi,bj)
409 CADJ & , KPPdiffKzT(:,:,:,bi,bj)
410 CADJ & , KPPdiffKzS(:,:,:,bi,bj)
411 CADJ & , KPPfrac (:,: ,bi,bj)
412 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
413 #endif /* ALLOW_AUTODIFF_TAMC */
414
415 #endif /* ALLOW_KPP */
416
417 #ifdef ALLOW_AUTODIFF_TAMC
418 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
419 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte
420 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
421 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
422 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
423 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
424 #ifdef ALLOW_PASSIVE_TRACER
425 CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
426 #endif
427 #endif /* ALLOW_AUTODIFF_TAMC */
428
429 #ifdef ALLOW_AIM
430 C AIM - atmospheric intermediate model, physics package code.
431 C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics
432 IF ( useAIM ) THEN
433 CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
434 CALL AIM_DO_ATMOS_PHYSICS( phiHyd, bi, bj, myTime, myThid )
435 CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid)
436 ENDIF
437 #endif /* ALLOW_AIM */
438
439
440 C-- Start of thermodynamics loop
441 DO k=Nr,1,-1
442 #ifdef ALLOW_AUTODIFF_TAMC
443 C? Patrick Is this formula correct?
444 cph Yes, but I rewrote it.
445 cph Also, the KappaR? need the index and subscript k!
446 kkey = (ikey-1)*Nr + k
447 #endif /* ALLOW_AUTODIFF_TAMC */
448
449 C-- km1 Points to level above k (=k-1)
450 C-- kup Cycles through 1,2 to point to layer above
451 C-- kDown Cycles through 2,1 to point to current layer
452
453 km1 = MAX(1,k-1)
454 kup = 1+MOD(k+1,2)
455 kDown= 1+MOD(k,2)
456
457 iMin = 1-OLx
458 iMax = sNx+OLx
459 jMin = 1-OLy
460 jMax = sNy+OLy
461
462 C-- Get temporary terms used by tendency routines
463 CALL CALC_COMMON_FACTORS (
464 I bi,bj,iMin,iMax,jMin,jMax,k,
465 O xA,yA,uTrans,vTrans,rTrans,maskUp,
466 I myThid)
467
468 #ifdef ALLOW_AUTODIFF_TAMC
469 CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
470 CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
471 #endif /* ALLOW_AUTODIFF_TAMC */
472
473 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
474 C-- Calculate the total vertical diffusivity
475 CALL CALC_DIFFUSIVITY(
476 I bi,bj,iMin,iMax,jMin,jMax,k,
477 I maskUp,
478 O KappaRT,KappaRS,
479 I myThid)
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 gTnm1, gSnm1, etc.
489 IF ( tempStepping ) THEN
490 CALL CALC_GT(
491 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
492 I xA,yA,uTrans,vTrans,rTrans,maskUp,
493 I KappaRT,
494 U fVerT,
495 I myTime, myThid)
496 tauAB = 0.5d0 + abEps
497 CALL TIMESTEP_TRACER(
498 I bi,bj,iMin,iMax,jMin,jMax,k,tauAB,
499 I theta, gT,
500 U gTnm1,
501 I myIter, myThid)
502 ENDIF
503 IF ( saltStepping ) THEN
504 CALL CALC_GS(
505 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
506 I xA,yA,uTrans,vTrans,rTrans,maskUp,
507 I KappaRS,
508 U fVerS,
509 I myTime, myThid)
510 tauAB = 0.5d0 + abEps
511 CALL TIMESTEP_TRACER(
512 I bi,bj,iMin,iMax,jMin,jMax,k,tauAB,
513 I salt, gS,
514 U gSnm1,
515 I myIter, myThid)
516 ENDIF
517 #ifdef ALLOW_PASSIVE_TRACER
518 IF ( tr1Stepping ) THEN
519 CALL CALC_GTR1(
520 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
521 I xA,yA,uTrans,vTrans,rTrans,maskUp,
522 I KappaRT,
523 U fVerTr1,
524 I myTime, myThid)
525 tauAB = 0.5d0 + abEps
526 CALL TIMESTEP_TRACER(
527 I bi,bj,iMin,iMax,jMin,jMax,k,tauAB,
528 I Tr1, gTr1,
529 U gTr1NM1,
530 I myIter, myThid)
531 ENDIF
532 #endif
533
534 #ifdef ALLOW_OBCS
535 C-- Apply open boundary conditions
536 IF (useOBCS) THEN
537 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
538 END IF
539 #endif /* ALLOW_OBCS */
540
541 C-- Freeze water
542 IF (allowFreezing) THEN
543 #ifdef ALLOW_AUTODIFF_TAMC
544 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k
545 CADJ & , key = kkey, byte = isbyte
546 #endif /* ALLOW_AUTODIFF_TAMC */
547 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
548 END IF
549
550 C-- end of thermodynamic k loop (Nr:1)
551 ENDDO
552
553
554 #ifdef ALLOW_AUTODIFF_TAMC
555 C? Patrick? What about this one?
556 cph Keys iikey and idkey don't seem to be needed
557 cph since storing occurs on different tape for each
558 cph impldiff call anyways.
559 cph Thus, common block comlev1_impl isn't needed either.
560 cph Storing below needed in the case useGMREDI.
561 iikey = (ikey-1)*maximpl
562 #endif /* ALLOW_AUTODIFF_TAMC */
563
564 C-- Implicit diffusion
565 IF (implicitDiffusion) THEN
566
567 IF (tempStepping) THEN
568 #ifdef ALLOW_AUTODIFF_TAMC
569 idkey = iikey + 1
570 CADJ STORE gTNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
571 #endif /* ALLOW_AUTODIFF_TAMC */
572 CALL IMPLDIFF(
573 I bi, bj, iMin, iMax, jMin, jMax,
574 I deltaTtracer, KappaRT, recip_HFacC,
575 U gTNm1,
576 I myThid )
577 ENDIF
578
579 IF (saltStepping) THEN
580 #ifdef ALLOW_AUTODIFF_TAMC
581 idkey = iikey + 2
582 CADJ STORE gSNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
583 #endif /* ALLOW_AUTODIFF_TAMC */
584 CALL IMPLDIFF(
585 I bi, bj, iMin, iMax, jMin, jMax,
586 I deltaTtracer, KappaRS, recip_HFacC,
587 U gSNm1,
588 I myThid )
589 ENDIF
590
591 #ifdef ALLOW_PASSIVE_TRACER
592 IF (tr1Stepping) THEN
593 #ifdef ALLOW_AUTODIFF_TAMC
594 CADJ STORE gTr1Nm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
595 #endif /* ALLOW_AUTODIFF_TAMC */
596 CALL IMPLDIFF(
597 I bi, bj, iMin, iMax, jMin, jMax,
598 I deltaTtracer, KappaRT, recip_HFacC,
599 U gTr1Nm1,
600 I myThid )
601 ENDIF
602 #endif
603
604 #ifdef ALLOW_OBCS
605 C-- Apply open boundary conditions
606 IF (useOBCS) THEN
607 DO K=1,Nr
608 CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid )
609 ENDDO
610 END IF
611 #endif /* ALLOW_OBCS */
612
613 C-- End If implicitDiffusion
614 ENDIF
615
616 Ccs-
617 ENDDO
618 ENDDO
619
620 #ifdef ALLOW_AIM
621 IF ( useAIM ) THEN
622 CALL AIM_AIM2DYN_EXCHANGES( myTime, myThid )
623 ENDIF
624 _EXCH_XYZ_R8(gTnm1,myThid)
625 _EXCH_XYZ_R8(gSnm1,myThid)
626 #else
627 IF (staggerTimeStep.AND.useCubedSphereExchange) THEN
628 _EXCH_XYZ_R8(gTnm1,myThid)
629 _EXCH_XYZ_R8(gSnm1,myThid)
630 ENDIF
631 #endif /* ALLOW_AIM */
632
633 RETURN
634 END

  ViewVC Help
Powered by ViewVC 1.1.22