/[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.20 - (show annotations) (download)
Tue May 18 18:01:12 1999 UTC (25 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint22
Changes since 1.19: +57 -1 lines
Modifications/additions for KPP mixing scheme. Instigated by Dimitri.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gt.F,v 1.19 1998/11/06 22:44:44 cnh 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 _RL hbl (1-OLx:sNx+OLx,1-OLy:sNy+OLy) ! used by KPP mixing scheme
101 _RL frac (1-OLx:sNx+OLx,1-OLy:sNy+OLy) ! used by KPP mixing scheme
102 integer jwtype ! index for Jerlov water type
103 #endif
104
105 afFacT = 1. _d 0
106 dfFacT = 1. _d 0
107 TOP_LAYER = K .EQ. 1
108
109 C--- Calculate advective and diffusive fluxes between cells.
110
111 C-- Zonal flux (fZon is at west face of "theta" cell)
112 #ifdef INCLUDE_T_ADVECTION_CODE
113 C o Advective component of zonal flux
114 DO j=jMin,jMax
115 DO i=iMin,iMax
116 af(i,j) =
117 & uTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i-1,j,k,bi,bj))*0.5 _d 0
118 ENDDO
119 ENDDO
120 #endif /* INCLUDE_T_ADVECTION_CODE */
121 #ifdef INCLUDE_T_DIFFUSION_CODE
122 C o Zonal tracer gradient
123 DO j=jMin,jMax
124 DO i=iMin,iMax
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 Diffusive component of zonal flux
130 DO j=jMin,jMax
131 DO i=iMin,iMax
132 df(i,j) = -(diffKhT+0.5*(KapGM(i,j)+KapGM(i-1,j)))*
133 & xA(i,j)*dTdx(i,j)
134 ENDDO
135 ENDDO
136 #endif /* INCLUDE_T_DIFFUSION_CODE */
137 C o Net zonal flux
138 DO j=jMin,jMax
139 DO i=iMin,iMax
140 fZon(i,j) = 0.
141 _ADT(& + afFacT*af(i,j) )
142 _LPT(& + dfFacT*df(i,j) )
143 ENDDO
144 ENDDO
145
146 C-- Meridional flux (fMer is at south face of "theta" cell)
147 #ifdef INCLUDE_T_ADVECTION_CODE
148 C o Advective component of meridional flux
149 DO j=jMin,jMax
150 DO i=iMin,iMax
151 af(i,j) =
152 & vTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j-1,k,bi,bj))*0.5 _d 0
153 ENDDO
154 ENDDO
155 #endif /* INCLUDE_T_ADVECTION_CODE */
156 #ifdef INCLUDE_T_DIFFUSION_CODE
157 C o Meridional tracer gradient
158 DO j=jMin,jMax
159 DO i=iMin,iMax
160 dTdy(i,j) = _recip_dyC(i,j,bi,bj)*
161 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
162 ENDDO
163 ENDDO
164 C o Diffusive component of meridional flux
165 DO j=jMin,jMax
166 DO i=iMin,iMax
167 df(i,j) = -(diffKhT+0.5*(KapGM(i,j)+KapGM(i,j-1)))*
168 & yA(i,j)*dTdy(i,j)
169 ENDDO
170 ENDDO
171 #endif /* INCLUDE_T_DIFFUSION_CODE */
172 C o Net meridional flux
173 DO j=jMin,jMax
174 DO i=iMin,iMax
175 fMer(i,j) = 0.
176 _ADT(& + afFacT*af(i,j) )
177 _LPT(& + dfFacT*df(i,j) )
178 ENDDO
179 ENDDO
180
181 #ifdef INCLUDE_T_DIFFUSION_CODE
182 C-- Terms that diffusion tensor projects onto z
183 DO j=jMin,jMax
184 DO i=iMin,iMax
185 dTdx(i,j) = 0.5*(
186 & +0.5*(_maskW(i+1,j,k,bi,bj)
187 & *_recip_dxC(i+1,j,bi,bj)*
188 & (theta(i+1,j,k,bi,bj)-theta(i,j,k,bi,bj))
189 & +_maskW(i,j,k,bi,bj)
190 & *_recip_dxC(i,j,bi,bj)*
191 & (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj)))
192 & +0.5*(_maskW(i+1,j,km1,bi,bj)
193 & *_recip_dxC(i+1,j,bi,bj)*
194 & (theta(i+1,j,km1,bi,bj)-theta(i,j,km1,bi,bj))
195 & +_maskW(i,j,km1,bi,bj)
196 & *_recip_dxC(i,j,bi,bj)*
197 & (theta(i,j,km1,bi,bj)-theta(i-1,j,km1,bi,bj)))
198 & )
199 ENDDO
200 ENDDO
201 DO j=jMin,jMax
202 DO i=iMin,iMax
203 dTdy(i,j) = 0.5*(
204 & +0.5*(_maskS(i,j,k,bi,bj)
205 & *_recip_dyC(i,j,bi,bj)*
206 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
207 & +_maskS(i,j+1,k,bi,bj)
208 & *_recip_dyC(i,j+1,bi,bj)*
209 & (theta(i,j+1,k,bi,bj)-theta(i,j,k,bi,bj)))
210 & +0.5*(_maskS(i,j,km1,bi,bj)
211 & *_recip_dyC(i,j,bi,bj)*
212 & (theta(i,j,km1,bi,bj)-theta(i,j-1,km1,bi,bj))
213 & +_maskS(i,j+1,km1,bi,bj)
214 & *_recip_dyC(i,j+1,bi,bj)*
215 & (theta(i,j+1,km1,bi,bj)-theta(i,j,km1,bi,bj)))
216 & )
217 ENDDO
218 ENDDO
219 #endif /* INCLUDE_T_DIFFUSION_CODE */
220
221 C-- Vertical flux ( fVerT(,,kUp) is at upper face of "theta" cell )
222 #ifdef INCLUDE_T_ADVECTION_CODE
223 C o Advective component of vertical flux
224 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
225 C (this plays the role of the free-surface correction)
226 DO j=jMin,jMax
227 DO i=iMin,iMax
228 af(i,j) =
229 & rTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0
230 ENDDO
231 ENDDO
232 #endif /* INCLUDE_T_ADVECTION_CODE */
233 #ifdef INCLUDE_T_DIFFUSION_CODE
234 C o Diffusive component of vertical flux
235 C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
236 C boundary condition.
237 DO j=jMin,jMax
238 DO i=iMin,iMax
239 df(i,j) = _rA(i,j,bi,bj)*(
240 & -KapGM(i,j)*K13(i,j,k)*dTdx(i,j)
241 & -KapGM(i,j)*K23(i,j,k)*dTdy(i,j)
242 & )
243 ENDDO
244 ENDDO
245 IF (.NOT.implicitDiffusion) THEN
246 DO j=jMin,jMax
247 DO i=iMin,iMax
248 df(i,j) = df(i,j) + _rA(i,j,bi,bj)*(
249 & -KappaRT(i,j,k)*recip_drC(k)
250 & *(theta(i,j,kM1,bi,bj)-theta(i,j,k,bi,bj))*rkFac
251 & )
252 ENDDO
253 ENDDO
254 ENDIF
255 #endif /* INCLUDE_T_DIFFUSION_CODE */
256
257 #ifdef ALLOW_KPP
258 IF (usingKPPmixing) THEN
259 C-- Compute fraction of solar short-wave flux penetrating to
260 C the bottom of the mixing layer
261 DO j=jMin,jMax
262 DO i=iMin,iMax
263 hbl(i,j) = KPPhbl(i,j,bi,bj)
264 ENDDO
265 ENDDO
266 j=(sNx+2*OLx)*(sNy+2*OLy)
267 jwtype = 3
268 CALL SWFRAC(
269 I j, -1., hbl, jwtype,
270 O frac )
271
272 C Add non local transport coefficient (ghat term) to right-hand-side
273 C The nonlocal transport term is noNrero only for scalars in unstable
274 C (convective) forcing conditions.
275 C Note: -[Qnet * delZ(1) + Qsw * (1-frac) / KPPhbl] * 4000 * rho
276 C is the total heat flux
277 C penetrating the mixed layer from the surface in (deg C / s)
278 IF ( TOP_LAYER ) THEN
279 DO j=jMin,jMax
280 DO i=iMin,iMax
281 df(i,j) = df(i,j) + _rA(i,j,bi,bj) *
282 & ( Qnet(i,j,bi,bj) * delZ(1) +
283 & Qsw(i,j,bi,bj) * (1.-frac(i,j))
284 & / KPPhbl(i,j,bi,bj) ) *
285 & ( KappaRT(i,j,k) * KPPghat(i,j,k, bi,bj) )
286 ENDDO
287 ENDDO
288 ELSE
289 DO j=jMin,jMax
290 DO i=iMin,iMax
291 df(i,j) = df(i,j) + _rA(i,j,bi,bj) *
292 & ( Qnet(i,j,bi,bj) * delZ(1) +
293 & Qsw(i,j,bi,bj) * (1.-frac(i,j))
294 & / KPPhbl(i,j,bi,bj) ) *
295 & ( KappaRT(i,j,k) * KPPghat(i,j,k, bi,bj)
296 & - KappaRT(i,j,k-1) * KPPghat(i,j,k-1,bi,bj) )
297 ENDDO
298 ENDDO
299 ENDIF
300 ENDIF
301 #endif /* ALLOW_KPP */
302
303 C o Net vertical flux
304 DO j=jMin,jMax
305 DO i=iMin,iMax
306 fVerT(i,j,kUp) = 0.
307 _ADT(& +afFacT*af(i,j)*maskUp(i,j) )
308 _LPT(& +dfFacT*df(i,j)*maskUp(i,j) )
309 ENDDO
310 ENDDO
311 #ifdef INCLUDE_T_ADVECTION_CODE
312 IF ( TOP_LAYER ) THEN
313 DO j=jMin,jMax
314 DO i=iMin,iMax
315 fVerT(i,j,kUp) = afFacT*af(i,j)*freeSurfFac
316 ENDDO
317 ENDDO
318 ENDIF
319 #endif /* INCLUDE_T_ADVECTION_CODE */
320
321 C-- Tendency is minus divergence of the fluxes.
322 C Note. Tendency terms will only be correct for range
323 C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
324 C will contain valid floating point numbers but
325 C they are not algorithmically correct. These points
326 C are not used.
327 DO j=jMin,jMax
328 DO i=iMin,iMax
329 #define _recip_VolT1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
330 #define _recip_VolT2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
331 gT(i,j,k,bi,bj)=
332 & -_recip_VolT1(i,j,k,bi,bj)
333 & _recip_VolT2(i,j,k,bi,bj)
334 & *(
335 & +( fZon(i+1,j)-fZon(i,j) )
336 & +( fMer(i,j+1)-fMer(i,j) )
337 & +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )*rkFac
338 & )
339 ENDDO
340 ENDDO
341
342 #ifdef INCLUDE_T_FORCING_CODE
343 C-- External thermal forcing term(s)
344 CALL EXTERNAL_FORCING_T(
345 I iMin,iMax,jMin,jMax,bi,bj,k,
346 I maskC,
347 I myCurrentTime,myThid)
348 #endif /* INCLUDE_T_FORCING_CODE */
349
350 #ifdef INCLUDE_LAT_CIRC_FFT_FILTER_CODE
351 C-- Zonal FFT filter of tendency
352 CALL FILTER_LATCIRCS_FFT_APPLY(
353 U gT,
354 I 1, sNy, k, k, bi, bj, 1, myThid)
355 #endif /* INCLUDE_LAT_CIRC_FFT_FILTER_CODE */
356
357
358 RETURN
359 END

  ViewVC Help
Powered by ViewVC 1.1.22