/[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.15 - (show annotations) (download)
Wed Oct 28 03:11:36 1998 UTC (25 years, 6 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint16
Changes since 1.14: +21 -11 lines
Changes to support
 - g77 compilation under Linux
 - LR(1) form of 64-bit is D or E for constants
 - Modified adjoint of exch with adjoint variables
   acuumulated.

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

  ViewVC Help
Powered by ViewVC 1.1.22