/[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.30 - (show annotations) (download)
Tue Feb 6 03:08:59 2001 UTC (23 years, 4 months ago) by cnh
Branch: MAIN
Changes since 1.29: +12 -2 lines
Changes to avoid references to unitialised variable.
Required to get exp0 through DEC compilers.

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

  ViewVC Help
Powered by ViewVC 1.1.22