/[MITgcm]/MITgcm/model/src/calc_gt.F
ViewVC logotype

Annotation of /MITgcm/model/src/calc_gt.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.19 - (hide annotations) (download)
Fri Nov 6 22:44:44 1998 UTC (25 years, 7 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint19, checkpoint18, checkpoint20, checkpoint21
Changes since 1.18: +55 -31 lines
Changes to allow for atmospheric integration builds of the code

1 cnh 1.19 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gt.F,v 1.18 1998/11/03 15:28:04 cnh Exp $
2 cnh 1.1
3 cnh 1.19 #include "CPP_OPTIONS.h"
4 cnh 1.1
5     CStartOfInterFace
6     SUBROUTINE CALC_GT(
7     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
8 cnh 1.14 I xA,yA,uTrans,vTrans,rTrans,maskup,maskC,
9     I K13,K23,KappaRT,KapGM,
10 cnh 1.1 U af,df,fZon,fMer,fVerT,
11 cnh 1.18 I myCurrentTime, myThid )
12 cnh 1.1 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 cnh 1.11 #include "FFIELDS.h"
46 cnh 1.1
47     C == Routine arguments ==
48     C fZon - Work array for flux of temperature in the east-west
49     C direction at the west face of a cell.
50     C fMer - Work array for flux of temperature in the north-south
51     C direction at the south face of a cell.
52     C fVerT - Flux of temperature (T) in the vertical
53     C direction at the upper(U) and lower(D) faces of a cell.
54     C maskUp - Land mask used to denote base of the domain.
55 adcroft 1.13 C maskC - Land mask for theta cells (used in TOP_LAYER only)
56 cnh 1.1 C xA - Tracer cell face area normal to X
57     C yA - Tracer cell face area normal to X
58     C uTrans - Zonal volume transport through cell face
59     C vTrans - Meridional volume transport through cell face
60 cnh 1.14 C rTrans - Vertical volume transport through cell face
61 cnh 1.1 C af - Advective flux component work array
62     C df - Diffusive flux component work array
63     C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
64     C results will be set.
65     C myThid - Instance number for this innvocation of CALC_GT
66     _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
67     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
68     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
69     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71     _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72     _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73 cnh 1.14 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 cnh 1.1 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 adcroft 1.13 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 cnh 1.16 _RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
77     _RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
78     _RL KappaRT(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
79 adcroft 1.3 _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 cnh 1.1 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81     _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 adcroft 1.3 INTEGER k,kUp,kDown,kM1
83 cnh 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
84     INTEGER myThid
85 cnh 1.18 _RL myCurrentTime
86 cnh 1.1 CEndOfInterface
87    
88     C == Local variables ==
89     C I, J, K - Loop counters
90 adcroft 1.3 INTEGER i,j
91 cnh 1.10 LOGICAL TOP_LAYER
92 adcroft 1.3 _RL afFacT, dfFacT
93     _RL dTdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94     _RL dTdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 cnh 1.1
96     afFacT = 1. _d 0
97     dfFacT = 1. _d 0
98 cnh 1.10 TOP_LAYER = K .EQ. 1
99 cnh 1.1
100     C--- Calculate advective and diffusive fluxes between cells.
101    
102     C-- Zonal flux (fZon is at west face of "theta" cell)
103 cnh 1.19 #ifdef INCLUDE_T_ADVECTION_CODE
104     C o Advective component of zonal flux
105 cnh 1.1 DO j=jMin,jMax
106     DO i=iMin,iMax
107     af(i,j) =
108     & uTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i-1,j,k,bi,bj))*0.5 _d 0
109     ENDDO
110     ENDDO
111 cnh 1.19 #endif /* INCLUDE_T_ADVECTION_CODE */
112     #ifdef INCLUDE_T_DIFFUSION_CODE
113     C o Zonal tracer gradient
114 adcroft 1.3 DO j=jMin,jMax
115     DO i=iMin,iMax
116 cnh 1.14 dTdx(i,j) = _recip_dxC(i,j,bi,bj)*
117 adcroft 1.3 & (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj))
118     ENDDO
119     ENDDO
120 cnh 1.19 C o Diffusive component of zonal flux
121 cnh 1.1 DO j=jMin,jMax
122     DO i=iMin,iMax
123 adcroft 1.3 df(i,j) = -(diffKhT+0.5*(KapGM(i,j)+KapGM(i-1,j)))*
124     & xA(i,j)*dTdx(i,j)
125 cnh 1.1 ENDDO
126     ENDDO
127 cnh 1.19 #endif /* INCLUDE_T_DIFFUSION_CODE */
128     C o Net zonal flux
129 cnh 1.1 DO j=jMin,jMax
130     DO i=iMin,iMax
131 cnh 1.19 fZon(i,j) = 0.
132     _ADT(& + afFacT*af(i,j) )
133     _LPT(& + dfFacT*df(i,j) )
134 cnh 1.1 ENDDO
135     ENDDO
136    
137     C-- Meridional flux (fMer is at south face of "theta" cell)
138 cnh 1.19 #ifdef INCLUDE_T_ADVECTION_CODE
139     C o Advective component of meridional flux
140 cnh 1.1 DO j=jMin,jMax
141     DO i=iMin,iMax
142     af(i,j) =
143     & vTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j-1,k,bi,bj))*0.5 _d 0
144     ENDDO
145     ENDDO
146 cnh 1.19 #endif /* INCLUDE_T_ADVECTION_CODE */
147     #ifdef INCLUDE_T_DIFFUSION_CODE
148     C o Meridional tracer gradient
149 adcroft 1.3 DO j=jMin,jMax
150     DO i=iMin,iMax
151 cnh 1.14 dTdy(i,j) = _recip_dyC(i,j,bi,bj)*
152 adcroft 1.3 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
153     ENDDO
154     ENDDO
155 cnh 1.19 C o Diffusive component of meridional flux
156 cnh 1.1 DO j=jMin,jMax
157     DO i=iMin,iMax
158 adcroft 1.3 df(i,j) = -(diffKhT+0.5*(KapGM(i,j)+KapGM(i,j-1)))*
159     & yA(i,j)*dTdy(i,j)
160 cnh 1.1 ENDDO
161     ENDDO
162 cnh 1.19 #endif /* INCLUDE_T_DIFFUSION_CODE */
163     C o Net meridional flux
164 cnh 1.1 DO j=jMin,jMax
165     DO i=iMin,iMax
166 cnh 1.19 fMer(i,j) = 0.
167     _ADT(& + afFacT*af(i,j) )
168     _LPT(& + dfFacT*df(i,j) )
169 cnh 1.1 ENDDO
170     ENDDO
171    
172 cnh 1.19 #ifdef INCLUDE_T_DIFFUSION_CODE
173     C-- Terms that diffusion tensor projects onto z
174 adcroft 1.3 DO j=jMin,jMax
175     DO i=iMin,iMax
176     dTdx(i,j) = 0.5*(
177 cnh 1.17 & +0.5*(_maskW(i+1,j,k,bi,bj)
178     & *_recip_dxC(i+1,j,bi,bj)*
179 adcroft 1.3 & (theta(i+1,j,k,bi,bj)-theta(i,j,k,bi,bj))
180 cnh 1.17 & +_maskW(i,j,k,bi,bj)
181     & *_recip_dxC(i,j,bi,bj)*
182 adcroft 1.3 & (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj)))
183 cnh 1.17 & +0.5*(_maskW(i+1,j,km1,bi,bj)
184     & *_recip_dxC(i+1,j,bi,bj)*
185 adcroft 1.3 & (theta(i+1,j,km1,bi,bj)-theta(i,j,km1,bi,bj))
186 cnh 1.17 & +_maskW(i,j,km1,bi,bj)
187     & *_recip_dxC(i,j,bi,bj)*
188 adcroft 1.3 & (theta(i,j,km1,bi,bj)-theta(i-1,j,km1,bi,bj)))
189     & )
190     ENDDO
191     ENDDO
192     DO j=jMin,jMax
193     DO i=iMin,iMax
194     dTdy(i,j) = 0.5*(
195 cnh 1.17 & +0.5*(_maskS(i,j,k,bi,bj)
196     & *_recip_dyC(i,j,bi,bj)*
197 adcroft 1.3 & (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
198 cnh 1.17 & +_maskS(i,j+1,k,bi,bj)
199     & *_recip_dyC(i,j+1,bi,bj)*
200 adcroft 1.3 & (theta(i,j+1,k,bi,bj)-theta(i,j,k,bi,bj)))
201 cnh 1.17 & +0.5*(_maskS(i,j,km1,bi,bj)
202     & *_recip_dyC(i,j,bi,bj)*
203 adcroft 1.3 & (theta(i,j,km1,bi,bj)-theta(i,j-1,km1,bi,bj))
204 cnh 1.17 & +_maskS(i,j+1,km1,bi,bj)
205     & *_recip_dyC(i,j+1,bi,bj)*
206 adcroft 1.3 & (theta(i,j+1,km1,bi,bj)-theta(i,j,km1,bi,bj)))
207     & )
208     ENDDO
209     ENDDO
210 cnh 1.19 #endif /* INCLUDE_T_DIFFUSION_CODE */
211 adcroft 1.3
212 cnh 1.19 C-- Vertical flux ( fVerT(,,kUp) is at upper face of "theta" cell )
213     #ifdef INCLUDE_T_ADVECTION_CODE
214     C o Advective component of vertical flux
215 adcroft 1.3 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
216     C (this plays the role of the free-surface correction)
217 cnh 1.1 DO j=jMin,jMax
218     DO i=iMin,iMax
219     af(i,j) =
220 cnh 1.14 & rTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0
221 cnh 1.1 ENDDO
222     ENDDO
223 cnh 1.19 #endif /* INCLUDE_T_ADVECTION_CODE */
224     #ifdef INCLUDE_T_DIFFUSION_CODE
225     C o Diffusive component of vertical flux
226     C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
227 adcroft 1.3 C boundary condition.
228 cnh 1.1 DO j=jMin,jMax
229     DO i=iMin,iMax
230 cnh 1.14 df(i,j) = _rA(i,j,bi,bj)*(
231 adcroft 1.3 & -KapGM(i,j)*K13(i,j,k)*dTdx(i,j)
232     & -KapGM(i,j)*K23(i,j,k)*dTdy(i,j)
233     & )
234 cnh 1.1 ENDDO
235     ENDDO
236 adcroft 1.9 IF (.NOT.implicitDiffusion) THEN
237     DO j=jMin,jMax
238     DO i=iMin,iMax
239 cnh 1.14 df(i,j) = df(i,j) + _rA(i,j,bi,bj)*(
240 cnh 1.16 & -KappaRT(i,j,k)*recip_drC(k)
241 cnh 1.15 & *(theta(i,j,kM1,bi,bj)-theta(i,j,k,bi,bj))*rkFac
242 adcroft 1.9 & )
243     ENDDO
244     ENDDO
245     ENDIF
246 cnh 1.19 #endif /* INCLUDE_T_DIFFUSION_CODE */
247     C o Net vertical flux
248 cnh 1.1 DO j=jMin,jMax
249     DO i=iMin,iMax
250 cnh 1.19 fVerT(i,j,kUp) = 0.
251     _ADT(& +afFacT*af(i,j)*maskUp(i,j) )
252     _LPT(& +dfFacT*df(i,j)*maskUp(i,j) )
253 cnh 1.1 ENDDO
254     ENDDO
255 cnh 1.19 #ifdef INCLUDE_T_ADVECTION_CODE
256 cnh 1.10 IF ( TOP_LAYER ) THEN
257     DO j=jMin,jMax
258     DO i=iMin,iMax
259     fVerT(i,j,kUp) = afFacT*af(i,j)*freeSurfFac
260     ENDDO
261     ENDDO
262     ENDIF
263 cnh 1.19 #endif /* INCLUDE_T_ADVECTION_CODE */
264 cnh 1.1
265     C-- Tendency is minus divergence of the fluxes.
266     C Note. Tendency terms will only be correct for range
267     C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
268     C will contain valid floating point numbers but
269     C they are not algorithmically correct. These points
270     C are not used.
271     DO j=jMin,jMax
272     DO i=iMin,iMax
273 cnh 1.17 #define _recip_VolT1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
274     #define _recip_VolT2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
275 cnh 1.1 gT(i,j,k,bi,bj)=
276 cnh 1.17 & -_recip_VolT1(i,j,k,bi,bj)
277     & _recip_VolT2(i,j,k,bi,bj)
278 cnh 1.1 & *(
279     & +( fZon(i+1,j)-fZon(i,j) )
280     & +( fMer(i,j+1)-fMer(i,j) )
281 cnh 1.14 & +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )*rkFac
282 cnh 1.1 & )
283     ENDDO
284     ENDDO
285    
286 cnh 1.19 #ifdef INCLUDE_T_FORCING_CODE
287 cnh 1.1 C-- External thermal forcing term(s)
288 cnh 1.19 CALL EXTERNAL_FORCING_T(
289     I iMin,iMax,jMin,jMax,bi,bj,k,
290     I maskC,
291     I myCurrentTime,myThid)
292     #endif /* INCLUDE_T_FORCING_CODE */
293    
294     #ifdef INCLUDE_LAT_CIRC_FFT_FILTER_CODE
295     C-- Zonal FFT filter of tendency
296     CALL FILTER_LATCIRCS_FFT_APPLY(
297     U gT,
298     I 1, sNy, k, k, bi, bj, 1, myThid)
299     #endif /* INCLUDE_LAT_CIRC_FFT_FILTER_CODE */
300    
301 cnh 1.1
302     RETURN
303     END

  ViewVC Help
Powered by ViewVC 1.1.22