1 |
mlosch |
1.16 |
C $Header: /u/gcmpack/MITgcm/model/src/find_rho.F,v 1.15 2001/09/26 18:09:15 cnh Exp $ |
2 |
cnh |
1.15 |
C $Name: $ |
3 |
cnh |
1.1 |
|
4 |
cnh |
1.9 |
#include "CPP_OPTIONS.h" |
5 |
adcroft |
1.5 |
#define USE_FACTORIZED_POLY |
6 |
cnh |
1.1 |
|
7 |
cnh |
1.15 |
CBOP |
8 |
|
|
C !ROUTINE: FIND_RHO |
9 |
|
|
C !INTERFACE: |
10 |
adcroft |
1.4 |
subroutine FIND_RHO( |
11 |
|
|
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, eqn, |
12 |
adcroft |
1.13 |
I tFld, sFld, |
13 |
adcroft |
1.4 |
O rholoc, |
14 |
|
|
I myThid ) |
15 |
cnh |
1.15 |
|
16 |
|
|
C !DESCRIPTION: \bv |
17 |
|
|
C *==========================================================* |
18 |
|
|
C | o SUBROUTINE FIND_RHO |
19 |
|
|
C | Calculates [rho(S,T,z)-Rhonil] of a slice |
20 |
|
|
C *==========================================================* |
21 |
|
|
C | |
22 |
|
|
C | k - is the Theta/Salt level |
23 |
|
|
C | kRef - determines pressure reference level |
24 |
|
|
C | (not used in 'LINEAR' mode) |
25 |
mlosch |
1.16 |
C | eqn - determines the eqn. of state: |
26 |
|
|
C | 'LINEAR', 'POLY3', 'UNESCO', 'JMD95Z', 'JMD95P' |
27 |
cnh |
1.15 |
C | |
28 |
|
|
C *==========================================================* |
29 |
|
|
C \ev |
30 |
|
|
|
31 |
|
|
C !USES: |
32 |
cnh |
1.1 |
implicit none |
33 |
heimbach |
1.12 |
C == Global variables == |
34 |
cnh |
1.1 |
#include "SIZE.h" |
35 |
cnh |
1.7 |
#include "EEPARAMS.h" |
36 |
cnh |
1.1 |
#include "PARAMS.h" |
37 |
mlosch |
1.16 |
#include "EOS.h" |
38 |
|
|
#include "GRID.h" |
39 |
heimbach |
1.12 |
|
40 |
cnh |
1.15 |
C !INPUT/OUTPUT PARAMETERS: |
41 |
heimbach |
1.12 |
C == Routine arguments == |
42 |
adcroft |
1.4 |
integer bi,bj,iMin,iMax,jMin,jMax |
43 |
|
|
integer k ! Level of Theta/Salt slice |
44 |
|
|
integer kRef ! Pressure reference level |
45 |
|
|
character*(*) eqn |
46 |
adcroft |
1.13 |
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
47 |
|
|
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
48 |
adcroft |
1.4 |
_RL rholoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
49 |
|
|
integer myThid |
50 |
heimbach |
1.12 |
|
51 |
cnh |
1.15 |
C !LOCAL VARIABLES: |
52 |
heimbach |
1.12 |
C == Local variables == |
53 |
adcroft |
1.4 |
integer i,j |
54 |
|
|
_RL refTemp,refSalt,sigRef,tP,sP,deltaSig |
55 |
mlosch |
1.16 |
_RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
56 |
|
|
_RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
57 |
adcroft |
1.10 |
character*(max_len_mbuf) msgbuf |
58 |
cnh |
1.15 |
CEOP |
59 |
heimbach |
1.11 |
|
60 |
|
|
#ifdef ALLOW_AUTODIFF_TAMC |
61 |
|
|
DO j=1-OLy,sNy+OLy |
62 |
|
|
DO i=1-OLx,sNx+OLx |
63 |
mlosch |
1.16 |
rholoc(i,j) = 0. _d 0 |
64 |
|
|
rhoP0(i,j) = 0. _d 0 |
65 |
|
|
bulkMod(i,j) = 0. _d 0 |
66 |
heimbach |
1.11 |
ENDDO |
67 |
|
|
ENDDO |
68 |
|
|
#endif |
69 |
cnh |
1.1 |
|
70 |
adcroft |
1.4 |
if (eqn.eq.'LINEAR') then |
71 |
|
|
|
72 |
|
|
C ***NOTE*** |
73 |
|
|
C In the linear EOS, to make the static stability calculation meaningful |
74 |
|
|
C we alway calculate the perturbation with respect to the surface level. |
75 |
|
|
C ********** |
76 |
adcroft |
1.6 |
refTemp=tRef(kRef) |
77 |
|
|
refSalt=sRef(kRef) |
78 |
adcroft |
1.4 |
|
79 |
cnh |
1.1 |
do j=jMin,jMax |
80 |
|
|
do i=iMin,iMax |
81 |
adcroft |
1.4 |
rholoc(i,j)=rhonil*( |
82 |
adcroft |
1.13 |
& sBeta*(sFld(i,j,k,bi,bj)-refSalt) |
83 |
|
|
& -tAlpha*(tFld(i,j,k,bi,bj)-refTemp) ) |
84 |
cnh |
1.1 |
enddo |
85 |
|
|
enddo |
86 |
cnh |
1.8 |
|
87 |
adcroft |
1.4 |
elseif (eqn.eq.'POLY3') then |
88 |
|
|
|
89 |
|
|
refTemp=eosRefT(kRef) |
90 |
|
|
refSalt=eosRefS(kRef) |
91 |
adcroft |
1.5 |
sigRef=eosSig0(kRef) + (1000.-Rhonil) |
92 |
adcroft |
1.4 |
|
93 |
|
|
do j=jMin,jMax |
94 |
|
|
do i=iMin,iMax |
95 |
adcroft |
1.13 |
tP=tFld(i,j,k,bi,bj)-refTemp |
96 |
|
|
sP=sFld(i,j,k,bi,bj)-refSalt |
97 |
adcroft |
1.5 |
#ifdef USE_FACTORIZED_POLY |
98 |
adcroft |
1.4 |
deltaSig= |
99 |
adcroft |
1.5 |
& (( eosC(9,kRef)*sP + eosC(5,kRef) )*sP + eosC(2,kRef) )*sP |
100 |
|
|
& + ( ( eosC(6,kRef) |
101 |
|
|
& *tP |
102 |
|
|
& +eosC(7,kRef)*sP + eosC(3,kRef) |
103 |
|
|
& )*tP |
104 |
|
|
& +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef) |
105 |
|
|
& )*tP |
106 |
|
|
#else |
107 |
|
|
deltaSig= |
108 |
|
|
& eosC(1,kRef)*tP |
109 |
|
|
& +eosC(2,kRef) *sP |
110 |
|
|
& +eosC(3,kRef)*tP*tP |
111 |
|
|
& +eosC(4,kRef)*tP *sP |
112 |
|
|
& +eosC(5,kRef) *sP*sP |
113 |
|
|
& +eosC(6,kRef)*tP*tP*tP |
114 |
|
|
& +eosC(7,kRef)*tP*tP *sP |
115 |
|
|
& +eosC(8,kRef)*tP *sP*sP |
116 |
|
|
& +eosC(9,kRef) *sP*sP*sP |
117 |
|
|
#endif |
118 |
|
|
rholoc(i,j)=sigRef+deltaSig |
119 |
adcroft |
1.4 |
enddo |
120 |
|
|
enddo |
121 |
|
|
|
122 |
mlosch |
1.16 |
elseif ( eqn(1:5).eq.'JMD95' .or. eqn.eq.'UNESCO' ) then |
123 |
|
|
C nonlinear equation of state in pressure coordinates |
124 |
|
|
|
125 |
|
|
CALL FIND_RHOP0( |
126 |
|
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
127 |
|
|
I tFld, sFld, |
128 |
|
|
O rhoP0, |
129 |
|
|
I myThid ) |
130 |
|
|
|
131 |
|
|
CALL FIND_BULKMOD( |
132 |
|
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
133 |
|
|
I tFld, sFld, |
134 |
|
|
O bulkMod, |
135 |
|
|
I myThid ) |
136 |
|
|
|
137 |
|
|
do j=jMin,jMax |
138 |
|
|
do i=iMin,iMax |
139 |
|
|
|
140 |
|
|
C density of sea water at pressure p |
141 |
|
|
rholoc(i,j) = rhoP0(i,j) |
142 |
|
|
& /(1. _d 0 - |
143 |
|
|
& pressure(i,j,k,bi,bj)*SItoBar/bulkMod(i,j)) |
144 |
|
|
& - rhonil |
145 |
|
|
|
146 |
|
|
end do |
147 |
|
|
end do |
148 |
|
|
|
149 |
adcroft |
1.4 |
else |
150 |
adcroft |
1.13 |
write(msgbuf,'(3a)') ' FIND_RHO: eqn = "',eqn,'"' |
151 |
adcroft |
1.10 |
call print_error( msgbuf, mythid ) |
152 |
|
|
stop 'ABNORMAL END: S/R FIND_RHO' |
153 |
adcroft |
1.4 |
endif |
154 |
cnh |
1.1 |
|
155 |
|
|
return |
156 |
|
|
end |
157 |
mlosch |
1.16 |
|
158 |
|
|
CBOP |
159 |
|
|
C !ROUTINE: FIND_RHOP0 |
160 |
|
|
C !INTERFACE: |
161 |
|
|
subroutine FIND_RHOP0( |
162 |
|
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
163 |
|
|
I tFld, sFld, |
164 |
|
|
O rhoP0, |
165 |
|
|
I myThid ) |
166 |
|
|
|
167 |
|
|
C !DESCRIPTION: \bv |
168 |
|
|
C *==========================================================* |
169 |
|
|
C | o SUBROUTINE FIND_RHOP0 |
170 |
|
|
C | Calculates rho(S,T,0) of a slice |
171 |
|
|
C *==========================================================* |
172 |
|
|
C | |
173 |
|
|
C | k - is the surface level |
174 |
|
|
C | |
175 |
|
|
C *==========================================================* |
176 |
|
|
C \ev |
177 |
|
|
|
178 |
|
|
C !USES: |
179 |
|
|
implicit none |
180 |
|
|
C == Global variables == |
181 |
|
|
#include "SIZE.h" |
182 |
|
|
#include "EEPARAMS.h" |
183 |
|
|
#include "PARAMS.h" |
184 |
|
|
#include "EOS.h" |
185 |
|
|
|
186 |
|
|
C !INPUT/OUTPUT PARAMETERS: |
187 |
|
|
C == Routine arguments == |
188 |
|
|
integer bi,bj,iMin,iMax,jMin,jMax |
189 |
|
|
integer k ! Level of surface slice |
190 |
|
|
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
191 |
|
|
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
192 |
|
|
_RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
193 |
|
|
_RL rfresh, rsalt |
194 |
|
|
integer myThid |
195 |
|
|
|
196 |
|
|
C !LOCAL VARIABLES: |
197 |
|
|
C == Local variables == |
198 |
|
|
integer i,j |
199 |
|
|
_RL t, t2, t3, t4, s, s3o2 |
200 |
|
|
CEOP |
201 |
|
|
|
202 |
|
|
do j=jMin,jMax |
203 |
|
|
do i=iMin,iMax |
204 |
|
|
C abbreviations |
205 |
|
|
t = tFld(i,j,k,bi,bj) |
206 |
|
|
t2 = t*t |
207 |
|
|
t3 = t2*t |
208 |
|
|
t4 = t3*t |
209 |
|
|
|
210 |
|
|
if ( sFld(i,j,k,bi,bj) .lt. 0. _d 0 ) then |
211 |
|
|
C issue a warning |
212 |
|
|
write(*,'(a,i3,a,i3,a,i3,a,e13.5)') |
213 |
|
|
& ' FIND_RHOP0: WARNING, s(',i,',',j,',',k,') = ', s |
214 |
|
|
s = 0. _d 0 |
215 |
|
|
else |
216 |
|
|
s = sFld(i,j,k,bi,bj) |
217 |
|
|
end if |
218 |
|
|
s3o2 = s*sqrt(s) |
219 |
|
|
|
220 |
|
|
C density of freshwater at the surface |
221 |
|
|
rfresh = |
222 |
|
|
& eosJMDCFw(1) |
223 |
|
|
& + eosJMDCFw(2)*t |
224 |
|
|
& + eosJMDCFw(3)*t2 |
225 |
|
|
& + eosJMDCFw(4)*t3 |
226 |
|
|
& + eosJMDCFw(5)*t4 |
227 |
|
|
& + eosJMDCFw(6)*t4*t |
228 |
|
|
C density of sea water at the surface |
229 |
|
|
rsalt = |
230 |
|
|
& s*( |
231 |
|
|
& eosJMDCSw(1) |
232 |
|
|
& + eosJMDCSw(2)*t |
233 |
|
|
& + eosJMDCSw(3)*t2 |
234 |
|
|
& + eosJMDCSw(4)*t3 |
235 |
|
|
& + eosJMDCSw(5)*t4 |
236 |
|
|
& ) |
237 |
|
|
& + s3o2*( |
238 |
|
|
& eosJMDCSw(6) |
239 |
|
|
& + eosJMDCSw(7)*t |
240 |
|
|
& + eosJMDCSw(8)*t2 |
241 |
|
|
& ) |
242 |
|
|
& + eosJMDCSw(9)*s*s |
243 |
|
|
|
244 |
|
|
rhoP0(i,j) = rfresh + rsalt |
245 |
|
|
|
246 |
|
|
end do |
247 |
|
|
end do |
248 |
|
|
|
249 |
|
|
return |
250 |
|
|
end |
251 |
|
|
|
252 |
|
|
C !ROUTINE: FIND_BULKMOD |
253 |
|
|
C !INTERFACE: |
254 |
|
|
subroutine FIND_BULKMOD( |
255 |
|
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
256 |
|
|
I tFld, sFld, |
257 |
|
|
O bulkMod, |
258 |
|
|
I myThid ) |
259 |
|
|
|
260 |
|
|
C !DESCRIPTION: \bv |
261 |
|
|
C *==========================================================* |
262 |
|
|
C | o SUBROUTINE FIND_BULKMOD |
263 |
|
|
C | Calculates the secant bulk modulus K(S,T,p) of a slice |
264 |
|
|
C *==========================================================* |
265 |
|
|
C | |
266 |
|
|
C | k - is the surface level |
267 |
|
|
C | |
268 |
|
|
C *==========================================================* |
269 |
|
|
C \ev |
270 |
|
|
|
271 |
|
|
C !USES: |
272 |
|
|
implicit none |
273 |
|
|
C == Global variables == |
274 |
|
|
#include "SIZE.h" |
275 |
|
|
#include "EEPARAMS.h" |
276 |
|
|
#include "PARAMS.h" |
277 |
|
|
#include "EOS.h" |
278 |
|
|
|
279 |
|
|
C !INPUT/OUTPUT PARAMETERS: |
280 |
|
|
C == Routine arguments == |
281 |
|
|
integer bi,bj,iMin,iMax,jMin,jMax |
282 |
|
|
integer k ! Level of surface slice |
283 |
|
|
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
284 |
|
|
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
285 |
|
|
_RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
286 |
|
|
_RL bMfresh, bMsalt, bMpres |
287 |
|
|
integer myThid |
288 |
|
|
|
289 |
|
|
C !LOCAL VARIABLES: |
290 |
|
|
C == Local variables == |
291 |
|
|
integer i,j |
292 |
|
|
_RL t, t2, t3, t4, s, s3o2, p, p2 |
293 |
|
|
CEOP |
294 |
|
|
|
295 |
|
|
do j=jMin,jMax |
296 |
|
|
do i=iMin,iMax |
297 |
|
|
C abbreviations |
298 |
|
|
t = tFld(i,j,k,bi,bj) |
299 |
|
|
t2 = t*t |
300 |
|
|
t3 = t2*t |
301 |
|
|
t4 = t3*t |
302 |
|
|
|
303 |
|
|
if ( sFld(i,j,k,bi,bj) .lt. 0. _d 0 ) then |
304 |
|
|
C issue a warning |
305 |
|
|
write(*,'(a,i3,a,i3,a,i3,a,e13.5)') |
306 |
|
|
& ' FIND_BULKMOD: WARNING, s(',i,',',j,',',k,') = ', s |
307 |
|
|
s = 0. _d 0 |
308 |
|
|
else |
309 |
|
|
s = sFld(i,j,k,bi,bj) |
310 |
|
|
end if |
311 |
|
|
s3o2 = s*sqrt(s) |
312 |
|
|
C |
313 |
|
|
p = pressure(i,j,k,bi,bj)*SItoBar |
314 |
|
|
p2 = p*p |
315 |
|
|
C secant bulk modulus of fresh water at the surface |
316 |
|
|
bMfresh = |
317 |
|
|
& eosJMDCKFw(1) |
318 |
|
|
& + eosJMDCKFw(2)*t |
319 |
|
|
& + eosJMDCKFw(3)*t2 |
320 |
|
|
& + eosJMDCKFw(4)*t3 |
321 |
|
|
& + eosJMDCKFw(5)*t4 |
322 |
|
|
C secant bulk modulus of sea water at the surface |
323 |
|
|
bMsalt = |
324 |
|
|
& s*( eosJMDCKSw(1) |
325 |
|
|
& + eosJMDCKSw(2)*t |
326 |
|
|
& + eosJMDCKSw(3)*t2 |
327 |
|
|
& + eosJMDCKSw(4)*t3 |
328 |
|
|
& ) |
329 |
|
|
& + s3o2*( eosJMDCKSw(5) |
330 |
|
|
& + eosJMDCKSw(6)*t |
331 |
|
|
& + eosJMDCKSw(7)*t2 |
332 |
|
|
& ) |
333 |
|
|
C secant bulk modulus of sea water at pressure p |
334 |
|
|
bMpres = |
335 |
|
|
& p*( eosJMDCKP(1) |
336 |
|
|
& + eosJMDCKP(2)*t |
337 |
|
|
& + eosJMDCKP(3)*t2 |
338 |
|
|
& + eosJMDCKP(4)*t3 |
339 |
|
|
& ) |
340 |
|
|
& + p*s*( eosJMDCKP(5) |
341 |
|
|
& + eosJMDCKP(6)*t |
342 |
|
|
& + eosJMDCKP(7)*t2 |
343 |
|
|
& ) |
344 |
|
|
& + p*s3o2*eosJMDCKP(8) |
345 |
|
|
& + p2*( eosJMDCKP(9) |
346 |
|
|
& + eosJMDCKP(10)*t |
347 |
|
|
& + eosJMDCKP(11)*t2 |
348 |
|
|
& ) |
349 |
|
|
& + p2*s*( eosJMDCKP(12) |
350 |
|
|
& + eosJMDCKP(13)*t |
351 |
|
|
& + eosJMDCKP(14)*t2 |
352 |
|
|
& ) |
353 |
|
|
|
354 |
|
|
bulkMod(i,j) = bMfresh + bMsalt + bMpres |
355 |
|
|
|
356 |
|
|
end do |
357 |
|
|
end do |
358 |
|
|
|
359 |
|
|
return |
360 |
|
|
end |
361 |
|
|
|