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

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

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


Revision 1.28 - (show annotations) (download)
Fri Feb 2 21:04:47 2001 UTC (23 years, 3 months ago) by adcroft
Branch: MAIN
Changes since 1.27: +6 -19 lines
Merged changes from branch "branch-atmos-merge" into MAIN (checkpoint34)
 - substantial modifications to algorithm sequence (dynamics.F)
 - packaged OBCS, Shapiro filter, Zonal filter, Atmospheric Physics

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gt.F,v 1.27.2.2 2001/01/24 21:46:21 adcroft Exp $
2
3 #include "CPP_OPTIONS.h"
4
5 CStartOfInterFace
6 SUBROUTINE CALC_GT(
7 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
8 I xA,yA,uTrans,vTrans,rTrans,maskup,maskC,
9 I KappaRT,
10 U fVerT,
11 I myCurrentTime, myThid )
12 C /==========================================================\
13 C | SUBROUTINE CALC_GT |
14 C | o Calculate the temperature tendency terms. |
15 C |==========================================================|
16 C | A procedure called EXTERNAL_FORCING_T is called from |
17 C | here. These procedures can be used to add per problem |
18 C | heat flux source terms. |
19 C | Note: Although it is slightly counter-intuitive the |
20 C | EXTERNAL_FORCING routine is not the place to put |
21 C | file I/O. Instead files that are required to |
22 C | calculate the external source terms are generally |
23 C | read during the model main loop. This makes the |
24 C | logisitics of multi-processing simpler and also |
25 C | makes the adjoint generation simpler. It also |
26 C | allows for I/O to overlap computation where that |
27 C | is supported by hardware. |
28 C | Aside from the problem specific term the code here |
29 C | forms the tendency terms due to advection and mixing |
30 C | The baseline implementation here uses a centered |
31 C | difference form for the advection term and a tensorial |
32 C | divergence of a flux form for the diffusive term. The |
33 C | diffusive term is formulated so that isopycnal mixing and|
34 C | GM-style subgrid-scale terms can be incorporated b simply|
35 C | setting the diffusion tensor terms appropriately. |
36 C \==========================================================/
37 IMPLICIT NONE
38
39 C == GLobal variables ==
40 #include "SIZE.h"
41 #include "DYNVARS.h"
42 #include "EEPARAMS.h"
43 #include "PARAMS.h"
44 #include "GRID.h"
45 #include "FFIELDS.h"
46 c #include "GM_ARRAYS.h"
47
48
49 C == Routine arguments ==
50 C fVerT - Flux of temperature (T) in the vertical
51 C direction at the upper(U) and lower(D) faces of a cell.
52 C maskUp - Land mask used to denote base of the domain.
53 C maskC - Land mask for theta cells (used in TOP_LAYER only)
54 C xA - Tracer cell face area normal to X
55 C yA - Tracer cell face area normal to X
56 C uTrans - Zonal volume transport through cell face
57 C vTrans - Meridional volume transport through cell face
58 C rTrans - Vertical volume transport through cell face
59 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
60 C results will be set.
61 C myThid - Instance number for this innvocation of CALC_GT
62 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
63 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
64 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
65 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
66 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
67 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
68 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
69 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70 _RL KappaRT(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
71 INTEGER k,kUp,kDown,kM1
72 INTEGER bi,bj,iMin,iMax,jMin,jMax
73 INTEGER myThid
74 _RL myCurrentTime
75 CEndOfInterface
76
77 C == Local variables ==
78 C I, J, K - Loop counters
79 INTEGER i,j
80 LOGICAL TOP_LAYER
81 _RL afFacT, dfFacT
82 _RL dTdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RL dTdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL df4 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88 _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89
90 #ifdef ALLOW_AUTODIFF_TAMC
91 C-- only the kUp part of fverT is set in this subroutine
92 C-- the kDown is still required
93
94 fVerT(1,1,kDown) = fVerT(1,1,kDown)
95 DO j=1-OLy,sNy+OLy
96 DO i=1-OLx,sNx+OLx
97 fZon(i,j) = 0.0
98 fMer(i,j) = 0.0
99 fVerT(i,j,kUp) = 0.0
100 ENDDO
101 ENDDO
102 #endif
103
104 afFacT = 1. _d 0
105 dfFacT = 1. _d 0
106 TOP_LAYER = K .EQ. 1
107
108 C--- Calculate advective and diffusive fluxes between cells.
109
110 #ifdef INCLUDE_T_DIFFUSION_CODE
111 C o Zonal tracer gradient
112 DO j=1-Oly,sNy+Oly
113 DO i=1-Olx+1,sNx+Olx
114 dTdx(i,j) = _recip_dxC(i,j,bi,bj)*
115 & (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj))
116 ENDDO
117 ENDDO
118 C o Meridional tracer gradient
119 DO j=1-Oly+1,sNy+Oly
120 DO i=1-Olx,sNx+Olx
121 dTdy(i,j) = _recip_dyC(i,j,bi,bj)*
122 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
123 ENDDO
124 ENDDO
125
126 C-- del^2 of T, needed for bi-harmonic (del^4) term
127 IF (diffK4T .NE. 0.) THEN
128 DO j=1-Oly+1,sNy+Oly-1
129 DO i=1-Olx+1,sNx+Olx-1
130 df4(i,j)= _recip_hFacC(i,j,k,bi,bj)
131 & *recip_drF(k)/_rA(i,j,bi,bj)
132 & *(
133 & +( xA(i+1,j)*dTdx(i+1,j)-xA(i,j)*dTdx(i,j) )
134 & +( yA(i,j+1)*dTdy(i,j+1)-yA(i,j)*dTdy(i,j) )
135 & )
136 ENDDO
137 ENDDO
138 ENDIF
139 #endif
140
141 C-- Zonal flux (fZon is at west face of "theta" cell)
142 #ifdef INCLUDE_T_ADVECTION_CODE
143 C o Advective component of zonal flux
144 DO j=jMin,jMax
145 DO i=iMin,iMax
146 af(i,j) =
147 & uTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i-1,j,k,bi,bj))*0.5 _d 0
148 ENDDO
149 ENDDO
150 #endif /* INCLUDE_T_ADVECTION_CODE */
151 #ifdef INCLUDE_T_DIFFUSION_CODE
152 C o Diffusive component of zonal flux
153 DO j=jMin,jMax
154 DO i=iMin,iMax
155 df(i,j) = -diffKhT*xA(i,j)*dTdx(i,j)
156 ENDDO
157 ENDDO
158 #ifdef ALLOW_GMREDI
159 IF (useGMRedi) CALL GMREDI_XTRANSPORT(
160 I iMin,iMax,jMin,jMax,bi,bj,K,
161 I xA,theta,
162 U df,
163 I myThid)
164 #endif
165 C o Add the bi-harmonic contribution
166 IF (diffK4T .NE. 0.) THEN
167 DO j=jMin,jMax
168 DO i=iMin,iMax
169 df(i,j) = df(i,j) + xA(i,j)*
170 & diffK4T*(df4(i,j)-df4(i-1,j))*_recip_dxC(i,j,bi,bj)
171 ENDDO
172 ENDDO
173 ENDIF
174 #endif /* INCLUDE_T_DIFFUSION_CODE */
175 C o Net zonal flux
176 DO j=jMin,jMax
177 DO i=iMin,iMax
178 fZon(i,j) = 0.
179 & _ADT( + afFacT*af(i,j) )
180 & _LPT( + dfFacT*df(i,j) )
181 ENDDO
182 ENDDO
183
184 C-- Meridional flux (fMer is at south face of "theta" cell)
185 #ifdef INCLUDE_T_ADVECTION_CODE
186 C o Advective component of meridional flux
187 DO j=jMin,jMax
188 DO i=iMin,iMax
189 af(i,j) =
190 & vTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j-1,k,bi,bj))*0.5 _d 0
191 ENDDO
192 ENDDO
193 #endif /* INCLUDE_T_ADVECTION_CODE */
194 #ifdef INCLUDE_T_DIFFUSION_CODE
195 C o Diffusive component of meridional flux
196 DO j=jMin,jMax
197 DO i=iMin,iMax
198 df(i,j) = -diffKhT*yA(i,j)*dTdy(i,j)
199 ENDDO
200 ENDDO
201 #ifdef ALLOW_GMREDI
202 IF (useGMRedi) CALL GMREDI_YTRANSPORT(
203 I iMin,iMax,jMin,jMax,bi,bj,K,
204 I yA,theta,
205 U df,
206 I myThid)
207 #endif
208 C o Add the bi-harmonic contribution
209 IF (diffK4T .NE. 0.) THEN
210 DO j=jMin,jMax
211 DO i=iMin,iMax
212 df(i,j) = df(i,j) + yA(i,j)*
213 & diffK4T*(df4(i,j)-df4(i,j-1))*_recip_dyC(i,j,bi,bj)
214 ENDDO
215 ENDDO
216 ENDIF
217 #endif /* INCLUDE_T_DIFFUSION_CODE */
218 C o Net meridional flux
219 DO j=jMin,jMax
220 DO i=iMin,iMax
221 fMer(i,j) = 0.
222 & _ADT( + afFacT*af(i,j) )
223 & _LPT( + dfFacT*df(i,j) )
224 ENDDO
225 ENDDO
226
227 #ifdef INCLUDE_T_DIFFUSION_CODE
228 C-- Terms that diffusion tensor projects onto z
229 DO j=jMin,jMax
230 DO i=iMin,iMax
231 dTdx(i,j) = 0.5*(
232 & +0.5*(_maskW(i+1,j,k,bi,bj)
233 & *_recip_dxC(i+1,j,bi,bj)*
234 & (theta(i+1,j,k,bi,bj)-theta(i,j,k,bi,bj))
235 & +_maskW(i,j,k,bi,bj)
236 & *_recip_dxC(i,j,bi,bj)*
237 & (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj)))
238 & +0.5*(_maskW(i+1,j,km1,bi,bj)
239 & *_recip_dxC(i+1,j,bi,bj)*
240 & (theta(i+1,j,km1,bi,bj)-theta(i,j,km1,bi,bj))
241 & +_maskW(i,j,km1,bi,bj)
242 & *_recip_dxC(i,j,bi,bj)*
243 & (theta(i,j,km1,bi,bj)-theta(i-1,j,km1,bi,bj)))
244 & )
245 ENDDO
246 ENDDO
247 DO j=jMin,jMax
248 DO i=iMin,iMax
249 dTdy(i,j) = 0.5*(
250 & +0.5*(_maskS(i,j,k,bi,bj)
251 & *_recip_dyC(i,j,bi,bj)*
252 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
253 & +_maskS(i,j+1,k,bi,bj)
254 & *_recip_dyC(i,j+1,bi,bj)*
255 & (theta(i,j+1,k,bi,bj)-theta(i,j,k,bi,bj)))
256 & +0.5*(_maskS(i,j,km1,bi,bj)
257 & *_recip_dyC(i,j,bi,bj)*
258 & (theta(i,j,km1,bi,bj)-theta(i,j-1,km1,bi,bj))
259 & +_maskS(i,j+1,km1,bi,bj)
260 & *_recip_dyC(i,j+1,bi,bj)*
261 & (theta(i,j+1,km1,bi,bj)-theta(i,j,km1,bi,bj)))
262 & )
263 ENDDO
264 ENDDO
265 #endif /* INCLUDE_T_DIFFUSION_CODE */
266
267 C-- Vertical flux ( fVerT(,,kUp) is at upper face of "theta" cell )
268 #ifdef INCLUDE_T_ADVECTION_CODE
269 C o Advective component of vertical flux
270 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
271 C (this plays the role of the free-surface correction)
272 DO j=jMin,jMax
273 DO i=iMin,iMax
274 af(i,j) =
275 & rTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0
276 ENDDO
277 ENDDO
278 #endif /* INCLUDE_T_ADVECTION_CODE */
279 #ifdef INCLUDE_T_DIFFUSION_CODE
280 C o Diffusive component of vertical flux
281 C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
282 C boundary condition.
283 IF (implicitDiffusion) THEN
284 DO j=jMin,jMax
285 DO i=iMin,iMax
286 df(i,j) = 0.
287 ENDDO
288 ENDDO
289 ELSE
290 DO j=jMin,jMax
291 DO i=iMin,iMax
292 df(i,j) = - _rA(i,j,bi,bj)*(
293 & KappaRT(i,j,k)*recip_drC(k)
294 & *(theta(i,j,kM1,bi,bj)-theta(i,j,k,bi,bj))*rkFac
295 & )
296 ENDDO
297 ENDDO
298 ENDIF
299 #endif /* INCLUDE_T_DIFFUSION_CODE */
300
301 #ifdef ALLOW_GMREDI
302 IF (useGMRedi) CALL GMREDI_RTRANSPORT(
303 I iMin,iMax,jMin,jMax,bi,bj,K,
304 I maskUp,theta,
305 U df,
306 I myThid)
307 #endif
308
309 #ifdef ALLOW_KPP
310 C-- Add non local KPP transport term (ghat) to diffusive T flux.
311 IF (useKPP) CALL KPP_TRANSPORT_T(
312 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
313 I maskC,KappaRT,
314 U df )
315 #endif
316
317 C o Net vertical flux
318 DO j=jMin,jMax
319 DO i=iMin,iMax
320 fVerT(i,j,kUp) = 0.
321 & _ADT( +afFacT*af(i,j)*maskUp(i,j) )
322 & _LPT( +dfFacT*df(i,j)*maskUp(i,j) )
323 ENDDO
324 ENDDO
325 #ifdef INCLUDE_T_ADVECTION_CODE
326 IF ( TOP_LAYER ) THEN
327 DO j=jMin,jMax
328 DO i=iMin,iMax
329 fVerT(i,j,kUp) = afFacT*af(i,j)*freeSurfFac
330 ENDDO
331 ENDDO
332 ENDIF
333 #endif /* INCLUDE_T_ADVECTION_CODE */
334
335 C-- Tendency is minus divergence of the fluxes.
336 C Note. Tendency terms will only be correct for range
337 C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
338 C will contain valid floating point numbers but
339 C they are not algorithmically correct. These points
340 C are not used.
341 DO j=jMin,jMax
342 DO i=iMin,iMax
343 #define _recip_VolT1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
344 #define _recip_VolT2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
345 gT(i,j,k,bi,bj)=
346 & -_recip_VolT1(i,j,k,bi,bj)
347 & _recip_VolT2(i,j,k,bi,bj)
348 & *(
349 & +( fZon(i+1,j)-fZon(i,j) )
350 & +( fMer(i,j+1)-fMer(i,j) )
351 & +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )*rkFac
352 & )
353 ENDDO
354 ENDDO
355
356 #ifdef INCLUDE_T_FORCING_CODE
357 C-- External thermal forcing term(s)
358 CALL EXTERNAL_FORCING_T(
359 I iMin,iMax,jMin,jMax,bi,bj,k,
360 I maskC,
361 I myCurrentTime,myThid)
362 #endif /* INCLUDE_T_FORCING_CODE */
363
364 RETURN
365 END

  ViewVC Help
Powered by ViewVC 1.1.22