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

Annotation of /MITgcm/model/src/temp_integrate.F

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


Revision 1.1 - (hide annotations) (download)
Tue Nov 19 16:58:38 2013 UTC (10 years, 6 months ago) by jmc
Branch: MAIN
- rename calc_gt.F to temp_integrate.F (includes call to TIMESTEP_TRACER)
 - rename calc_gs.F to salt_integrate.F (includes call to TIMESTEP_TRACER)
 - remove k from thermodynamics.F and move it to temp_integrate.F
   salt_integrate.F and ptracers_integrate.F (now done inside the tracer
   loop).
 - compute locally (in thermodynamics.F) 3-D velocity field that is used to
   advect tracers; pass it as argument to GAD_ADVECTION, GAD_SOM_ADVECT,
   PTRACERS_ADVECTION, TEMP_INTEGRATE, SALT_INTEGRATE, PTRACERS_INTEGRATE,
   GAD_IMPLICIT_R and PTRACERS_IMPLICIT

1 jmc 1.1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gt.F,v 1.60 2013/02/19 13:42:19 jmc Exp $
2     C $Name: $
3    
4     #include "PACKAGES_CONFIG.h"
5     #include "CPP_OPTIONS.h"
6    
7     CBOP
8     C !ROUTINE: TEMP_INTEGRATE
9     C !INTERFACE:
10     SUBROUTINE TEMP_INTEGRATE(
11     I bi, bj, iMin, iMax, jMin, jMax,
12     I uFld, vFld, wFld, KappaRk,
13     I myTime, myIter, myThid )
14     C !DESCRIPTION: \bv
15     C *==========================================================*
16     C | SUBROUTINE TEMP_INTEGRATE
17     C | o Calculate tendency for temperature
18     C | and integrates forward in time.
19     C *==========================================================*
20     C | A procedure called EXTERNAL_FORCING_T is called from
21     C | here. These procedures can be used to add per problem
22     C | heat flux source terms.
23     C | Note: Although it is slightly counter-intuitive the
24     C | EXTERNAL_FORCING routine is not the place to put
25     C | file I/O. Instead files that are required to
26     C | calculate the external source terms are generally
27     C | read during the model main loop. This makes the
28     C | logistics of multi-processing simpler and also
29     C | makes the adjoint generation simpler. It also
30     C | allows for I/O to overlap computation where that
31     C | is supported by hardware.
32     C | Aside from the problem specific term the code here
33     C | forms the tendency terms due to advection and mixing
34     C | The baseline implementation here uses a centered
35     C | difference form for the advection term and a tensorial
36     C | divergence of a flux form for the diffusive term. The
37     C | diffusive term is formulated so that isopycnal mixing and
38     C | GM-style subgrid-scale terms can be incorporated b simply
39     C | setting the diffusion tensor terms appropriately.
40     C *==========================================================*
41     C \ev
42    
43     C !USES:
44     IMPLICIT NONE
45     C == GLobal variables ==
46     #include "SIZE.h"
47     #include "DYNVARS.h"
48     #include "EEPARAMS.h"
49     #include "PARAMS.h"
50     #include "RESTART.h"
51     #ifdef ALLOW_GENERIC_ADVDIFF
52     #include "GAD.h"
53     #endif
54     #ifdef ALLOW_AUTODIFF_TAMC
55     # include "tamc.h"
56     # include "tamc_keys.h"
57     #endif
58    
59     C !INPUT/OUTPUT PARAMETERS:
60     C == Routine arguments ==
61     C bi, bj, :: tile indices
62     C iMin,iMax :: loop range for called routines
63     C jMin,jMax :: loop range for called routines
64     C uFld,vFld :: Local copy of horizontal velocity field
65     C wFld :: Local copy of vertical velocity field
66     C KappaRk :: Vertical diffusion for Tempertature
67     C myTime :: current time
68     C myIter :: current iteration number
69     C myThid :: my Thread Id. number
70     INTEGER bi, bj, iMin, iMax, jMin, jMax
71     _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
72     _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
73     _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
74     _RL KappaRk(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
75     _RL myTime
76     INTEGER myIter
77     INTEGER myThid
78     CEOP
79    
80     #ifdef ALLOW_GENERIC_ADVDIFF
81     C !LOCAL VARIABLES:
82     C k :: vertical index
83     C kM1 :: =k-1 for k>1, =1 for k=1
84     C kUp :: index into 2 1/2D array, toggles between 1|2
85     C kDown :: index into 2 1/2D array, toggles between 2|1
86     C xA :: Tracer cell face area normal to X
87     C yA :: Tracer cell face area normal to X
88     C maskUp :: Land/water mask for Wvel points (interface k)
89     C uTrans :: Zonal volume transport through cell face
90     C vTrans :: Meridional volume transport through cell face
91     C rTrans :: Vertical volume transport at interface k
92     C rTransKp :: Vertical volume transport at inteface k+1
93     C fVerT :: Flux of temperature (T) in the vertical direction
94     C at the upper(U) and lower(D) faces of a cell.
95     INTEGER i, j, k
96     INTEGER kUp, kDown, kM1
97     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103     _RL rTransKp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
105     _RL gt_AB (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106     LOGICAL calcAdvection
107     INTEGER iterNb
108     #ifdef ALLOW_ADAMSBASHFORTH_3
109     INTEGER m1, m2
110     #endif
111    
112     #ifdef ALLOW_AUTODIFF_TAMC
113     act1 = bi - myBxLo(myThid)
114     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
115     act2 = bj - myByLo(myThid)
116     max2 = myByHi(myThid) - myByLo(myThid) + 1
117     act3 = myThid - 1
118     max3 = nTx*nTy
119     act4 = ikey_dynamics - 1
120     itdkey = (act1 + 1) + act2*max1
121     & + act3*max1*max2
122     & + act4*max1*max2*max3
123     #endif /* ALLOW_AUTODIFF_TAMC */
124    
125     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
126    
127     calcAdvection = tempAdvection .AND. .NOT.tempMultiDimAdvec
128     iterNb = myIter
129     IF (staggerTimeStep) iterNb = myIter -1
130    
131     DO j=1-OLy,sNy+OLy
132     DO i=1-OLx,sNx+OLx
133     fVerT(i,j,1) = 0. _d 0
134     fVerT(i,j,2) = 0. _d 0
135     ENDDO
136     ENDDO
137    
138     DO k=Nr,1,-1
139     #ifdef ALLOW_AUTODIFF_TAMC
140     kkey = (itdkey-1)*Nr + k
141     #endif /* ALLOW_AUTODIFF_TAMC */
142     kM1 = MAX(1,k-1)
143     kUp = 1+MOD(k+1,2)
144     kDown= 1+MOD(k,2)
145    
146     #ifdef ALLOW_AUTODIFF_TAMC
147     CADJ STORE rtrans(:,:) = comlev1_bibj_k, key=kkey,
148     CADJ & byte=isbyte, kind = isbyte
149     CADJ STORE fVerT(:,:,:) = comlev1_bibj_k, key=kkey,
150     CADJ & byte=isbyte, kind = isbyte
151     CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey,
152     CADJ & byte=isbyte, kind = isbyte
153     # ifdef ALLOW_ADAMSBASHFORTH_3
154     CADJ STORE gtNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey,
155     CADJ & byte=isbyte, kind = isbyte
156     CADJ STORE gtNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey,
157     CADJ & byte=isbyte, kind = isbyte
158     # else
159     CADJ STORE gtNm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey,
160     CADJ & byte=isbyte, kind = isbyte
161     # endif
162     #endif /* ALLOW_AUTODIFF_TAMC */
163     CALL CALC_ADV_FLOW(
164     I uFld, vFld, wFld,
165     U rTrans,
166     O uTrans, vTrans, rTransKp,
167     O maskUp, xA, yA,
168     I k, bi, bj, myThid )
169    
170     #ifdef ALLOW_ADAMSBASHFORTH_3
171     m1 = 1 + MOD(iterNb+1,2)
172     m2 = 1 + MOD( iterNb ,2)
173     CALL GAD_CALC_RHS(
174     I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown,
175     I xA, yA, maskUp, uFld(1-OLx,1-OLy,k),
176     I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k),
177     I uTrans, vTrans, rTrans, rTransKp,
178     I diffKhT, diffK4T, KappaRk(1-OLx,1-OLy,k), diffKr4T,
179     I gtNm(1-OLx,1-OLy,1,1,1,m2), theta, dTtracerLev,
180     I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme,
181     I calcAdvection, tempImplVertAdv, AdamsBashforth_T,
182     I tempVertDiff4, useGMRedi, useKPP,
183     U fVerT, gT,
184     I myTime, myIter, myThid )
185     #else /* ALLOW_ADAMSBASHFORTH_3 */
186     CALL GAD_CALC_RHS(
187     I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown,
188     I xA, yA, maskUp, uFld(1-OLx,1-OLy,k),
189     I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k),
190     I uTrans, vTrans, rTrans, rTransKp,
191     I diffKhT, diffK4T, KappaRk(1-OLx,1-OLy,k), diffKr4T,
192     I gtNm1, theta, dTtracerLev,
193     I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme,
194     I calcAdvection, tempImplVertAdv, AdamsBashforth_T,
195     I tempVertDiff4, useGMRedi, useKPP,
196     U fVerT, gT,
197     I myTime, myIter, myThid )
198     #endif
199    
200     C-- External thermal forcing term(s) inside Adams-Bashforth:
201     IF ( tempForcing .AND. tracForcingOutAB.NE.1 )
202     & CALL EXTERNAL_FORCING_T(
203     I iMin, iMax, jMin, jMax, bi, bj, k,
204     I myTime, myThid )
205    
206     IF ( AdamsBashforthGt ) THEN
207     #ifdef ALLOW_ADAMSBASHFORTH_3
208     CALL ADAMS_BASHFORTH3(
209     I bi, bj, k, Nr,
210     U gT, gtNm, gt_AB,
211     I tempStartAB, iterNb, myThid )
212     #else
213     CALL ADAMS_BASHFORTH2(
214     I bi, bj, k, Nr,
215     U gT, gtNm1, gt_AB,
216     I tempStartAB, iterNb, myThid )
217     #endif
218     #ifdef ALLOW_DIAGNOSTICS
219     IF ( useDiagnostics ) THEN
220     CALL DIAGNOSTICS_FILL(gt_AB,'AB_gT ',k,1,2,bi,bj,myThid)
221     ENDIF
222     #endif /* ALLOW_DIAGNOSTICS */
223     ENDIF
224    
225     C-- External thermal forcing term(s) outside Adams-Bashforth:
226     IF ( tempForcing .AND. tracForcingOutAB.EQ.1 )
227     & CALL EXTERNAL_FORCING_T(
228     I iMin, iMax, jMin, jMax, bi, bj, k,
229     I myTime, myThid )
230    
231     #ifdef NONLIN_FRSURF
232     IF (nonlinFreeSurf.GT.0) THEN
233     CALL FREESURF_RESCALE_G(
234     I bi, bj, k,
235     U gT,
236     I myThid )
237     IF ( AdamsBashforthGt ) THEN
238     #ifdef ALLOW_ADAMSBASHFORTH_3
239     # ifdef ALLOW_AUTODIFF_TAMC
240     CADJ STORE gtNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey,
241     CADJ & byte=isbyte, kind = isbyte
242     CADJ STORE gtNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey,
243     CADJ & byte=isbyte, kind = isbyte
244     # endif
245     CALL FREESURF_RESCALE_G(
246     I bi, bj, k,
247     U gtNm(1-OLx,1-OLy,1,1,1,1),
248     I myThid )
249     CALL FREESURF_RESCALE_G(
250     I bi, bj, k,
251     U gtNm(1-OLx,1-OLy,1,1,1,2),
252     I myThid )
253     #else
254     CALL FREESURF_RESCALE_G(
255     I bi, bj, k,
256     U gtNm1,
257     I myThid )
258     #endif
259     ENDIF
260     ENDIF
261     #endif /* NONLIN_FRSURF */
262    
263     #ifdef ALLOW_ADAMSBASHFORTH_3
264     IF ( AdamsBashforth_T ) THEN
265     CALL TIMESTEP_TRACER(
266     I bi, bj, k, dTtracerLev(k),
267     I gtNm(1-OLx,1-OLy,1,1,1,m2),
268     U gT,
269     I myIter, myThid )
270     ELSE
271     #endif
272     CALL TIMESTEP_TRACER(
273     I bi, bj, k, dTtracerLev(k),
274     I theta,
275     U gT,
276     I myIter, myThid )
277     #ifdef ALLOW_ADAMSBASHFORTH_3
278     ENDIF
279     #endif
280    
281     C- end of vertical index (k) loop (Nr:1)
282     ENDDO
283    
284     #endif /* ALLOW_GENERIC_ADVDIFF */
285    
286     RETURN
287     END

  ViewVC Help
Powered by ViewVC 1.1.22