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

Contents of /MITgcm/model/src/calc_gs.F

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


Revision 1.17 - (show annotations) (download)
Fri Nov 6 22:44:44 1998 UTC (25 years, 6 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint19, checkpoint18, checkpoint20, checkpoint21
Changes since 1.16: +4 -14 lines
Changes to allow for atmospheric integration builds of the code

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gs.F,v 1.16 1998/11/03 15:28:04 cnh Exp $
2
3 #include "CPP_OPTIONS.h"
4
5 CStartOfInterFace
6 SUBROUTINE CALC_GS(
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,KappaRS,KapGM,
10 U af,df,fZon,fMer,fVerS,
11 I myCurrentTime, myThid )
12 C /==========================================================\
13 C | SUBROUTINE CALC_GS |
14 C | o Calculate the salt tendency terms. |
15 C |==========================================================|
16 C | A procedure called EXTERNAL_FORCING_S is called from |
17 C | here. These procedures can be used to add per problem |
18 C | E-P 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
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 fVerS - Flux of salt (S) 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 C maskC - Land mask for salt cells (used in TOP_LAYER only)
56 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 C wTrans - Vertical volume transport through cell face
61 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 fVerS (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 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _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 KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
79 _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 INTEGER k,kUp,kDown,kM1
83 INTEGER bi,bj,iMin,iMax,jMin,jMax
84 INTEGER myThid
85 _RL myCurrentTime
86 CEndOfInterface
87
88 C == Local variables ==
89 C I, J, K - Loop counters
90 INTEGER i,j
91 LOGICAL TOP_LAYER
92 _RL afFacS, dfFacS
93 _RL dSdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94 _RL dSdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95
96 afFacS = 1. _d 0
97 dfFacS = 1. _d 0
98 TOP_LAYER = K .EQ. 1
99
100 C--- Calculate advective and diffusive fluxes between cells.
101
102 C-- Zonal flux (fZon is at west face of "salt" cell)
103 C Advective component of zonal flux
104 DO j=jMin,jMax
105 DO i=iMin,iMax
106 af(i,j) =
107 & uTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i-1,j,k,bi,bj))*0.5 _d 0
108 ENDDO
109 ENDDO
110 C Zonal tracer gradient
111 DO j=jMin,jMax
112 DO i=iMin,iMax
113 dSdx(i,j) = _recip_dxC(i,j,bi,bj)*
114 & (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj))
115 ENDDO
116 ENDDO
117 C Diffusive component of zonal flux
118 DO j=jMin,jMax
119 DO i=iMin,iMax
120 df(i,j) = -(diffKhS+0.5*(KapGM(i,j)+KapGM(i-1,j)))*
121 & xA(i,j)*dSdx(i,j)
122 ENDDO
123 ENDDO
124 C Net zonal flux
125 DO j=jMin,jMax
126 DO i=iMin,iMax
127 fZon(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)
128 ENDDO
129 ENDDO
130
131 C-- Meridional flux (fMer is at south face of "salt" cell)
132 C Advective component of meridional flux
133 DO j=jMin,jMax
134 DO i=iMin,iMax
135 C Advective component of meridional flux
136 af(i,j) =
137 & vTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j-1,k,bi,bj))*0.5 _d 0
138 ENDDO
139 ENDDO
140 C Zonal tracer gradient
141 DO j=jMin,jMax
142 DO i=iMin,iMax
143 dSdy(i,j) = _recip_dyC(i,j,bi,bj)*
144 & (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))
145 ENDDO
146 ENDDO
147 C Diffusive component of meridional flux
148 DO j=jMin,jMax
149 DO i=iMin,iMax
150 df(i,j) = -(diffKhS+0.5*(KapGM(i,j)+KapGM(i,j-1)))*
151 & yA(i,j)*dSdy(i,j)
152 ENDDO
153 ENDDO
154 C Net meridional flux
155 DO j=jMin,jMax
156 DO i=iMin,iMax
157 fMer(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)
158 ENDDO
159 ENDDO
160
161 C-- Interpolate terms for Redi/GM scheme
162 DO j=jMin,jMax
163 DO i=iMin,iMax
164 dSdx(i,j) = 0.5*(
165 & +0.5*(_maskW(i+1,j,k,bi,bj)
166 & *_recip_dxC(i+1,j,bi,bj)*
167 & (salt(i+1,j,k,bi,bj)-salt(i,j,k,bi,bj))
168 & +_maskW(i,j,k,bi,bj)
169 & *_recip_dxC(i,j,bi,bj)*
170 & (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj)))
171 & +0.5*(_maskW(i+1,j,km1,bi,bj)
172 & *_recip_dxC(i+1,j,bi,bj)*
173 & (salt(i+1,j,km1,bi,bj)-salt(i,j,km1,bi,bj))
174 & +_maskW(i,j,km1,bi,bj)
175 & *_recip_dxC(i,j,bi,bj)*
176 & (salt(i,j,km1,bi,bj)-salt(i-1,j,km1,bi,bj)))
177 & )
178 ENDDO
179 ENDDO
180 DO j=jMin,jMax
181 DO i=iMin,iMax
182 dSdy(i,j) = 0.5*(
183 & +0.5*(_maskS(i,j,k,bi,bj)
184 & *_recip_dyC(i,j,bi,bj)*
185 & (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))
186 & +_maskS(i,j+1,k,bi,bj)
187 & *_recip_dyC(i,j+1,bi,bj)*
188 & (salt(i,j+1,k,bi,bj)-salt(i,j,k,bi,bj)))
189 & +0.5*(_maskS(i,j,km1,bi,bj)
190 & *_recip_dyC(i,j,bi,bj)*
191 & (salt(i,j,km1,bi,bj)-salt(i,j-1,km1,bi,bj))
192 & +_maskS(i,j+1,km1,bi,bj)
193 & *_recip_dyC(i,j+1,bi,bj)*
194 & (salt(i,j+1,km1,bi,bj)-salt(i,j,km1,bi,bj)))
195 & )
196 ENDDO
197 ENDDO
198
199 C-- Vertical flux (fVerS) above
200 C Advective component of vertical flux
201 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
202 C (this plays the role of the free-surface correction)
203 DO j=jMin,jMax
204 DO i=iMin,iMax
205 af(i,j) =
206 & rTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j,kM1,bi,bj))*0.5 _d 0
207 ENDDO
208 ENDDO
209 C Diffusive component of vertical flux
210 C Note: For K=1 then KM1=1 this gives a dS/dz = 0 upper
211 C boundary condition.
212 DO j=jMin,jMax
213 DO i=iMin,iMax
214 df(i,j) = _rA(i,j,bi,bj)*(
215 & -KapGM(i,j)*K13(i,j,k)*dSdx(i,j)
216 & -KapGM(i,j)*K23(i,j,k)*dSdy(i,j)
217 & )
218 ENDDO
219 ENDDO
220 IF (.NOT.implicitDiffusion) THEN
221 DO j=jMin,jMax
222 DO i=iMin,iMax
223 df(i,j) = df(i,j) + _rA(i,j,bi,bj)*(
224 & -KappaRS(i,j,k)*recip_drC(k)
225 & *(salt(i,j,kM1,bi,bj)-salt(i,j,k,bi,bj))*rkFac
226 & )
227 ENDDO
228 ENDDO
229 ENDIF
230 C Net vertical flux
231 DO j=jMin,jMax
232 DO i=iMin,iMax
233 fVerS(i,j,kUp) = ( afFacS*af(i,j)+ dfFacS*df(i,j) )*maskUp(i,j)
234 ENDDO
235 ENDDO
236 IF ( TOP_LAYER ) THEN
237 DO j=jMin,jMax
238 DO i=iMin,iMax
239 fVerS(i,j,kUp) = afFacS*af(i,j)*freeSurfFac
240 ENDDO
241 ENDDO
242 ENDIF
243
244 C-- Tendency is minus divergence of the fluxes.
245 C Note. Tendency terms will only be correct for range
246 C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
247 C will contain valid floating point numbers but
248 C they are not algorithmically correct. These points
249 C are not used.
250 DO j=jMin,jMax
251 DO i=iMin,iMax
252 #define _recip_VolS1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
253 #define _recip_VolS2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
254 gS(i,j,k,bi,bj)=
255 & -_recip_VolS1(i,j,k,bi,bj)
256 & _recip_VolS2(i,j,k,bi,bj)
257 & *(
258 & +( fZon(i+1,j)-fZon(i,j) )
259 & +( fMer(i,j+1)-fMer(i,j) )
260 & +( fVerS(i,j,kUp)-fVerS(i,j,kDown) )*rkFac
261 & )
262 ENDDO
263 ENDDO
264
265 C-- External forcing term(s)
266 CALL EXTERNAL_FORCING_S(
267 I iMin,iMax,jMin,jMax,bi,bj,k,
268 I maskC,
269 I myCurrentTime,myThid)
270
271 #ifdef INCLUDE_LAT_CIRC_FFT_FILTER_CODE
272 C--
273 CALL FILTER_LATCIRCS_FFT_APPLY( gS, 1, sNy, k, k, bi, bj, 1, myThid)
274 #endif
275
276 RETURN
277 END

  ViewVC Help
Powered by ViewVC 1.1.22