/[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.21 - (show annotations) (download)
Mon May 24 14:24:24 1999 UTC (25 years ago) by adcroft
Branch: MAIN
Changes since 1.20: +6 -4 lines
Fixed up types for arguments to SWFRAC().

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gt.F,v 1.20 1999/05/18 18:01:12 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 K13,K23,KappaRT,KapGM,
10 U af,df,fZon,fMer,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 #ifdef ALLOW_KPP
47 #include "KPPMIX.h"
48 #endif
49
50
51 C == Routine arguments ==
52 C fZon - Work array for flux of temperature in the east-west
53 C direction at the west face of a cell.
54 C fMer - Work array for flux of temperature in the north-south
55 C direction at the south face of a cell.
56 C fVerT - Flux of temperature (T) in the vertical
57 C direction at the upper(U) and lower(D) faces of a cell.
58 C maskUp - Land mask used to denote base of the domain.
59 C maskC - Land mask for theta cells (used in TOP_LAYER only)
60 C xA - Tracer cell face area normal to X
61 C yA - Tracer cell face area normal to X
62 C uTrans - Zonal volume transport through cell face
63 C vTrans - Meridional volume transport through cell face
64 C rTrans - Vertical volume transport through cell face
65 C af - Advective flux component work array
66 C df - Diffusive flux component work array
67 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
68 C results will be set.
69 C myThid - Instance number for this innvocation of CALC_GT
70 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71 _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
73 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
81 _RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
82 _RL KappaRT(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
83 _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 INTEGER k,kUp,kDown,kM1
87 INTEGER bi,bj,iMin,iMax,jMin,jMax
88 INTEGER myThid
89 _RL myCurrentTime
90 CEndOfInterface
91
92 C == Local variables ==
93 C I, J, K - Loop counters
94 INTEGER i,j
95 LOGICAL TOP_LAYER
96 _RL afFacT, dfFacT
97 _RL dTdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RL dTdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 #ifdef ALLOW_KPP
100 _RS hbl (1-OLx:sNx+OLx,1-OLy:sNy+OLy) ! used by KPP mixing scheme
101 _RS frac (1-OLx:sNx+OLx,1-OLy:sNy+OLy) ! used by KPP mixing scheme
102 _RS negone ! used as argument to SWFRAC
103 integer jwtype ! index for Jerlov water type
104 #endif
105
106 afFacT = 1. _d 0
107 dfFacT = 1. _d 0
108 TOP_LAYER = K .EQ. 1
109
110 C--- Calculate advective and diffusive fluxes between cells.
111
112 C-- Zonal flux (fZon is at west face of "theta" cell)
113 #ifdef INCLUDE_T_ADVECTION_CODE
114 C o Advective component of zonal flux
115 DO j=jMin,jMax
116 DO i=iMin,iMax
117 af(i,j) =
118 & uTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i-1,j,k,bi,bj))*0.5 _d 0
119 ENDDO
120 ENDDO
121 #endif /* INCLUDE_T_ADVECTION_CODE */
122 #ifdef INCLUDE_T_DIFFUSION_CODE
123 C o Zonal tracer gradient
124 DO j=jMin,jMax
125 DO i=iMin,iMax
126 dTdx(i,j) = _recip_dxC(i,j,bi,bj)*
127 & (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj))
128 ENDDO
129 ENDDO
130 C o Diffusive component of zonal flux
131 DO j=jMin,jMax
132 DO i=iMin,iMax
133 df(i,j) = -(diffKhT+0.5*(KapGM(i,j)+KapGM(i-1,j)))*
134 & xA(i,j)*dTdx(i,j)
135 ENDDO
136 ENDDO
137 #endif /* INCLUDE_T_DIFFUSION_CODE */
138 C o Net zonal flux
139 DO j=jMin,jMax
140 DO i=iMin,iMax
141 fZon(i,j) = 0.
142 _ADT(& + afFacT*af(i,j) )
143 _LPT(& + dfFacT*df(i,j) )
144 ENDDO
145 ENDDO
146
147 C-- Meridional flux (fMer is at south face of "theta" cell)
148 #ifdef INCLUDE_T_ADVECTION_CODE
149 C o Advective component of meridional flux
150 DO j=jMin,jMax
151 DO i=iMin,iMax
152 af(i,j) =
153 & vTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j-1,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 Meridional tracer gradient
159 DO j=jMin,jMax
160 DO i=iMin,iMax
161 dTdy(i,j) = _recip_dyC(i,j,bi,bj)*
162 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
163 ENDDO
164 ENDDO
165 C o Diffusive component of meridional flux
166 DO j=jMin,jMax
167 DO i=iMin,iMax
168 df(i,j) = -(diffKhT+0.5*(KapGM(i,j)+KapGM(i,j-1)))*
169 & yA(i,j)*dTdy(i,j)
170 ENDDO
171 ENDDO
172 #endif /* INCLUDE_T_DIFFUSION_CODE */
173 C o Net meridional flux
174 DO j=jMin,jMax
175 DO i=iMin,iMax
176 fMer(i,j) = 0.
177 _ADT(& + afFacT*af(i,j) )
178 _LPT(& + dfFacT*df(i,j) )
179 ENDDO
180 ENDDO
181
182 #ifdef INCLUDE_T_DIFFUSION_CODE
183 C-- Terms that diffusion tensor projects onto z
184 DO j=jMin,jMax
185 DO i=iMin,iMax
186 dTdx(i,j) = 0.5*(
187 & +0.5*(_maskW(i+1,j,k,bi,bj)
188 & *_recip_dxC(i+1,j,bi,bj)*
189 & (theta(i+1,j,k,bi,bj)-theta(i,j,k,bi,bj))
190 & +_maskW(i,j,k,bi,bj)
191 & *_recip_dxC(i,j,bi,bj)*
192 & (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj)))
193 & +0.5*(_maskW(i+1,j,km1,bi,bj)
194 & *_recip_dxC(i+1,j,bi,bj)*
195 & (theta(i+1,j,km1,bi,bj)-theta(i,j,km1,bi,bj))
196 & +_maskW(i,j,km1,bi,bj)
197 & *_recip_dxC(i,j,bi,bj)*
198 & (theta(i,j,km1,bi,bj)-theta(i-1,j,km1,bi,bj)))
199 & )
200 ENDDO
201 ENDDO
202 DO j=jMin,jMax
203 DO i=iMin,iMax
204 dTdy(i,j) = 0.5*(
205 & +0.5*(_maskS(i,j,k,bi,bj)
206 & *_recip_dyC(i,j,bi,bj)*
207 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
208 & +_maskS(i,j+1,k,bi,bj)
209 & *_recip_dyC(i,j+1,bi,bj)*
210 & (theta(i,j+1,k,bi,bj)-theta(i,j,k,bi,bj)))
211 & +0.5*(_maskS(i,j,km1,bi,bj)
212 & *_recip_dyC(i,j,bi,bj)*
213 & (theta(i,j,km1,bi,bj)-theta(i,j-1,km1,bi,bj))
214 & +_maskS(i,j+1,km1,bi,bj)
215 & *_recip_dyC(i,j+1,bi,bj)*
216 & (theta(i,j+1,km1,bi,bj)-theta(i,j,km1,bi,bj)))
217 & )
218 ENDDO
219 ENDDO
220 #endif /* INCLUDE_T_DIFFUSION_CODE */
221
222 C-- Vertical flux ( fVerT(,,kUp) is at upper face of "theta" cell )
223 #ifdef INCLUDE_T_ADVECTION_CODE
224 C o Advective component of vertical flux
225 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
226 C (this plays the role of the free-surface correction)
227 DO j=jMin,jMax
228 DO i=iMin,iMax
229 af(i,j) =
230 & rTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0
231 ENDDO
232 ENDDO
233 #endif /* INCLUDE_T_ADVECTION_CODE */
234 #ifdef INCLUDE_T_DIFFUSION_CODE
235 C o Diffusive component of vertical flux
236 C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
237 C boundary condition.
238 DO j=jMin,jMax
239 DO i=iMin,iMax
240 df(i,j) = _rA(i,j,bi,bj)*(
241 & -KapGM(i,j)*K13(i,j,k)*dTdx(i,j)
242 & -KapGM(i,j)*K23(i,j,k)*dTdy(i,j)
243 & )
244 ENDDO
245 ENDDO
246 IF (.NOT.implicitDiffusion) THEN
247 DO j=jMin,jMax
248 DO i=iMin,iMax
249 df(i,j) = df(i,j) + _rA(i,j,bi,bj)*(
250 & -KappaRT(i,j,k)*recip_drC(k)
251 & *(theta(i,j,kM1,bi,bj)-theta(i,j,k,bi,bj))*rkFac
252 & )
253 ENDDO
254 ENDDO
255 ENDIF
256 #endif /* INCLUDE_T_DIFFUSION_CODE */
257
258 #ifdef ALLOW_KPP
259 IF (usingKPPmixing) THEN
260 C-- Compute fraction of solar short-wave flux penetrating to
261 C the bottom of the mixing layer
262 DO j=jMin,jMax
263 DO i=iMin,iMax
264 hbl(i,j) = KPPhbl(i,j,bi,bj)
265 ENDDO
266 ENDDO
267 j=(sNx+2*OLx)*(sNy+2*OLy)
268 jwtype = 3
269 negone = -1.
270 CALL SWFRAC(
271 I j, negone, hbl, jwtype,
272 O frac )
273
274 C Add non local transport coefficient (ghat term) to right-hand-side
275 C The nonlocal transport term is noNrero only for scalars in unstable
276 C (convective) forcing conditions.
277 C Note: -[Qnet * delZ(1) + Qsw * (1-frac) / KPPhbl] * 4000 * rho
278 C is the total heat flux
279 C penetrating the mixed layer from the surface in (deg C / s)
280 IF ( TOP_LAYER ) THEN
281 DO j=jMin,jMax
282 DO i=iMin,iMax
283 df(i,j) = df(i,j) + _rA(i,j,bi,bj) *
284 & ( Qnet(i,j,bi,bj) * delZ(1) +
285 & Qsw(i,j,bi,bj) * (1.-frac(i,j))
286 & / KPPhbl(i,j,bi,bj) ) *
287 & ( KappaRT(i,j,k) * KPPghat(i,j,k, bi,bj) )
288 ENDDO
289 ENDDO
290 ELSE
291 DO j=jMin,jMax
292 DO i=iMin,iMax
293 df(i,j) = df(i,j) + _rA(i,j,bi,bj) *
294 & ( Qnet(i,j,bi,bj) * delZ(1) +
295 & Qsw(i,j,bi,bj) * (1.-frac(i,j))
296 & / KPPhbl(i,j,bi,bj) ) *
297 & ( KappaRT(i,j,k) * KPPghat(i,j,k, bi,bj)
298 & - KappaRT(i,j,k-1) * KPPghat(i,j,k-1,bi,bj) )
299 ENDDO
300 ENDDO
301 ENDIF
302 ENDIF
303 #endif /* ALLOW_KPP */
304
305 C o Net vertical flux
306 DO j=jMin,jMax
307 DO i=iMin,iMax
308 fVerT(i,j,kUp) = 0.
309 _ADT(& +afFacT*af(i,j)*maskUp(i,j) )
310 _LPT(& +dfFacT*df(i,j)*maskUp(i,j) )
311 ENDDO
312 ENDDO
313 #ifdef INCLUDE_T_ADVECTION_CODE
314 IF ( TOP_LAYER ) THEN
315 DO j=jMin,jMax
316 DO i=iMin,iMax
317 fVerT(i,j,kUp) = afFacT*af(i,j)*freeSurfFac
318 ENDDO
319 ENDDO
320 ENDIF
321 #endif /* INCLUDE_T_ADVECTION_CODE */
322
323 C-- Tendency is minus divergence of the fluxes.
324 C Note. Tendency terms will only be correct for range
325 C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
326 C will contain valid floating point numbers but
327 C they are not algorithmically correct. These points
328 C are not used.
329 DO j=jMin,jMax
330 DO i=iMin,iMax
331 #define _recip_VolT1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
332 #define _recip_VolT2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
333 gT(i,j,k,bi,bj)=
334 & -_recip_VolT1(i,j,k,bi,bj)
335 & _recip_VolT2(i,j,k,bi,bj)
336 & *(
337 & +( fZon(i+1,j)-fZon(i,j) )
338 & +( fMer(i,j+1)-fMer(i,j) )
339 & +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )*rkFac
340 & )
341 ENDDO
342 ENDDO
343
344 #ifdef INCLUDE_T_FORCING_CODE
345 C-- External thermal forcing term(s)
346 CALL EXTERNAL_FORCING_T(
347 I iMin,iMax,jMin,jMax,bi,bj,k,
348 I maskC,
349 I myCurrentTime,myThid)
350 #endif /* INCLUDE_T_FORCING_CODE */
351
352 #ifdef INCLUDE_LAT_CIRC_FFT_FILTER_CODE
353 C-- Zonal FFT filter of tendency
354 CALL FILTER_LATCIRCS_FFT_APPLY(
355 U gT,
356 I 1, sNy, k, k, bi, bj, 1, myThid)
357 #endif /* INCLUDE_LAT_CIRC_FFT_FILTER_CODE */
358
359
360 RETURN
361 END

  ViewVC Help
Powered by ViewVC 1.1.22