/[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.27 - (show annotations) (download)
Wed Feb 7 16:27:29 2001 UTC (23 years, 3 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint38, c37_adj, checkpoint39, checkpoint37, checkpoint36, checkpoint35
Branch point for: pre38
Changes since 1.26: +2 -3 lines
Applied same fix for DEC in both calc_gs and calc_gt. Keeping these
to routines close together will make it easier to replace them
with a generic routine later.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gs.F,v 1.26 2001/02/06 05:12:46 cnh Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CStartOfInterFace
7 SUBROUTINE CALC_GS(
8 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
9 I xA,yA,uTrans,vTrans,rTrans,maskup,maskC,
10 I KappaRS,
11 U fVerS,
12 I myCurrentTime, myThid )
13 C /==========================================================\
14 C | SUBROUTINE CALC_GS |
15 C | o Calculate the salt tendency terms. |
16 C |==========================================================|
17 C | A procedure called EXTERNAL_FORCING_S is called from |
18 C | here. These procedures can be used to add per problem |
19 C | E-P flux source terms. |
20 C | Note: Although it is slightly counter-intuitive the |
21 C | EXTERNAL_FORCING routine is not the place to put |
22 C | file I/O. Instead files that are required to |
23 C | calculate the external source terms are generally |
24 C | read during the model main loop. This makes the |
25 C | logisitics of multi-processing simpler and also |
26 C | makes the adjoint generation simpler. It also |
27 C | allows for I/O to overlap computation where that |
28 C | is supported by hardware. |
29 C | Aside from the problem specific term the code here |
30 C | forms the tendency terms due to advection and mixing |
31 C | The baseline implementation here uses a centered |
32 C | difference form for the advection term and a tensorial |
33 C | divergence of a flux form for the diffusive term. The |
34 C | diffusive term is formulated so that isopycnal mixing and|
35 C | GM-style subgrid-scale terms can be incorporated b simply|
36 C | setting the diffusion tensor terms appropriately. |
37 C \==========================================================/
38 IMPLICIT NONE
39
40 C == GLobal variables ==
41 #include "SIZE.h"
42 #include "DYNVARS.h"
43 #include "EEPARAMS.h"
44 #include "PARAMS.h"
45 #include "GRID.h"
46 #include "FFIELDS.h"
47
48 C == Routine arguments ==
49 C fVerS - Flux of salt (S) in the vertical
50 C direction at the upper(U) and lower(D) faces of a cell.
51 C maskUp - Land mask used to denote base of the domain.
52 C maskC - Land mask for salt cells (used in TOP_LAYER only)
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 rTrans - Vertical volume transport through cell face
58 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
59 C results will be set.
60 C myThid - Instance number for this innvocation of CALC_GT
61 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
62 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
63 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
64 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
65 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
66 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
67 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
68 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
69 _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
70 INTEGER k,kUp,kDown,kM1
71 INTEGER bi,bj,iMin,iMax,jMin,jMax
72 _RL myCurrentTime
73 INTEGER myThid
74 CEndOfInterface
75
76 C == Local variables ==
77 C I, J, K - Loop counters
78 INTEGER i,j
79 LOGICAL TOP_LAYER
80 _RL afFacS, dfFacS
81 _RL dSdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RL dSdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RL df4 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87 _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88
89 #ifdef ALLOW_AUTODIFF_TAMC
90 C-- only the kUp part of fverS is set in this subroutine
91 C-- the kDown is still required
92 fVerS(1,1,kDown) = fVerS(1,1,kDown)
93 #endif
94 DO j=1-OLy,sNy+OLy
95 DO i=1-OLx,sNx+OLx
96 fZon(i,j) = 0.0
97 fMer(i,j) = 0.0
98 fVerS(i,j,kUp) = 0.0
99 ENDDO
100 ENDDO
101
102 afFacS = 1. _d 0
103 dfFacS = 1. _d 0
104 TOP_LAYER = K .EQ. 1
105
106 C--- Calculate advective and diffusive fluxes between cells.
107
108 #ifdef INCLUDE_T_DIFFUSION_CODE
109 C o Zonal tracer gradient
110 DO j=1-Oly,sNy+Oly
111 DO i=1-Olx+1,sNx+Olx
112 dSdx(i,j) = _recip_dxC(i,j,bi,bj)*
113 & (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj))
114 ENDDO
115 ENDDO
116 C o Meridional tracer gradient
117 DO j=1-Oly+1,sNy+Oly
118 DO i=1-Olx,sNx+Olx
119 dSdy(i,j) = _recip_dyC(i,j,bi,bj)*
120 & (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))
121 ENDDO
122 ENDDO
123
124 C-- del^2 of S, needed for bi-harmonic (del^4) term
125 IF (diffK4S .NE. 0.) THEN
126 DO j=1-Oly+1,sNy+Oly-1
127 DO i=1-Olx+1,sNx+Olx-1
128 df4(i,j)= _recip_hFacC(i,j,k,bi,bj)
129 & *recip_drF(k)/_rA(i,j,bi,bj)
130 & *(
131 & +( xA(i+1,j)*dSdx(i+1,j)-xA(i,j)*dSdx(i,j) )
132 & +( yA(i,j+1)*dSdy(i,j+1)-yA(i,j)*dSdy(i,j) )
133 & )
134 ENDDO
135 ENDDO
136 ENDIF
137 #endif
138
139 C-- Zonal flux (fZon is at west face of "salt" cell)
140 C Advective component of zonal flux
141 DO j=jMin,jMax
142 DO i=iMin,iMax
143 af(i,j) =
144 & uTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i-1,j,k,bi,bj))*0.5 _d 0
145 ENDDO
146 ENDDO
147 C o Diffusive component of zonal flux
148 DO j=jMin,jMax
149 DO i=iMin,iMax
150 df(i,j) = -diffKhS*xA(i,j)*dSdx(i,j)
151 ENDDO
152 ENDDO
153 #ifdef ALLOW_GMREDI
154 IF (useGMRedi) CALL GMREDI_XTRANSPORT(
155 I iMin,iMax,jMin,jMax,bi,bj,K,
156 I xA,salt,
157 U df,
158 I myThid)
159 #endif
160 C o Add the bi-harmonic contribution
161 IF (diffK4S .NE. 0.) THEN
162 DO j=jMin,jMax
163 DO i=iMin,iMax
164 df(i,j) = df(i,j) + xA(i,j)*
165 & diffK4S*(df4(i,j)-df4(i-1,j))*_recip_dxC(i,j,bi,bj)
166 ENDDO
167 ENDDO
168 ENDIF
169 C Net zonal flux
170 DO j=jMin,jMax
171 DO i=iMin,iMax
172 fZon(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)
173 ENDDO
174 ENDDO
175
176 C-- Meridional flux (fMer is at south face of "salt" cell)
177 C Advective component of meridional flux
178 DO j=jMin,jMax
179 DO i=iMin,iMax
180 C Advective component of meridional flux
181 af(i,j) =
182 & vTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j-1,k,bi,bj))*0.5 _d 0
183 ENDDO
184 ENDDO
185 C Diffusive component of meridional flux
186 DO j=jMin,jMax
187 DO i=iMin,iMax
188 df(i,j) = -diffKhS*yA(i,j)*dSdy(i,j)
189 ENDDO
190 ENDDO
191 #ifdef ALLOW_GMREDI
192 IF (useGMRedi) CALL GMREDI_YTRANSPORT(
193 I iMin,iMax,jMin,jMax,bi,bj,K,
194 I yA,salt,
195 U df,
196 I myThid)
197 #endif
198 C o Add the bi-harmonic contribution
199 IF (diffK4S .NE. 0.) THEN
200 DO j=jMin,jMax
201 DO i=iMin,iMax
202 df(i,j) = df(i,j) + yA(i,j)*
203 & diffK4S*(df4(i,j)-df4(i,j-1))*_recip_dyC(i,j,bi,bj)
204 ENDDO
205 ENDDO
206 ENDIF
207
208 C Net meridional flux
209 DO j=jMin,jMax
210 DO i=iMin,iMax
211 fMer(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)
212 ENDDO
213 ENDDO
214
215 C-- Vertical flux (fVerS) above
216 C Advective component of vertical flux
217 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
218 C (this plays the role of the free-surface correction)
219 DO j=jMin,jMax
220 DO i=iMin,iMax
221 af(i,j) =
222 & rTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j,kM1,bi,bj))*0.5 _d 0
223 ENDDO
224 ENDDO
225 C o Diffusive component of vertical flux
226 C Note: For K=1 then KM1=1 and this gives a dS/dr = 0 upper
227 C boundary condition.
228 IF (implicitDiffusion) THEN
229 DO j=jMin,jMax
230 DO i=iMin,iMax
231 df(i,j) = 0.
232 ENDDO
233 ENDDO
234 ELSE
235 DO j=jMin,jMax
236 DO i=iMin,iMax
237 df(i,j) = - _rA(i,j,bi,bj)*(
238 & KappaRS(i,j,k)*recip_drC(k)
239 & *(salt(i,j,kM1,bi,bj)-salt(i,j,k,bi,bj))*rkFac
240 & )
241 ENDDO
242 ENDDO
243 ENDIF
244
245 #ifdef ALLOW_GMREDI
246 IF (useGMRedi) CALL GMREDI_RTRANSPORT(
247 I iMin,iMax,jMin,jMax,bi,bj,K,
248 I maskUp,salt,
249 U df,
250 I myThid)
251 #endif
252
253 #ifdef ALLOW_KPP
254 C-- Add non-local KPP transport term (ghat) to diffusive salt flux.
255 IF (useKPP) CALL KPP_TRANSPORT_S(
256 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
257 I maskC,KappaRS,
258 U df )
259 #endif
260
261 C Net vertical flux
262 DO j=jMin,jMax
263 DO i=iMin,iMax
264 fVerS(i,j,kUp) = ( afFacS*af(i,j)+ dfFacS*df(i,j) )*maskUp(i,j)
265 ENDDO
266 ENDDO
267 IF ( TOP_LAYER ) THEN
268 DO j=jMin,jMax
269 DO i=iMin,iMax
270 fVerS(i,j,kUp) = afFacS*af(i,j)*freeSurfFac
271 ENDDO
272 ENDDO
273 ENDIF
274
275 C-- Tendency is minus divergence of the fluxes.
276 C Note. Tendency terms will only be correct for range
277 C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
278 C will contain valid floating point numbers but
279 C they are not algorithmically correct. These points
280 C are not used.
281 DO j=jMin,jMax-1
282 DO i=iMin,iMax-1
283 C DO j=1-2,OLy+2
284 C DO i=1-2,OLx+2
285 #define _recip_VolS1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
286 #define _recip_VolS2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
287 gS(i,j,k,bi,bj)=
288 & -_recip_VolS1(i,j,k,bi,bj)
289 & _recip_VolS2(i,j,k,bi,bj)
290 & *(
291 & +( fZon(i+1,j)-fZon(i,j) )
292 & +( fMer(i,j+1)-fMer(i,j) )
293 & +( fVerS(i,j,kUp)-fVerS(i,j,kDown) )*rkFac
294 & )
295 ENDDO
296 ENDDO
297
298 C-- External forcing term(s)
299 CALL EXTERNAL_FORCING_S(
300 I iMin,iMax,jMin,jMax,bi,bj,k,
301 I maskC,
302 I myCurrentTime,myThid)
303
304 RETURN
305 END

  ViewVC Help
Powered by ViewVC 1.1.22