/[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.18 - (hide annotations) (download)
Tue May 18 18:01:12 1999 UTC (25 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint22
Changes since 1.17: +32 -3 lines
Modifications/additions for KPP mixing scheme. Instigated by Dimitri.

1 adcroft 1.18 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gs.F,v 1.17 1998/11/06 22:44:44 cnh Exp $
2 cnh 1.1
3 cnh 1.17 #include "CPP_OPTIONS.h"
4 cnh 1.1
5     CStartOfInterFace
6     SUBROUTINE CALC_GS(
7     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
8 cnh 1.12 I xA,yA,uTrans,vTrans,rTrans,maskup,maskC,
9     I K13,K23,KappaRS,KapGM,
10 adcroft 1.7 U af,df,fZon,fMer,fVerS,
11 cnh 1.16 I myCurrentTime, myThid )
12 cnh 1.1 C /==========================================================\
13     C | SUBROUTINE CALC_GS |
14 adcroft 1.7 C | o Calculate the salt tendency terms. |
15 cnh 1.1 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 adcroft 1.7 C | E-P flux source terms. |
19 cnh 1.1 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 cnh 1.8 #include "FFIELDS.h"
46 adcroft 1.18 #ifdef ALLOW_KPP
47     #include "KPPMIX.h"
48     #endif
49 cnh 1.1
50     C == Routine arguments ==
51     C fZon - Work array for flux of temperature in the east-west
52     C direction at the west face of a cell.
53     C fMer - Work array for flux of temperature in the north-south
54     C direction at the south face of a cell.
55 adcroft 1.7 C fVerS - Flux of salt (S) in the vertical
56 cnh 1.1 C direction at the upper(U) and lower(D) faces of a cell.
57     C maskUp - Land mask used to denote base of the domain.
58 adcroft 1.10 C maskC - Land mask for salt cells (used in TOP_LAYER only)
59 cnh 1.1 C xA - Tracer cell face area normal to X
60     C yA - Tracer cell face area normal to X
61     C uTrans - Zonal volume transport through cell face
62     C vTrans - Meridional volume transport through cell face
63 adcroft 1.18 C rTrans - Vertical volume transport through cell face
64 cnh 1.1 C af - Advective flux component work array
65     C df - Diffusive flux component work array
66     C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
67     C results will be set.
68 adcroft 1.7 C myThid - Instance number for this innvocation of CALC_GT
69 cnh 1.1 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71     _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
72     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74     _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75     _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 cnh 1.12 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 cnh 1.1 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 adcroft 1.10 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 cnh 1.14 _RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
80     _RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
81     _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
82 adcroft 1.7 _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 cnh 1.1 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 adcroft 1.7 INTEGER k,kUp,kDown,kM1
86 cnh 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
87 adcroft 1.18 _RL myCurrentTime
88 cnh 1.1 INTEGER myThid
89     CEndOfInterface
90    
91     C == Local variables ==
92     C I, J, K - Loop counters
93 adcroft 1.7 INTEGER i,j
94     LOGICAL TOP_LAYER
95     _RL afFacS, dfFacS
96     _RL dSdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97     _RL dSdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 cnh 1.1
99     afFacS = 1. _d 0
100     dfFacS = 1. _d 0
101 adcroft 1.7 TOP_LAYER = K .EQ. 1
102 cnh 1.1
103     C--- Calculate advective and diffusive fluxes between cells.
104    
105     C-- Zonal flux (fZon is at west face of "salt" cell)
106     C Advective component of zonal flux
107     DO j=jMin,jMax
108     DO i=iMin,iMax
109     af(i,j) =
110     & uTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i-1,j,k,bi,bj))*0.5 _d 0
111     ENDDO
112     ENDDO
113 adcroft 1.7 C Zonal tracer gradient
114     DO j=jMin,jMax
115     DO i=iMin,iMax
116 cnh 1.12 dSdx(i,j) = _recip_dxC(i,j,bi,bj)*
117 adcroft 1.7 & (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj))
118     ENDDO
119     ENDDO
120 cnh 1.1 C Diffusive component of zonal flux
121     DO j=jMin,jMax
122     DO i=iMin,iMax
123 adcroft 1.7 df(i,j) = -(diffKhS+0.5*(KapGM(i,j)+KapGM(i-1,j)))*
124     & xA(i,j)*dSdx(i,j)
125 cnh 1.1 ENDDO
126     ENDDO
127     C Net zonal flux
128     DO j=jMin,jMax
129     DO i=iMin,iMax
130     fZon(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)
131     ENDDO
132     ENDDO
133    
134     C-- Meridional flux (fMer is at south face of "salt" cell)
135     C Advective component of meridional flux
136     DO j=jMin,jMax
137     DO i=iMin,iMax
138     C Advective component of meridional flux
139     af(i,j) =
140     & vTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j-1,k,bi,bj))*0.5 _d 0
141     ENDDO
142     ENDDO
143 adcroft 1.7 C Zonal tracer gradient
144     DO j=jMin,jMax
145     DO i=iMin,iMax
146 cnh 1.12 dSdy(i,j) = _recip_dyC(i,j,bi,bj)*
147 adcroft 1.7 & (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))
148     ENDDO
149     ENDDO
150 cnh 1.1 C Diffusive component of meridional flux
151     DO j=jMin,jMax
152     DO i=iMin,iMax
153 adcroft 1.7 df(i,j) = -(diffKhS+0.5*(KapGM(i,j)+KapGM(i,j-1)))*
154     & yA(i,j)*dSdy(i,j)
155 cnh 1.1 ENDDO
156     ENDDO
157     C Net meridional flux
158     DO j=jMin,jMax
159     DO i=iMin,iMax
160     fMer(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)
161     ENDDO
162     ENDDO
163    
164 adcroft 1.7 C-- Interpolate terms for Redi/GM scheme
165     DO j=jMin,jMax
166     DO i=iMin,iMax
167     dSdx(i,j) = 0.5*(
168 cnh 1.15 & +0.5*(_maskW(i+1,j,k,bi,bj)
169     & *_recip_dxC(i+1,j,bi,bj)*
170 adcroft 1.7 & (salt(i+1,j,k,bi,bj)-salt(i,j,k,bi,bj))
171 cnh 1.15 & +_maskW(i,j,k,bi,bj)
172     & *_recip_dxC(i,j,bi,bj)*
173 adcroft 1.7 & (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj)))
174 cnh 1.15 & +0.5*(_maskW(i+1,j,km1,bi,bj)
175     & *_recip_dxC(i+1,j,bi,bj)*
176 adcroft 1.7 & (salt(i+1,j,km1,bi,bj)-salt(i,j,km1,bi,bj))
177 cnh 1.15 & +_maskW(i,j,km1,bi,bj)
178     & *_recip_dxC(i,j,bi,bj)*
179 adcroft 1.7 & (salt(i,j,km1,bi,bj)-salt(i-1,j,km1,bi,bj)))
180     & )
181     ENDDO
182     ENDDO
183     DO j=jMin,jMax
184     DO i=iMin,iMax
185     dSdy(i,j) = 0.5*(
186 cnh 1.15 & +0.5*(_maskS(i,j,k,bi,bj)
187     & *_recip_dyC(i,j,bi,bj)*
188 adcroft 1.7 & (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))
189 cnh 1.15 & +_maskS(i,j+1,k,bi,bj)
190     & *_recip_dyC(i,j+1,bi,bj)*
191 adcroft 1.7 & (salt(i,j+1,k,bi,bj)-salt(i,j,k,bi,bj)))
192 cnh 1.15 & +0.5*(_maskS(i,j,km1,bi,bj)
193     & *_recip_dyC(i,j,bi,bj)*
194 adcroft 1.7 & (salt(i,j,km1,bi,bj)-salt(i,j-1,km1,bi,bj))
195 cnh 1.15 & +_maskS(i,j+1,km1,bi,bj)
196     & *_recip_dyC(i,j+1,bi,bj)*
197 adcroft 1.7 & (salt(i,j+1,km1,bi,bj)-salt(i,j,km1,bi,bj)))
198     & )
199     ENDDO
200     ENDDO
201    
202 cnh 1.1 C-- Vertical flux (fVerS) above
203     C Advective component of vertical flux
204 adcroft 1.7 C Note: For K=1 then KM1=1 this gives a barZ(T) = T
205     C (this plays the role of the free-surface correction)
206 cnh 1.1 DO j=jMin,jMax
207     DO i=iMin,iMax
208     af(i,j) =
209 cnh 1.12 & rTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j,kM1,bi,bj))*0.5 _d 0
210 cnh 1.1 ENDDO
211     ENDDO
212     C Diffusive component of vertical flux
213 adcroft 1.7 C Note: For K=1 then KM1=1 this gives a dS/dz = 0 upper
214     C boundary condition.
215 cnh 1.1 DO j=jMin,jMax
216     DO i=iMin,iMax
217 cnh 1.13 df(i,j) = _rA(i,j,bi,bj)*(
218 adcroft 1.7 & -KapGM(i,j)*K13(i,j,k)*dSdx(i,j)
219     & -KapGM(i,j)*K23(i,j,k)*dSdy(i,j)
220     & )
221 cnh 1.1 ENDDO
222     ENDDO
223 adcroft 1.7 IF (.NOT.implicitDiffusion) THEN
224     DO j=jMin,jMax
225     DO i=iMin,iMax
226 cnh 1.12 df(i,j) = df(i,j) + _rA(i,j,bi,bj)*(
227     & -KappaRS(i,j,k)*recip_drC(k)
228 cnh 1.13 & *(salt(i,j,kM1,bi,bj)-salt(i,j,k,bi,bj))*rkFac
229 adcroft 1.7 & )
230     ENDDO
231     ENDDO
232     ENDIF
233 adcroft 1.18 #ifdef ALLOW_KPP
234     IF (usingKPPmixing) THEN
235     C-- Add non local transport coefficient (ghat term) to right-hand-side
236     C The nonlocal transport term is noNrero only for scalars in unstable
237     C (convective) forcing conditions.
238     IF ( TOP_LAYER ) THEN
239     DO j=jMin,jMax
240     DO i=iMin,iMax
241     df(i,j) = df(i,j) - _rA(i,j,bi,bj) *
242     & EmPmR(i,j,bi,bj) * delZ(1) *
243     & ( KappaRS(i,j,k) * KPPghat(i,j,k,bi,bj) )
244     ENDDO
245     ENDDO
246     ELSE
247     DO j=jMin,jMax
248     DO i=iMin,iMax
249     df(i,j) = df(i,j) - _rA(i,j,bi,bj) *
250     & EmPmR(i,j,bi,bj) * delZ(1) *
251     & ( KappaRS(i,j,k) * KPPghat(i,j,k,bi,bj)
252     & - KappaRS(i,j,k-1) * KPPghat(i,j,k-1,bi,bj) )
253     ENDDO
254     ENDDO
255     ENDIF
256     ENDIF
257     #endif /* ALLOW_KPP */
258    
259 cnh 1.1 C Net vertical flux
260     DO j=jMin,jMax
261     DO i=iMin,iMax
262 adcroft 1.7 fVerS(i,j,kUp) = ( afFacS*af(i,j)+ dfFacS*df(i,j) )*maskUp(i,j)
263 cnh 1.1 ENDDO
264     ENDDO
265 adcroft 1.7 IF ( TOP_LAYER ) THEN
266     DO j=jMin,jMax
267     DO i=iMin,iMax
268     fVerS(i,j,kUp) = afFacS*af(i,j)*freeSurfFac
269     ENDDO
270     ENDDO
271     ENDIF
272 cnh 1.1
273     C-- Tendency is minus divergence of the fluxes.
274     C Note. Tendency terms will only be correct for range
275     C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points
276     C will contain valid floating point numbers but
277     C they are not algorithmically correct. These points
278     C are not used.
279     DO j=jMin,jMax
280     DO i=iMin,iMax
281 cnh 1.15 #define _recip_VolS1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
282     #define _recip_VolS2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
283 cnh 1.1 gS(i,j,k,bi,bj)=
284 cnh 1.15 & -_recip_VolS1(i,j,k,bi,bj)
285     & _recip_VolS2(i,j,k,bi,bj)
286 cnh 1.1 & *(
287     & +( fZon(i+1,j)-fZon(i,j) )
288     & +( fMer(i,j+1)-fMer(i,j) )
289 cnh 1.12 & +( fVerS(i,j,kUp)-fVerS(i,j,kDown) )*rkFac
290 cnh 1.1 & )
291     ENDDO
292     ENDDO
293    
294 cnh 1.16 C-- External forcing term(s)
295     CALL EXTERNAL_FORCING_S(
296     I iMin,iMax,jMin,jMax,bi,bj,k,
297 cnh 1.17 I maskC,
298 cnh 1.16 I myCurrentTime,myThid)
299    
300 cnh 1.17 #ifdef INCLUDE_LAT_CIRC_FFT_FILTER_CODE
301 cnh 1.16 C--
302     CALL FILTER_LATCIRCS_FFT_APPLY( gS, 1, sNy, k, k, bi, bj, 1, myThid)
303     #endif
304 cnh 1.1
305     RETURN
306     END

  ViewVC Help
Powered by ViewVC 1.1.22