/[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.32 - (show annotations) (download)
Tue May 29 14:01:36 2001 UTC (23 years, 1 month ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre3, checkpoint40pre1, checkpoint40pre7, checkpoint40pre6, checkpoint40pre2, checkpoint40pre4, checkpoint40pre5
Changes since 1.31: +70 -78 lines
Merge from branch pre38:
 o essential mods for cubed sphere
 o debugged atmosphere, dynamcis + physics (aim)
 o new packages (mom_vecinv, mom_fluxform, ...)

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gt.F,v 1.31.2.4 2001/04/09 18:52:44 adcroft Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 #define COSINEMETH_III
7 #undef ISOTROPIC_COS_SCALING
8
9 CStartOfInterFace
10 SUBROUTINE CALC_GT(
11 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
12 I xA,yA,uTrans,vTrans,rTrans,maskUp,
13 I KappaRT,
14 U fVerT,
15 I myCurrentTime, myThid )
16 C /==========================================================\
17 C | SUBROUTINE CALC_GT |
18 C | o Calculate the temperature tendency terms. |
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 | logisitics 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 IMPLICIT NONE
42
43 C == GLobal variables ==
44 #include "SIZE.h"
45 #include "DYNVARS.h"
46 #include "EEPARAMS.h"
47 #include "PARAMS.h"
48 #include "GRID.h"
49 #include "FFIELDS.h"
50 c #include "GM_ARRAYS.h"
51
52
53 C == Routine arguments ==
54 C fVerT - Flux of temperature (T) in the vertical
55 C direction at the upper(U) and lower(D) faces of a cell.
56 C maskUp - Land mask used to denote base of the domain.
57 C xA - Tracer cell face area normal to X
58 C yA - Tracer cell face area normal to X
59 C uTrans - Zonal volume transport through cell face
60 C vTrans - Meridional volume transport through cell face
61 C rTrans - Vertical volume transport through cell face
62 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
63 C results will be set.
64 C myThid - Instance number for this innvocation of CALC_GT
65 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
66 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
67 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
68 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
69 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72 _RL KappaRT(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
73 INTEGER k,kUp,kDown,kM1
74 INTEGER bi,bj,iMin,iMax,jMin,jMax
75 INTEGER myThid
76 _RL myCurrentTime
77 CEndOfInterface
78
79 C == Local variables ==
80 C I, J, K - Loop counters
81 INTEGER i,j
82 LOGICAL TOP_LAYER
83 _RL afFacT, dfFacT
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 fVerT(1,1,kDown) = fVerT(1,1,kDown)
94 #endif
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
103 afFacT = 1. _d 0
104 dfFacT = 1. _d 0
105 TOP_LAYER = K .EQ. 1
106
107 C--- Calculate advective and diffusive fluxes between cells.
108
109 #ifdef INCLUDE_T_DIFFUSION_CODE
110 C o Zonal tracer gradient
111 DO j=1-Oly,sNy+Oly
112 DO i=1-Olx+1,sNx+Olx
113 fZon(i,j) = _recip_dxC(i,j,bi,bj)*xA(i,j)
114 & *(theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj))
115 #ifdef COSINEMETH_III
116 & *sqCosFacU(j,bi,bj)
117 #endif
118 ENDDO
119 ENDDO
120 C o Meridional tracer gradient
121 DO j=1-Oly+1,sNy+Oly
122 DO i=1-Olx,sNx+Olx
123 fMer(i,j) = _recip_dyC(i,j,bi,bj)*yA(i,j)
124 & *(theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
125 #ifdef ISOTROPIC_COS_SCALING
126 #ifdef COSINEMETH_III
127 & *sqCosFacV(j,bi,bj)
128 #endif
129 #endif
130 ENDDO
131 ENDDO
132 C-- del^2 of T, needed for bi-harmonic (del^4) term
133 IF (diffK4T .NE. 0.) THEN
134 DO j=1-Oly+1,sNy+Oly-1
135 DO i=1-Olx+1,sNx+Olx-1
136 df4(i,j)= _recip_hFacC(i,j,k,bi,bj)
137 & *recip_drF(k)/_rA(i,j,bi,bj)
138 & *(
139 & +( fZon(i+1,j)-fZon(i,j) )
140 & +( fMer(i,j+1)-fMer(i,j) )
141 & )
142 ENDDO
143 ENDDO
144 ENDIF
145 #endif
146
147 C-- Zonal flux (fZon is at west face of "theta" cell)
148 #ifdef INCLUDE_T_ADVECTION_CODE
149 C o Advective component of zonal flux
150 DO j=jMin,jMax
151 DO i=iMin,iMax
152 af(i,j) =
153 & uTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i-1,j,k,bi,bj))*0.5 _d 0
154 ENDDO
155 ENDDO
156 #endif /* INCLUDE_T_ADVECTION_CODE */
157 #ifdef INCLUDE_T_DIFFUSION_CODE
158 C o Diffusive component of zonal flux
159 DO j=jMin,jMax
160 DO i=iMin,iMax
161 df(i,j) = -diffKhT*xA(i,j)*_recip_dxC(i,j,bi,bj)*
162 & (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj))
163 & *CosFacU(j,bi,bj)
164 ENDDO
165 ENDDO
166 #ifdef ALLOW_GMREDI
167 IF (useGMRedi) CALL GMREDI_XTRANSPORT(
168 I iMin,iMax,jMin,jMax,bi,bj,K,
169 I xA,theta,
170 U df,
171 I myThid)
172 #endif
173 C o Add the bi-harmonic contribution
174 IF (diffK4T .NE. 0.) THEN
175 DO j=jMin,jMax
176 DO i=iMin,iMax
177 df(i,j) = df(i,j) + xA(i,j)*
178 & diffK4T*(df4(i,j)-df4(i-1,j))*_recip_dxC(i,j,bi,bj)
179 #ifdef COSINEMETH_III
180 & *sqCosFacU(j,bi,bj)
181 #else
182 & *CosFacU(j,bi,bj)
183 #endif
184 ENDDO
185 ENDDO
186 ENDIF
187 #endif /* INCLUDE_T_DIFFUSION_CODE */
188 C o Net zonal flux
189 DO j=jMin,jMax
190 DO i=iMin,iMax
191 fZon(i,j) = 0.
192 & _ADT( + afFacT*af(i,j) )
193 & _LPT( + dfFacT*df(i,j) )
194 ENDDO
195 ENDDO
196
197 C-- Meridional flux (fMer is at south face of "theta" cell)
198 #ifdef INCLUDE_T_ADVECTION_CODE
199 C o Advective component of meridional flux
200 DO j=jMin,jMax
201 DO i=iMin,iMax
202 af(i,j) =
203 & vTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j-1,k,bi,bj))*0.5 _d 0
204 ENDDO
205 ENDDO
206 #endif /* INCLUDE_T_ADVECTION_CODE */
207 #ifdef INCLUDE_T_DIFFUSION_CODE
208 C o Diffusive component of meridional flux
209 DO j=jMin,jMax
210 DO i=iMin,iMax
211 df(i,j) = -diffKhT*yA(i,j)*_recip_dyC(i,j,bi,bj)*
212 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
213 #ifdef ISOTROPIC_COS_SCALING
214 & *CosFacV(j,bi,bj)
215 #endif
216 ENDDO
217 ENDDO
218 #ifdef ALLOW_GMREDI
219 IF (useGMRedi) CALL GMREDI_YTRANSPORT(
220 I iMin,iMax,jMin,jMax,bi,bj,K,
221 I yA,theta,
222 U df,
223 I myThid)
224 #endif
225 C o Add the bi-harmonic contribution
226 IF (diffK4T .NE. 0.) THEN
227 DO j=jMin,jMax
228 DO i=iMin,iMax
229 df(i,j) = df(i,j) + yA(i,j)*
230 & diffK4T*(df4(i,j)-df4(i,j-1))*_recip_dyC(i,j,bi,bj)
231 #ifdef ISOTROPIC_COS_SCALING
232 #ifdef COSINEMETH_III
233 & *sqCosFacV(j,bi,bj)
234 #else
235 & *CosFacV(j,bi,bj)
236 #endif
237 #endif
238 ENDDO
239 ENDDO
240 ENDIF
241 #endif /* INCLUDE_T_DIFFUSION_CODE */
242 C o Net meridional flux
243 DO j=jMin,jMax
244 DO i=iMin,iMax
245 fMer(i,j) = 0.
246 & _ADT( + afFacT*af(i,j) )
247 & _LPT( + dfFacT*df(i,j) )
248 ENDDO
249 ENDDO
250
251 C-- Vertical flux ( fVerT(,,kUp) is at upper face of "Tracer" cell )
252 #ifdef INCLUDE_T_ADVECTION_CODE
253 C o Advective component of vertical flux : assume W_bottom=0 (mask)
254 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
255 C (this plays the role of the free-surface correction)
256 IF ( rigidLid .AND. TOP_LAYER) THEN
257 DO j=jMin,jMax
258 DO i=iMin,iMax
259 af(i,j) = 0.
260 ENDDO
261 ENDDO
262 ELSEIF ( rigidLid ) THEN
263 DO j=jMin,jMax
264 DO i=iMin,iMax
265 af(i,j) = rTrans(i,j)*
266 & (theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0
267 ENDDO
268 ENDDO
269 ELSE
270 C- include "free-surface correction" :
271 DO j=jMin,jMax
272 DO i=iMin,iMax
273 af(i,j) = rTrans(i,j)*(
274 & maskC(i,j,kM1,bi,bj)*
275 & (theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0
276 & +(maskC(i,j,k,bi,bj)-maskC(i,j,kM1,bi,bj))*
277 & theta(i,j,k,bi,bj) )
278 ENDDO
279 ENDDO
280 ENDIF
281 #endif /* INCLUDE_T_ADVECTION_CODE */
282 #ifdef INCLUDE_T_DIFFUSION_CODE
283 C o Diffusive component of vertical flux
284 C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
285 C boundary condition.
286 IF (implicitDiffusion) THEN
287 DO j=jMin,jMax
288 DO i=iMin,iMax
289 df(i,j) = 0.
290 ENDDO
291 ENDDO
292 ELSE
293 DO j=jMin,jMax
294 DO i=iMin,iMax
295 df(i,j) = - _rA(i,j,bi,bj)*(
296 & KappaRT(i,j,k)*recip_drC(k)
297 & *(theta(i,j,kM1,bi,bj)-theta(i,j,k,bi,bj))*rkFac
298 & )
299 ENDDO
300 ENDDO
301 ENDIF
302 #endif /* INCLUDE_T_DIFFUSION_CODE */
303
304 #ifdef ALLOW_GMREDI
305 IF (useGMRedi) CALL GMREDI_RTRANSPORT(
306 I iMin,iMax,jMin,jMax,bi,bj,K,
307 I maskUp,theta,
308 U df,
309 I myThid)
310 #endif
311
312 #ifdef ALLOW_KPP
313 C-- Add non local KPP transport term (ghat) to diffusive T flux.
314 IF (useKPP) CALL KPP_TRANSPORT_T(
315 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
316 I KappaRT,
317 U df )
318 #endif
319
320 C o Net vertical flux
321 DO j=jMin,jMax
322 DO i=iMin,iMax
323 c fVerT(i,j,kUp) = afFacT*af(i,j) + dfFacT*df(i,j)*maskUp(i,j)
324 fVerT(i,j,kUp) = 0.
325 & _ADT( +afFacT*af(i,j) )
326 & _LPT( +dfFacT*df(i,j)*maskUp(i,j) )
327 ENDDO
328 ENDDO
329
330 C-- Tendency is minus divergence of the fluxes.
331 C Note. Tendency terms will only be correct for range
332 C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
333 C will contain valid floating point numbers but
334 C they are not algorithmically correct. These points
335 C are not used.
336 DO j=jMin,jMax
337 DO i=iMin,iMax
338 gT(i,j,k,bi,bj)=
339 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
340 & *recip_rA(i,j,bi,bj)
341 & *(
342 & +( fZon(i+1,j)-fZon(i,j) )
343 & +( fMer(i,j+1)-fMer(i,j) )
344 & +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )*rkFac
345 & )
346 ENDDO
347 ENDDO
348
349 #ifdef INCLUDE_T_FORCING_CODE
350 C-- External thermal forcing term(s)
351 CALL EXTERNAL_FORCING_T(
352 I iMin,iMax,jMin,jMax,bi,bj,k,
353 I myCurrentTime,myThid)
354 #endif /* INCLUDE_T_FORCING_CODE */
355
356 RETURN
357 END

  ViewVC Help
Powered by ViewVC 1.1.22