1 |
cnh |
1.8 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gs.F,v 1.7 1998/06/10 16:05:39 adcroft Exp $ |
2 |
cnh |
1.1 |
|
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,wTrans,maskup, |
9 |
adcroft |
1.7 |
I K13,K23,KappaZS,KapGM, |
10 |
|
|
U af,df,fZon,fMer,fVerS, |
11 |
cnh |
1.1 |
I myThid ) |
12 |
|
|
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 |
cnh |
1.1 |
|
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 |
adcroft |
1.7 |
C fVerS - Flux of salt (S) in the vertical |
53 |
cnh |
1.1 |
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 xA - Tracer cell face area normal to X |
56 |
|
|
C yA - Tracer cell face area normal to X |
57 |
|
|
C uTrans - Zonal volume transport through cell face |
58 |
|
|
C vTrans - Meridional volume transport through cell face |
59 |
|
|
C wTrans - Vertical volume transport through cell face |
60 |
|
|
C af - Advective flux component work array |
61 |
|
|
C df - Diffusive flux component work array |
62 |
|
|
C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation |
63 |
|
|
C results will be set. |
64 |
adcroft |
1.7 |
C myThid - Instance number for this innvocation of CALC_GT |
65 |
cnh |
1.1 |
_RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
66 |
|
|
_RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
67 |
|
|
_RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
68 |
|
|
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
69 |
|
|
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
70 |
|
|
_RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
71 |
|
|
_RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
72 |
|
|
_RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
73 |
|
|
_RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
74 |
adcroft |
1.7 |
_RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz) |
75 |
|
|
_RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz) |
76 |
|
|
_RL KappaZS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz) |
77 |
|
|
_RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
78 |
cnh |
1.1 |
_RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
79 |
|
|
_RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
80 |
adcroft |
1.7 |
INTEGER k,kUp,kDown,kM1 |
81 |
cnh |
1.1 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
82 |
|
|
INTEGER myThid |
83 |
|
|
CEndOfInterface |
84 |
|
|
|
85 |
|
|
C == Local variables == |
86 |
|
|
C I, J, K - Loop counters |
87 |
adcroft |
1.7 |
INTEGER i,j |
88 |
|
|
LOGICAL TOP_LAYER |
89 |
|
|
_RL afFacS, dfFacS |
90 |
|
|
_RL dSdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
91 |
|
|
_RL dSdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
92 |
cnh |
1.1 |
|
93 |
|
|
afFacS = 1. _d 0 |
94 |
|
|
dfFacS = 1. _d 0 |
95 |
adcroft |
1.7 |
TOP_LAYER = K .EQ. 1 |
96 |
cnh |
1.1 |
|
97 |
|
|
C--- Calculate advective and diffusive fluxes between cells. |
98 |
|
|
|
99 |
|
|
C-- Zonal flux (fZon is at west face of "salt" cell) |
100 |
|
|
C Advective component of zonal flux |
101 |
|
|
DO j=jMin,jMax |
102 |
|
|
DO i=iMin,iMax |
103 |
|
|
af(i,j) = |
104 |
|
|
& uTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i-1,j,k,bi,bj))*0.5 _d 0 |
105 |
|
|
ENDDO |
106 |
|
|
ENDDO |
107 |
adcroft |
1.7 |
C Zonal tracer gradient |
108 |
|
|
DO j=jMin,jMax |
109 |
|
|
DO i=iMin,iMax |
110 |
|
|
dSdx(i,j) = _rdxC(i,j,bi,bj)* |
111 |
|
|
& (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj)) |
112 |
|
|
ENDDO |
113 |
|
|
ENDDO |
114 |
cnh |
1.1 |
C Diffusive component of zonal flux |
115 |
|
|
DO j=jMin,jMax |
116 |
|
|
DO i=iMin,iMax |
117 |
adcroft |
1.7 |
df(i,j) = -(diffKhS+0.5*(KapGM(i,j)+KapGM(i-1,j)))* |
118 |
|
|
& xA(i,j)*dSdx(i,j) |
119 |
cnh |
1.1 |
ENDDO |
120 |
|
|
ENDDO |
121 |
|
|
C Net zonal flux |
122 |
|
|
DO j=jMin,jMax |
123 |
|
|
DO i=iMin,iMax |
124 |
|
|
fZon(i,j) = afFacS*af(i,j) + dfFacS*df(i,j) |
125 |
|
|
ENDDO |
126 |
|
|
ENDDO |
127 |
|
|
|
128 |
|
|
C-- Meridional flux (fMer is at south face of "salt" cell) |
129 |
|
|
C Advective component of meridional flux |
130 |
|
|
DO j=jMin,jMax |
131 |
|
|
DO i=iMin,iMax |
132 |
|
|
C Advective component of meridional flux |
133 |
|
|
af(i,j) = |
134 |
|
|
& vTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j-1,k,bi,bj))*0.5 _d 0 |
135 |
|
|
ENDDO |
136 |
|
|
ENDDO |
137 |
adcroft |
1.7 |
C Zonal tracer gradient |
138 |
|
|
DO j=jMin,jMax |
139 |
|
|
DO i=iMin,iMax |
140 |
|
|
dSdy(i,j) = _rdyC(i,j,bi,bj)* |
141 |
|
|
& (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj)) |
142 |
|
|
ENDDO |
143 |
|
|
ENDDO |
144 |
cnh |
1.1 |
C Diffusive component of meridional flux |
145 |
|
|
DO j=jMin,jMax |
146 |
|
|
DO i=iMin,iMax |
147 |
adcroft |
1.7 |
df(i,j) = -(diffKhS+0.5*(KapGM(i,j)+KapGM(i,j-1)))* |
148 |
|
|
& yA(i,j)*dSdy(i,j) |
149 |
cnh |
1.1 |
ENDDO |
150 |
|
|
ENDDO |
151 |
|
|
C Net meridional flux |
152 |
|
|
DO j=jMin,jMax |
153 |
|
|
DO i=iMin,iMax |
154 |
|
|
fMer(i,j) = afFacS*af(i,j) + dfFacS*df(i,j) |
155 |
|
|
ENDDO |
156 |
|
|
ENDDO |
157 |
|
|
|
158 |
adcroft |
1.7 |
C-- Interpolate terms for Redi/GM scheme |
159 |
|
|
DO j=jMin,jMax |
160 |
|
|
DO i=iMin,iMax |
161 |
|
|
dSdx(i,j) = 0.5*( |
162 |
|
|
& +0.5*(_maskW(i+1,j,k,bi,bj)*_rdxC(i+1,j,bi,bj)* |
163 |
|
|
& (salt(i+1,j,k,bi,bj)-salt(i,j,k,bi,bj)) |
164 |
|
|
& +_maskW(i,j,k,bi,bj)*_rdxC(i,j,bi,bj)* |
165 |
|
|
& (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj))) |
166 |
|
|
& +0.5*(_maskW(i+1,j,km1,bi,bj)*_rdxC(i+1,j,bi,bj)* |
167 |
|
|
& (salt(i+1,j,km1,bi,bj)-salt(i,j,km1,bi,bj)) |
168 |
|
|
& +_maskW(i,j,km1,bi,bj)*_rdxC(i,j,bi,bj)* |
169 |
|
|
& (salt(i,j,km1,bi,bj)-salt(i-1,j,km1,bi,bj))) |
170 |
|
|
& ) |
171 |
|
|
ENDDO |
172 |
|
|
ENDDO |
173 |
|
|
DO j=jMin,jMax |
174 |
|
|
DO i=iMin,iMax |
175 |
|
|
dSdy(i,j) = 0.5*( |
176 |
|
|
& +0.5*(_maskS(i,j,k,bi,bj)*_rdyC(i,j,bi,bj)* |
177 |
|
|
& (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj)) |
178 |
|
|
& +_maskS(i,j+1,k,bi,bj)*_rdyC(i,j+1,bi,bj)* |
179 |
|
|
& (salt(i,j+1,k,bi,bj)-salt(i,j,k,bi,bj))) |
180 |
|
|
& +0.5*(_maskS(i,j,km1,bi,bj)*_rdyC(i,j,bi,bj)* |
181 |
|
|
& (salt(i,j,km1,bi,bj)-salt(i,j-1,km1,bi,bj)) |
182 |
|
|
& +_maskS(i,j+1,km1,bi,bj)*_rdyC(i,j+1,bi,bj)* |
183 |
|
|
& (salt(i,j+1,km1,bi,bj)-salt(i,j,km1,bi,bj))) |
184 |
|
|
& ) |
185 |
|
|
ENDDO |
186 |
|
|
ENDDO |
187 |
|
|
|
188 |
cnh |
1.1 |
C-- Vertical flux (fVerS) above |
189 |
|
|
C Advective component of vertical flux |
190 |
adcroft |
1.7 |
C Note: For K=1 then KM1=1 this gives a barZ(T) = T |
191 |
|
|
C (this plays the role of the free-surface correction) |
192 |
cnh |
1.1 |
DO j=jMin,jMax |
193 |
|
|
DO i=iMin,iMax |
194 |
|
|
af(i,j) = |
195 |
|
|
& wTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j,kM1,bi,bj))*0.5 _d 0 |
196 |
|
|
ENDDO |
197 |
|
|
ENDDO |
198 |
|
|
C Diffusive component of vertical flux |
199 |
adcroft |
1.7 |
C Note: For K=1 then KM1=1 this gives a dS/dz = 0 upper |
200 |
|
|
C boundary condition. |
201 |
cnh |
1.1 |
DO j=jMin,jMax |
202 |
|
|
DO i=iMin,iMax |
203 |
adcroft |
1.7 |
df(i,j) = _zA(i,j,bi,bj)*( |
204 |
|
|
& -KapGM(i,j)*K13(i,j,k)*dSdx(i,j) |
205 |
|
|
& -KapGM(i,j)*K23(i,j,k)*dSdy(i,j) |
206 |
|
|
& ) |
207 |
cnh |
1.1 |
ENDDO |
208 |
|
|
ENDDO |
209 |
adcroft |
1.7 |
IF (.NOT.implicitDiffusion) THEN |
210 |
|
|
DO j=jMin,jMax |
211 |
|
|
DO i=iMin,iMax |
212 |
|
|
df(i,j) = df(i,j) + _zA(i,j,bi,bj)*( |
213 |
|
|
& -KappaZS(i,j,k)*rdzC(k) |
214 |
|
|
& *(salt(i,j,kM1,bi,bj)-salt(i,j,k,bi,bj)) |
215 |
|
|
& ) |
216 |
|
|
ENDDO |
217 |
|
|
ENDDO |
218 |
|
|
ENDIF |
219 |
cnh |
1.1 |
C Net vertical flux |
220 |
|
|
DO j=jMin,jMax |
221 |
|
|
DO i=iMin,iMax |
222 |
adcroft |
1.7 |
fVerS(i,j,kUp) = ( afFacS*af(i,j)+ dfFacS*df(i,j) )*maskUp(i,j) |
223 |
cnh |
1.1 |
ENDDO |
224 |
|
|
ENDDO |
225 |
adcroft |
1.7 |
IF ( TOP_LAYER ) THEN |
226 |
|
|
DO j=jMin,jMax |
227 |
|
|
DO i=iMin,iMax |
228 |
|
|
fVerS(i,j,kUp) = afFacS*af(i,j)*freeSurfFac |
229 |
|
|
ENDDO |
230 |
|
|
ENDDO |
231 |
|
|
ENDIF |
232 |
cnh |
1.1 |
|
233 |
|
|
C-- Tendency is minus divergence of the fluxes. |
234 |
|
|
C Note. Tendency terms will only be correct for range |
235 |
|
|
C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points |
236 |
|
|
C will contain valid floating point numbers but |
237 |
|
|
C they are not algorithmically correct. These points |
238 |
|
|
C are not used. |
239 |
|
|
DO j=jMin,jMax |
240 |
|
|
DO i=iMin,iMax |
241 |
adcroft |
1.7 |
C & -_rhFacC(i,j,k,bi,bj)*rdzF(k)*_rdxF(i,j,bi,bj)*_rdyF(i,j,bi,bj) |
242 |
|
|
C & -_rhFacC(i,j,k,bi,bj)*rdzF(k)/_zA(i,j,bi,bj) |
243 |
|
|
C #define _rVolS(i,j,k,bi,bj) _rhFacC(i,j,k,bi,bj)*rdzF(k)*_rdxF(i,j,bi,bj)*_rdyF(i,j,bi,bj) |
244 |
|
|
#define _rVolS(i,j,k,bi,bj) _rhFacC(i,j,k,bi,bj)*rdzF(k)/_zA(i,j,bi,bj) |
245 |
cnh |
1.1 |
gS(i,j,k,bi,bj)= |
246 |
adcroft |
1.7 |
& -_rVolS(i,j,k,bi,bj) |
247 |
cnh |
1.1 |
& *( |
248 |
|
|
& +( fZon(i+1,j)-fZon(i,j) ) |
249 |
|
|
& +( fMer(i,j+1)-fMer(i,j) ) |
250 |
|
|
& +( fVerS(i,j,kUp)-fVerS(i,j,kDown) ) |
251 |
|
|
& ) |
252 |
|
|
ENDDO |
253 |
|
|
ENDDO |
254 |
|
|
|
255 |
adcroft |
1.7 |
C-- External P-E forcing term(s) |
256 |
cnh |
1.8 |
C o Surface relaxation term |
257 |
|
|
IF ( TOP_LAYER ) THEN |
258 |
|
|
DO j=jMin,jMax |
259 |
|
|
DO i=iMin,iMax |
260 |
|
|
gS(i,j,k,bi,bj)=gS(i,j,k,bi,bj) |
261 |
|
|
& -lambdaSaltClimRelax*(salt(i,j,k,bi,bj)-saltClim(i,j,k,bi,bj)) |
262 |
|
|
ENDDO |
263 |
|
|
ENDDO |
264 |
|
|
ENDIF |
265 |
|
|
|
266 |
cnh |
1.1 |
|
267 |
|
|
RETURN |
268 |
|
|
END |