/[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.25 - (show annotations) (download)
Tue Feb 6 04:47:51 2001 UTC (23 years, 3 months ago) by cnh
Branch: MAIN
Changes since 1.24: +6 -4 lines
More references to unset array element fixes

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gs.F,v 1.24 2001/02/04 14:38:45 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
93 fVerS(1,1,kDown) = fVerS(1,1,kDown)
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 #endif
102
103 afFacS = 1. _d 0
104 dfFacS = 1. _d 0
105 TOP_LAYER = K .EQ. 1
106
107 C--- Calculate advective and diffusive fluxes between cells.
108
109 #ifdef INCLUDE_T_DIFFUSION_CODE
110 C o Zonal tracer gradient
111 DO j=1-Oly,sNy+Oly
112 DO i=1-Olx+1,sNx+Olx
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 o Meridional tracer gradient
118 DO j=1-Oly+1,sNy+Oly
119 DO i=1-Olx,sNx+Olx
120 dSdy(i,j) = _recip_dyC(i,j,bi,bj)*
121 & (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))
122 ENDDO
123 ENDDO
124
125 C-- del^2 of S, needed for bi-harmonic (del^4) term
126 IF (diffK4S .NE. 0.) THEN
127 DO j=1-Oly+1,sNy+Oly-1
128 DO i=1-Olx+1,sNx+Olx-1
129 df4(i,j)= _recip_hFacC(i,j,k,bi,bj)
130 & *recip_drF(k)/_rA(i,j,bi,bj)
131 & *(
132 & +( xA(i+1,j)*dSdx(i+1,j)-xA(i,j)*dSdx(i,j) )
133 & +( yA(i,j+1)*dSdy(i,j+1)-yA(i,j)*dSdy(i,j) )
134 & )
135 ENDDO
136 ENDDO
137 ENDIF
138 #endif
139
140 C-- Zonal flux (fZon is at west face of "salt" cell)
141 C Advective component of zonal flux
142 DO j=jMin,jMax
143 DO i=iMin,iMax
144 af(i,j) =
145 & uTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i-1,j,k,bi,bj))*0.5 _d 0
146 ENDDO
147 ENDDO
148 C o Diffusive component of zonal flux
149 DO j=jMin,jMax
150 DO i=iMin,iMax
151 df(i,j) = -diffKhS*xA(i,j)*dSdx(i,j)
152 ENDDO
153 ENDDO
154 #ifdef ALLOW_GMREDI
155 IF (useGMRedi) CALL GMREDI_XTRANSPORT(
156 I iMin,iMax,jMin,jMax,bi,bj,K,
157 I xA,salt,
158 U df,
159 I myThid)
160 #endif
161 C o Add the bi-harmonic contribution
162 IF (diffK4S .NE. 0.) THEN
163 DO j=jMin,jMax
164 DO i=iMin,iMax
165 df(i,j) = df(i,j) + xA(i,j)*
166 & diffK4S*(df4(i,j)-df4(i-1,j))*_recip_dxC(i,j,bi,bj)
167 ENDDO
168 ENDDO
169 ENDIF
170 C Net zonal flux
171 DO j=jMin,jMax
172 DO i=iMin,iMax
173 fZon(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)
174 ENDDO
175 ENDDO
176
177 C-- Meridional flux (fMer is at south face of "salt" cell)
178 C Advective component of meridional flux
179 DO j=jMin,jMax
180 DO i=iMin,iMax
181 C Advective component of meridional flux
182 af(i,j) =
183 & vTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j-1,k,bi,bj))*0.5 _d 0
184 ENDDO
185 ENDDO
186 C Diffusive component of meridional flux
187 DO j=jMin,jMax
188 DO i=iMin,iMax
189 df(i,j) = -diffKhS*yA(i,j)*dSdy(i,j)
190 ENDDO
191 ENDDO
192 #ifdef ALLOW_GMREDI
193 IF (useGMRedi) CALL GMREDI_YTRANSPORT(
194 I iMin,iMax,jMin,jMax,bi,bj,K,
195 I yA,salt,
196 U df,
197 I myThid)
198 #endif
199 C o Add the bi-harmonic contribution
200 IF (diffK4S .NE. 0.) THEN
201 DO j=jMin,jMax
202 DO i=iMin,iMax
203 df(i,j) = df(i,j) + yA(i,j)*
204 & diffK4S*(df4(i,j)-df4(i,j-1))*_recip_dyC(i,j,bi,bj)
205 ENDDO
206 ENDDO
207 ENDIF
208
209 C Net meridional flux
210 DO j=jMin,jMax
211 DO i=iMin,iMax
212 fMer(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)
213 ENDDO
214 ENDDO
215
216 C-- Vertical flux (fVerS) above
217 C Advective component of vertical flux
218 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
219 C (this plays the role of the free-surface correction)
220 DO j=jMin,jMax
221 DO i=iMin,iMax
222 af(i,j) =
223 & rTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j,kM1,bi,bj))*0.5 _d 0
224 ENDDO
225 ENDDO
226 C o Diffusive component of vertical flux
227 C Note: For K=1 then KM1=1 and this gives a dS/dr = 0 upper
228 C boundary condition.
229 IF (implicitDiffusion) THEN
230 DO j=jMin,jMax
231 DO i=iMin,iMax
232 df(i,j) = 0.
233 ENDDO
234 ENDDO
235 ELSE
236 DO j=jMin,jMax
237 DO i=iMin,iMax
238 df(i,j) = - _rA(i,j,bi,bj)*(
239 & KappaRS(i,j,k)*recip_drC(k)
240 & *(salt(i,j,kM1,bi,bj)-salt(i,j,k,bi,bj))*rkFac
241 & )
242 ENDDO
243 ENDDO
244 ENDIF
245
246 #ifdef ALLOW_GMREDI
247 IF (useGMRedi) CALL GMREDI_RTRANSPORT(
248 I iMin,iMax,jMin,jMax,bi,bj,K,
249 I maskUp,salt,
250 U df,
251 I myThid)
252 #endif
253
254 #ifdef ALLOW_KPP
255 C-- Add non-local KPP transport term (ghat) to diffusive salt flux.
256 IF (useKPP) CALL KPP_TRANSPORT_S(
257 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
258 I maskC,KappaRS,
259 U df )
260 #endif
261
262 C Net vertical flux
263 DO j=jMin,jMax
264 DO i=iMin,iMax
265 fVerS(i,j,kUp) = ( afFacS*af(i,j)+ dfFacS*df(i,j) )*maskUp(i,j)
266 ENDDO
267 ENDDO
268 IF ( TOP_LAYER ) THEN
269 DO j=jMin,jMax
270 DO i=iMin,iMax
271 fVerS(i,j,kUp) = afFacS*af(i,j)*freeSurfFac
272 ENDDO
273 ENDDO
274 ENDIF
275
276 C-- Tendency is minus divergence of the fluxes.
277 C Note. Tendency terms will only be correct for range
278 C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
279 C will contain valid floating point numbers but
280 C they are not algorithmically correct. These points
281 C are not used.
282 C DO j=jMin,jMax
283 C DO i=iMin,iMax
284 DO j=1-1,OLy+1
285 DO i=1-1,OLx+1
286 #define _recip_VolS1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
287 #define _recip_VolS2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
288 gS(i,j,k,bi,bj)=
289 & -_recip_VolS1(i,j,k,bi,bj)
290 & _recip_VolS2(i,j,k,bi,bj)
291 & *(
292 & +( fZon(i+1,j)-fZon(i,j) )
293 & +( fMer(i,j+1)-fMer(i,j) )
294 & +( fVerS(i,j,kUp)-fVerS(i,j,kDown) )*rkFac
295 & )
296 ENDDO
297 ENDDO
298
299 C-- External forcing term(s)
300 CALL EXTERNAL_FORCING_S(
301 I iMin,iMax,jMin,jMax,bi,bj,k,
302 I maskC,
303 I myCurrentTime,myThid)
304
305 RETURN
306 END

  ViewVC Help
Powered by ViewVC 1.1.22