/[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.22 - (show annotations) (download)
Wed May 26 20:26:42 1999 UTC (25 years, 1 month ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint23, checkpoint24, checkpoint25
Changes since 1.21: +51 -15 lines
Added bi-harmonic diffusion to calc_gs (Salt) and calc_gt (Temperature).

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

  ViewVC Help
Powered by ViewVC 1.1.22