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

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

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


Revision 1.24 - (hide annotations) (download)
Sun Feb 4 14:38:45 2001 UTC (23 years, 3 months ago) by cnh
Branch: MAIN
Changes since 1.23: +2 -1 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

1 cnh 1.24 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gs.F,v 1.23 2001/02/02 21:04:47 adcroft Exp $
2     C $Name: $
3 cnh 1.1
4 cnh 1.17 #include "CPP_OPTIONS.h"
5 cnh 1.1
6     CStartOfInterFace
7     SUBROUTINE CALC_GS(
8     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
9 cnh 1.12 I xA,yA,uTrans,vTrans,rTrans,maskup,maskC,
10 adcroft 1.21 I KappaRS,
11 adcroft 1.23 U fVerS,
12 cnh 1.16 I myCurrentTime, myThid )
13 cnh 1.1 C /==========================================================\
14     C | SUBROUTINE CALC_GS |
15 adcroft 1.7 C | o Calculate the salt tendency terms. |
16 cnh 1.1 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 adcroft 1.7 C | E-P flux source terms. |
20 cnh 1.1 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 cnh 1.8 #include "FFIELDS.h"
47 cnh 1.1
48     C == Routine arguments ==
49 adcroft 1.7 C fVerS - Flux of salt (S) in the vertical
50 cnh 1.1 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 adcroft 1.10 C maskC - Land mask for salt cells (used in TOP_LAYER only)
53 cnh 1.1 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 adcroft 1.18 C rTrans - Vertical volume transport through cell face
58 cnh 1.1 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
59     C results will be set.
60 adcroft 1.7 C myThid - Instance number for this innvocation of CALC_GT
61 cnh 1.1 _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 cnh 1.12 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
67 cnh 1.1 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
68 adcroft 1.10 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
69 cnh 1.14 _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
70 adcroft 1.7 INTEGER k,kUp,kDown,kM1
71 cnh 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
72 adcroft 1.18 _RL myCurrentTime
73 cnh 1.1 INTEGER myThid
74     CEndOfInterface
75    
76     C == Local variables ==
77     C I, J, K - Loop counters
78 adcroft 1.7 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 adcroft 1.19 _RL df4 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 adcroft 1.23 _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 heimbach 1.20
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 cnh 1.1
103     afFacS = 1. _d 0
104     dfFacS = 1. _d 0
105 adcroft 1.7 TOP_LAYER = K .EQ. 1
106 cnh 1.1
107     C--- Calculate advective and diffusive fluxes between cells.
108    
109 adcroft 1.19 #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 cnh 1.1 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 adcroft 1.19 C o Diffusive component of zonal flux
149 cnh 1.1 DO j=jMin,jMax
150     DO i=iMin,iMax
151 adcroft 1.21 df(i,j) = -diffKhS*xA(i,j)*dSdx(i,j)
152 cnh 1.1 ENDDO
153     ENDDO
154 adcroft 1.21 #ifdef ALLOW_GMREDI
155 heimbach 1.22 IF (useGMRedi) CALL GMREDI_XTRANSPORT(
156 adcroft 1.21 I iMin,iMax,jMin,jMax,bi,bj,K,
157     I xA,salt,
158     U df,
159     I myThid)
160     #endif
161 adcroft 1.19 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 cnh 1.1 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 adcroft 1.21 df(i,j) = -diffKhS*yA(i,j)*dSdy(i,j)
190 cnh 1.1 ENDDO
191     ENDDO
192 adcroft 1.21 #ifdef ALLOW_GMREDI
193 heimbach 1.22 IF (useGMRedi) CALL GMREDI_YTRANSPORT(
194 adcroft 1.21 I iMin,iMax,jMin,jMax,bi,bj,K,
195     I yA,salt,
196     U df,
197     I myThid)
198     #endif
199 adcroft 1.19 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 cnh 1.1 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 adcroft 1.7 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 cnh 1.1 DO j=jMin,jMax
221     DO i=iMin,iMax
222     af(i,j) =
223 cnh 1.12 & rTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j,kM1,bi,bj))*0.5 _d 0
224 cnh 1.1 ENDDO
225     ENDDO
226 adcroft 1.21 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 adcroft 1.7 C boundary condition.
229 adcroft 1.21 IF (implicitDiffusion) THEN
230     DO j=jMin,jMax
231     DO i=iMin,iMax
232     df(i,j) = 0.
233     ENDDO
234 cnh 1.1 ENDDO
235 adcroft 1.21 ELSE
236 adcroft 1.7 DO j=jMin,jMax
237     DO i=iMin,iMax
238 adcroft 1.21 df(i,j) = - _rA(i,j,bi,bj)*(
239     & KappaRS(i,j,k)*recip_drC(k)
240 cnh 1.13 & *(salt(i,j,kM1,bi,bj)-salt(i,j,k,bi,bj))*rkFac
241 adcroft 1.7 & )
242     ENDDO
243     ENDDO
244     ENDIF
245 adcroft 1.21
246     #ifdef ALLOW_GMREDI
247 heimbach 1.22 IF (useGMRedi) CALL GMREDI_RTRANSPORT(
248 adcroft 1.21 I iMin,iMax,jMin,jMax,bi,bj,K,
249     I maskUp,salt,
250     U df,
251     I myThid)
252     #endif
253    
254 adcroft 1.18 #ifdef ALLOW_KPP
255 adcroft 1.21 C-- Add non-local KPP transport term (ghat) to diffusive salt flux.
256 heimbach 1.22 IF (useKPP) CALL KPP_TRANSPORT_S(
257 adcroft 1.21 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
258     I maskC,KappaRS,
259     U df )
260     #endif
261 adcroft 1.18
262 cnh 1.1 C Net vertical flux
263     DO j=jMin,jMax
264     DO i=iMin,iMax
265 adcroft 1.7 fVerS(i,j,kUp) = ( afFacS*af(i,j)+ dfFacS*df(i,j) )*maskUp(i,j)
266 cnh 1.1 ENDDO
267     ENDDO
268 adcroft 1.7 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 cnh 1.1
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     DO j=jMin,jMax
283     DO i=iMin,iMax
284 cnh 1.15 #define _recip_VolS1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
285     #define _recip_VolS2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
286 cnh 1.1 gS(i,j,k,bi,bj)=
287 cnh 1.15 & -_recip_VolS1(i,j,k,bi,bj)
288     & _recip_VolS2(i,j,k,bi,bj)
289 cnh 1.1 & *(
290     & +( fZon(i+1,j)-fZon(i,j) )
291     & +( fMer(i,j+1)-fMer(i,j) )
292 cnh 1.12 & +( fVerS(i,j,kUp)-fVerS(i,j,kDown) )*rkFac
293 cnh 1.1 & )
294     ENDDO
295     ENDDO
296    
297 cnh 1.16 C-- External forcing term(s)
298     CALL EXTERNAL_FORCING_S(
299     I iMin,iMax,jMin,jMax,bi,bj,k,
300 cnh 1.17 I maskC,
301 cnh 1.16 I myCurrentTime,myThid)
302 cnh 1.1
303     RETURN
304     END

  ViewVC Help
Powered by ViewVC 1.1.22