1 |
cnh |
1.4 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gs.F,v 1.3 1998/05/27 21:01:47 cnh 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 |
|
|
U af,df,fZon,fMer, fVerS, |
10 |
|
|
I myThid ) |
11 |
|
|
C /==========================================================\ |
12 |
|
|
C | SUBROUTINE CALC_GS | |
13 |
|
|
C | o Calculate the salinity tendency terms. | |
14 |
|
|
C |==========================================================| |
15 |
|
|
C | A procedure called EXTERNAL_FORCING_S is called from | |
16 |
|
|
C | here. These procedures can be used to add per problem | |
17 |
|
|
C | fresh water flux source terms. | |
18 |
|
|
C | Note: Although it is slightly counter-intuitive the | |
19 |
|
|
C | EXTERNAL_FORCING routine is not the place to put | |
20 |
|
|
C | file I/O. Instead files that are required to | |
21 |
|
|
C | calculate the external source terms are generally | |
22 |
|
|
C | read during the model main loop. This makes the | |
23 |
|
|
C | logisitics of multi-processing simpler and also | |
24 |
|
|
C | makes the adjoint generation simpler. It also | |
25 |
|
|
C | allows for I/O to overlap computation where that | |
26 |
|
|
C | is supported by hardware. | |
27 |
|
|
C | Aside from the problem specific term the code here | |
28 |
|
|
C | forms the tendency terms due to advection and mixing | |
29 |
|
|
C | The baseline implementation here uses a centered | |
30 |
|
|
C | difference form for the advection term and a tensorial | |
31 |
|
|
C | divergence of a flux form for the diffusive term. The | |
32 |
|
|
C | diffusive term is formulated so that isopycnal mixing and| |
33 |
|
|
C | GM-style subgrid-scale terms can be incorporated b simply| |
34 |
|
|
C | setting the diffusion tensor terms appropriately. | |
35 |
|
|
C \==========================================================/ |
36 |
|
|
IMPLICIT NONE |
37 |
|
|
|
38 |
|
|
C == GLobal variables == |
39 |
|
|
#include "SIZE.h" |
40 |
|
|
#include "DYNVARS.h" |
41 |
|
|
#include "EEPARAMS.h" |
42 |
|
|
#include "PARAMS.h" |
43 |
|
|
#include "GRID.h" |
44 |
|
|
|
45 |
|
|
C == Routine arguments == |
46 |
|
|
C fZon - Work array for flux of temperature in the east-west |
47 |
|
|
C direction at the west face of a cell. |
48 |
|
|
C fMer - Work array for flux of temperature in the north-south |
49 |
|
|
C direction at the south face of a cell. |
50 |
|
|
C fVerS - Flux of salinity (S) in the vertical |
51 |
|
|
C direction at the upper(U) and lower(D) faces of a cell. |
52 |
|
|
C maskUp - Land mask used to denote base of the domain. |
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 wTrans - Vertical volume transport through cell face |
58 |
|
|
C af - Advective flux component work array |
59 |
|
|
C df - Diffusive flux component work array |
60 |
|
|
C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation |
61 |
|
|
C results will be set. |
62 |
|
|
C myThid - Instance number for this innvocation of CALC_GS |
63 |
|
|
_RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
64 |
|
|
_RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
65 |
|
|
_RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
66 |
|
|
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
67 |
|
|
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
68 |
|
|
_RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
69 |
|
|
_RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
70 |
|
|
_RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
71 |
|
|
_RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
72 |
|
|
_RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
73 |
|
|
_RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
74 |
|
|
INTEGER kUp,kDown,kM1 |
75 |
|
|
INTEGER bi,bj,iMin,iMax,jMin,jMax |
76 |
|
|
INTEGER myThid |
77 |
|
|
CEndOfInterface |
78 |
|
|
|
79 |
|
|
C == Local variables == |
80 |
|
|
C I, J, K - Loop counters |
81 |
|
|
INTEGER i,j,k |
82 |
|
|
INTEGER afFacS, dfFacS |
83 |
|
|
|
84 |
|
|
afFacS = 1. _d 0 |
85 |
|
|
dfFacS = 1. _d 0 |
86 |
|
|
|
87 |
|
|
C--- |
88 |
|
|
C--- Calculate advective and diffusive fluxes between cells. |
89 |
|
|
C--- |
90 |
|
|
|
91 |
|
|
C-- Zonal flux (fZon is at west face of "salt" cell) |
92 |
|
|
C Advective component of zonal flux |
93 |
|
|
DO j=jMin,jMax |
94 |
|
|
DO i=iMin,iMax |
95 |
|
|
af(i,j) = |
96 |
|
|
& uTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i-1,j,k,bi,bj))*0.5 _d 0 |
97 |
|
|
ENDDO |
98 |
|
|
ENDDO |
99 |
|
|
C Diffusive component of zonal flux |
100 |
|
|
DO j=jMin,jMax |
101 |
|
|
DO i=iMin,iMax |
102 |
|
|
df(i,j) = |
103 |
cnh |
1.3 |
& -diffKhS*xA(i,j)*_rdxC(i,j,bi,bj) |
104 |
cnh |
1.1 |
& *(salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj)) |
105 |
|
|
ENDDO |
106 |
|
|
ENDDO |
107 |
|
|
C Net zonal flux |
108 |
|
|
DO j=jMin,jMax |
109 |
|
|
DO i=iMin,iMax |
110 |
|
|
fZon(i,j) = afFacS*af(i,j) + dfFacS*df(i,j) |
111 |
|
|
ENDDO |
112 |
|
|
ENDDO |
113 |
|
|
|
114 |
|
|
C-- Meridional flux (fMer is at south face of "salt" cell) |
115 |
|
|
C Advective component of meridional flux |
116 |
|
|
DO j=jMin,jMax |
117 |
|
|
DO i=iMin,iMax |
118 |
|
|
C Advective component of meridional flux |
119 |
|
|
af(i,j) = |
120 |
|
|
& vTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j-1,k,bi,bj))*0.5 _d 0 |
121 |
|
|
ENDDO |
122 |
|
|
ENDDO |
123 |
|
|
C Diffusive component of meridional flux |
124 |
|
|
DO j=jMin,jMax |
125 |
|
|
DO i=iMin,iMax |
126 |
|
|
df(i,j) = |
127 |
cnh |
1.4 |
& -diffKhS*yA(i,j)*_rdyC(i,j,bi,bj) |
128 |
cnh |
1.1 |
& *(salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj)) |
129 |
|
|
ENDDO |
130 |
|
|
ENDDO |
131 |
|
|
C Net meridional flux |
132 |
|
|
DO j=jMin,jMax |
133 |
|
|
DO i=iMin,iMax |
134 |
|
|
fMer(i,j) = afFacS*af(i,j) + dfFacS*df(i,j) |
135 |
|
|
ENDDO |
136 |
|
|
ENDDO |
137 |
|
|
|
138 |
|
|
C-- Vertical flux (fVerS) above |
139 |
|
|
C Note: For K=1 then KM1=1 this gives a dS/dz = 0 upper |
140 |
|
|
C boundary condition. |
141 |
|
|
C Advective component of vertical flux |
142 |
|
|
DO j=jMin,jMax |
143 |
|
|
DO i=iMin,iMax |
144 |
|
|
af(i,j) = |
145 |
|
|
& wTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j,kM1,bi,bj))*0.5 _d 0 |
146 |
|
|
ENDDO |
147 |
|
|
ENDDO |
148 |
|
|
C Diffusive component of vertical flux |
149 |
|
|
DO j=jMin,jMax |
150 |
|
|
DO i=iMin,iMax |
151 |
|
|
df(i,j) = |
152 |
|
|
& -diffKzS*zA(i,j,bi,bj)*rdzC(k) |
153 |
|
|
& *(salt(i,j,kM1,bi,bj)-salt(i,j,k,bi,bj)) |
154 |
|
|
ENDDO |
155 |
|
|
ENDDO |
156 |
|
|
C Net vertical flux |
157 |
|
|
DO j=jMin,jMax |
158 |
|
|
DO i=iMin,iMax |
159 |
|
|
fVerS(i,j,kUp) = (afFacS*af(i,j) + dfFacS*df(i,j))*maskUp(i,j) |
160 |
|
|
ENDDO |
161 |
|
|
ENDDO |
162 |
|
|
|
163 |
|
|
C-- Tendency is minus divergence of the fluxes. |
164 |
|
|
C Note. Tendency terms will only be correct for range |
165 |
|
|
C i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points |
166 |
|
|
C will contain valid floating point numbers but |
167 |
|
|
C they are not algorithmically correct. These points |
168 |
|
|
C are not used. |
169 |
|
|
DO j=jMin,jMax |
170 |
|
|
DO i=iMin,iMax |
171 |
|
|
gS(i,j,k,bi,bj)= |
172 |
cnh |
1.4 |
& -rHFacC(i,j,k,bi,bj)*rdzF(k)*_rdxF(i,j,bi,bj)*_rdyF(i,j,bi,bj) |
173 |
cnh |
1.1 |
& *( |
174 |
|
|
& +( fZon(i+1,j)-fZon(i,j) ) |
175 |
|
|
& +( fMer(i,j+1)-fMer(i,j) ) |
176 |
|
|
& +( fVerS(i,j,kUp)-fVerS(i,j,kDown) ) |
177 |
|
|
& ) |
178 |
|
|
ENDDO |
179 |
|
|
ENDDO |
180 |
|
|
|
181 |
|
|
C-- External haline forcing term(s) |
182 |
|
|
|
183 |
|
|
RETURN |
184 |
|
|
END |