/[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.2 - (hide annotations) (download)
Fri Apr 24 02:05:40 1998 UTC (26 years, 1 month ago) by cnh
Branch: MAIN
CVS Tags: checkpoint1, kloop1, kloop2
Changes since 1.1: +1 -1 lines
Further $Id to $Header conversions

1 cnh 1.2 C $Header: calc_gt.F,v 1.1.1.1 1998/04/22 19:15:30 cnh Exp $
2 cnh 1.1
3     #include "CPP_EEOPTIONS.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,wTrans,maskup,
9     U af,df,fZon,fMer,fVerT,
10     I myThid )
11     C /==========================================================\
12     C | SUBROUTINE CALC_GT |
13     C | o Calculate the temperature tendency terms. |
14     C |==========================================================|
15     C | A procedure called EXTERNAL_FORCING_T is called from |
16     C | here. These procedures can be used to add per problem |
17     C | heat flux source terms. |
18     C | Note: Although it is slightly counter-intuitive the |
19     C | EXTERNAL_FORCING routine is not the place to put |
20     C | file I/O. Instead files that are required to |
21     C | calculate the external source terms are generally |
22     C | read during the model main loop. This makes the |
23     C | logisitics of multi-processing simpler and also |
24     C | makes the adjoint generation simpler. It also |
25     C | allows for I/O to overlap computation where that |
26     C | is supported by hardware. |
27     C | Aside from the problem specific term the code here |
28     C | forms the tendency terms due to advection and mixing |
29     C | The baseline implementation here uses a centered |
30     C | difference form for the advection term and a tensorial |
31     C | divergence of a flux form for the diffusive term. The |
32     C | diffusive term is formulated so that isopycnal mixing and|
33     C | GM-style subgrid-scale terms can be incorporated b simply|
34     C | setting the diffusion tensor terms appropriately. |
35     C \==========================================================/
36     IMPLICIT NONE
37    
38     C == GLobal variables ==
39     #include "SIZE.h"
40     #include "DYNVARS.h"
41     #include "EEPARAMS.h"
42     #include "PARAMS.h"
43     #include "GRID.h"
44    
45     C == Routine arguments ==
46     C fZon - Work array for flux of temperature in the east-west
47     C direction at the west face of a cell.
48     C fMer - Work array for flux of temperature in the north-south
49     C direction at the south face of a cell.
50     C fVerT - Flux of temperature (T) in the vertical
51     C direction at the upper(U) and lower(D) faces of a cell.
52     C maskUp - Land mask used to denote base of the domain.
53     C xA - Tracer cell face area normal to X
54     C yA - Tracer cell face area normal to X
55     C uTrans - Zonal volume transport through cell face
56     C vTrans - Meridional volume transport through cell face
57     C wTrans - Vertical volume transport through cell face
58     C af - Advective flux component work array
59     C df - Diffusive flux component work array
60     C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
61     C results will be set.
62     C myThid - Instance number for this innvocation of CALC_GT
63     _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
64     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
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 wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71     _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72     _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73     _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74     INTEGER kUp,kDown,kM1
75     INTEGER bi,bj,iMin,iMax,jMin,jMax
76     INTEGER myThid
77     CEndOfInterface
78    
79     C == Local variables ==
80     C I, J, K - Loop counters
81     INTEGER i,j,k
82     REAL afFacT, dfFacT
83     REAL dutdxFac
84    
85     afFacT = 1. _d 0
86     dfFacT = 1. _d 0
87     dutdxFac = afFacT
88    
89     C---
90     C--- Calculate advective and diffusive fluxes between cells.
91     C---
92    
93     C-- Zonal flux (fZon is at west face of "theta" cell)
94     C Advective component of zonal flux
95     DO j=jMin,jMax
96     DO i=iMin,iMax
97     af(i,j) =
98     & uTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i-1,j,k,bi,bj))*0.5 _d 0
99     ENDDO
100     ENDDO
101     C Diffusive component of zonal flux
102     DO j=jMin,jMax
103     DO i=iMin,iMax
104     df(i,j) =
105     & -diffKhT*xA(i,j)*rdxC(i,j,bi,bj)
106     & *(theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj))
107     ENDDO
108     ENDDO
109     C Net zonal flux
110     DO j=jMin,jMax
111     DO i=iMin,iMax
112     fZon(i,j) = afFacT*af(i,j) + dfFacT*df(i,j)
113     ENDDO
114     ENDDO
115    
116     C-- Meridional flux (fMer is at south face of "theta" cell)
117     C Advective component of meridional flux
118     DO j=jMin,jMax
119     DO i=iMin,iMax
120     C Advective component of meridional flux
121     af(i,j) =
122     & vTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j-1,k,bi,bj))*0.5 _d 0
123     ENDDO
124     ENDDO
125     C Diffusive component of meridional flux
126     DO j=jMin,jMax
127     DO i=iMin,iMax
128     df(i,j) =
129     & -diffKhT*yA(i,j)*rdyC(i,j,bi,bj)
130     C & -1.D3*rdyC(i,j,bi,bj)*dZF(K)*delX(1)*hFacC(i,j,k,bi,bj)*
131     C & hFacC(i,j-1,k,bi,bj)
132     & *(theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
133     ENDDO
134     ENDDO
135     C Net meridional flux
136     DO j=jMin,jMax
137     DO i=iMin,iMax
138     fMer(i,j) = afFacT*af(i,j) + dfFacT*df(i,j)
139     ENDDO
140     ENDDO
141    
142     C-- Vertical flux (fVerT) above
143     C Note: For K=1 then KM1=1 this gives a dT/dz = 0 upper
144     C boundary condition.
145     C Advective component of vertical flux
146     DO j=jMin,jMax
147     DO i=iMin,iMax
148     af(i,j) =
149     & wTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0
150     ENDDO
151     ENDDO
152     C Diffusive component of vertical flux
153     DO j=jMin,jMax
154     DO i=iMin,iMax
155     df(i,j) =
156     & -diffKzT*zA(i,j,bi,bj)*rdzC(k)
157     & *(theta(i,j,kM1,bi,bj)-theta(i,j,k,bi,bj))
158     ENDDO
159     ENDDO
160     C Net vertical flux
161     DO j=jMin,jMax
162     DO i=iMin,iMax
163     fVerT(i,j,kUp) = (afFacT*af(i,j) + dfFacT*df(i,j))*maskUp(i,j)
164     ENDDO
165     ENDDO
166    
167     C-- Tendency is minus divergence of the fluxes.
168     C Note. Tendency terms will only be correct for range
169     C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
170     C will contain valid floating point numbers but
171     C they are not algorithmically correct. These points
172     C are not used.
173     DO j=jMin,jMax
174     DO i=iMin,iMax
175     gT(i,j,k,bi,bj)=
176     & -rHFacC(i,j,k,bi,bj)*rdzF(k)*rDxF(i,j,bi,bj)*rDyF(i,j,bi,bj)
177     & *(
178     & +( fZon(i+1,j)-fZon(i,j) )
179     & +( fMer(i,j+1)-fMer(i,j) )
180     & +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )
181     & )
182     ENDDO
183     ENDDO
184    
185     C-- External thermal forcing term(s)
186    
187     RETURN
188     END

  ViewVC Help
Powered by ViewVC 1.1.22