1 |
C $Header: /u/gcmpack/MITgcm/model/src/find_rho.F,v 1.19 2002/09/05 20:49:33 mlosch Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
#define USE_FACTORIZED_POLY |
6 |
|
7 |
CBOP |
8 |
C !ROUTINE: FIND_RHO |
9 |
C !INTERFACE: |
10 |
subroutine FIND_RHO( |
11 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
12 |
I tFld, sFld, |
13 |
O rholoc, |
14 |
I myThid ) |
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 |
C | |
26 |
C *==========================================================* |
27 |
C \ev |
28 |
|
29 |
C !USES: |
30 |
implicit none |
31 |
C == Global variables == |
32 |
#include "SIZE.h" |
33 |
#include "EEPARAMS.h" |
34 |
#include "PARAMS.h" |
35 |
#include "EOS.h" |
36 |
#include "GRID.h" |
37 |
|
38 |
C !INPUT/OUTPUT PARAMETERS: |
39 |
C == Routine arguments == |
40 |
integer bi,bj,iMin,iMax,jMin,jMax |
41 |
integer k ! Level of Theta/Salt slice |
42 |
integer kRef ! Pressure reference level |
43 |
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
44 |
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
45 |
_RL rholoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
46 |
integer myThid |
47 |
|
48 |
C !LOCAL VARIABLES: |
49 |
C == Local variables == |
50 |
integer i,j |
51 |
_RL refTemp,refSalt,sigRef,tP,sP,deltaSig |
52 |
_RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
53 |
_RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
54 |
_RL rhoNum(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
55 |
_RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
56 |
character*(max_len_mbuf) msgbuf |
57 |
CEOP |
58 |
|
59 |
#ifdef ALLOW_AUTODIFF_TAMC |
60 |
DO j=1-OLy,sNy+OLy |
61 |
DO i=1-OLx,sNx+OLx |
62 |
rholoc(i,j) = 0. _d 0 |
63 |
rhoP0(i,j) = 0. _d 0 |
64 |
bulkMod(i,j) = 0. _d 0 |
65 |
ENDDO |
66 |
ENDDO |
67 |
#endif |
68 |
|
69 |
if (equationOfState.eq.'LINEAR') then |
70 |
|
71 |
C ***NOTE*** |
72 |
C In the linear EOS, to make the static stability calculation meaningful |
73 |
C we alway calculate the perturbation with respect to the surface level. |
74 |
C ********** |
75 |
refTemp=tRef(kRef) |
76 |
refSalt=sRef(kRef) |
77 |
|
78 |
do j=jMin,jMax |
79 |
do i=iMin,iMax |
80 |
rholoc(i,j)=rhonil*( |
81 |
& sBeta*(sFld(i,j,k,bi,bj)-refSalt) |
82 |
& -tAlpha*(tFld(i,j,k,bi,bj)-refTemp) ) |
83 |
enddo |
84 |
enddo |
85 |
|
86 |
elseif (equationOfState.eq.'POLY3') then |
87 |
|
88 |
refTemp=eosRefT(kRef) |
89 |
refSalt=eosRefS(kRef) |
90 |
sigRef=eosSig0(kRef) + (1000.-Rhonil) |
91 |
|
92 |
do j=jMin,jMax |
93 |
do i=iMin,iMax |
94 |
tP=tFld(i,j,k,bi,bj)-refTemp |
95 |
sP=sFld(i,j,k,bi,bj)-refSalt |
96 |
#ifdef USE_FACTORIZED_POLY |
97 |
deltaSig= |
98 |
& (( eosC(9,kRef)*sP + eosC(5,kRef) )*sP + eosC(2,kRef) )*sP |
99 |
& + ( ( eosC(6,kRef) |
100 |
& *tP |
101 |
& +eosC(7,kRef)*sP + eosC(3,kRef) |
102 |
& )*tP |
103 |
& +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef) |
104 |
& )*tP |
105 |
#else |
106 |
deltaSig= |
107 |
& eosC(1,kRef)*tP |
108 |
& +eosC(2,kRef) *sP |
109 |
& +eosC(3,kRef)*tP*tP |
110 |
& +eosC(4,kRef)*tP *sP |
111 |
& +eosC(5,kRef) *sP*sP |
112 |
& +eosC(6,kRef)*tP*tP*tP |
113 |
& +eosC(7,kRef)*tP*tP *sP |
114 |
& +eosC(8,kRef)*tP *sP*sP |
115 |
& +eosC(9,kRef) *sP*sP*sP |
116 |
#endif |
117 |
rholoc(i,j)=sigRef+deltaSig |
118 |
enddo |
119 |
enddo |
120 |
|
121 |
elseif ( equationOfState(1:5).eq.'JMD95' |
122 |
& .or. equationOfState.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, kRef, |
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,kRef,bi,bj)*SItoBar/bulkMod(i,j)) |
144 |
& - rhonil |
145 |
|
146 |
end do |
147 |
end do |
148 |
|
149 |
elseif ( equationOfState.eq.'MDJWF' ) then |
150 |
|
151 |
CALL FIND_RHONUM( bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
152 |
& tFld, sFld, rhoNum, myThid ) |
153 |
|
154 |
CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
155 |
& tFld, sFld, rhoDen, myThid ) |
156 |
|
157 |
do j=jMin,jMax |
158 |
do i=iMin,iMax |
159 |
|
160 |
rholoc(i,j) = rhoNum(i,j)*rhoDen(i,j) - rhoNil |
161 |
|
162 |
end do |
163 |
end do |
164 |
elseif( equationOfState .eq. 'IDEALG' ) then |
165 |
C |
166 |
else |
167 |
write(msgbuf,'(3a)') |
168 |
& ' FIND_RHO: equationOfState = "',equationOfState,'"' |
169 |
call print_error( msgbuf, mythid ) |
170 |
stop 'ABNORMAL END: S/R FIND_RHO' |
171 |
endif |
172 |
|
173 |
return |
174 |
end |
175 |
|
176 |
CBOP |
177 |
C !ROUTINE: FIND_RHOP0 |
178 |
C !INTERFACE: |
179 |
subroutine FIND_RHOP0( |
180 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
181 |
I tFld, sFld, |
182 |
O rhoP0, |
183 |
I myThid ) |
184 |
|
185 |
C !DESCRIPTION: \bv |
186 |
C *==========================================================* |
187 |
C | o SUBROUTINE FIND_RHOP0 |
188 |
C | Calculates rho(S,T,0) of a slice |
189 |
C *==========================================================* |
190 |
C | |
191 |
C | k - is the surface level |
192 |
C | |
193 |
C *==========================================================* |
194 |
C \ev |
195 |
|
196 |
C !USES: |
197 |
implicit none |
198 |
C == Global variables == |
199 |
#include "SIZE.h" |
200 |
#include "EEPARAMS.h" |
201 |
#include "PARAMS.h" |
202 |
#include "EOS.h" |
203 |
|
204 |
C !INPUT/OUTPUT PARAMETERS: |
205 |
C == Routine arguments == |
206 |
integer bi,bj,iMin,iMax,jMin,jMax |
207 |
integer k ! Level of surface slice |
208 |
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
209 |
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
210 |
_RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
211 |
_RL rfresh, rsalt |
212 |
integer myThid |
213 |
|
214 |
C !LOCAL VARIABLES: |
215 |
C == Local variables == |
216 |
integer i,j |
217 |
_RL t, t2, t3, t4, s, s3o2 |
218 |
CEOP |
219 |
|
220 |
do j=jMin,jMax |
221 |
do i=iMin,iMax |
222 |
C abbreviations |
223 |
t = tFld(i,j,k,bi,bj) |
224 |
t2 = t*t |
225 |
t3 = t2*t |
226 |
t4 = t3*t |
227 |
|
228 |
s = sFld(i,j,k,bi,bj) |
229 |
if ( s .lt. 0. _d 0 ) then |
230 |
C issue a warning |
231 |
write(*,'(a,i3,a,i3,a,i3,a,i3,a,i3,a,e13.5)') |
232 |
& ' FIND_RHOP0: WARNING, salinity(', |
233 |
& i,',',j,',',k,',',bi,',',bj,') = ', s |
234 |
s = 0. _d 0 |
235 |
end if |
236 |
s3o2 = s*sqrt(s) |
237 |
|
238 |
C density of freshwater at the surface |
239 |
rfresh = |
240 |
& eosJMDCFw(1) |
241 |
& + eosJMDCFw(2)*t |
242 |
& + eosJMDCFw(3)*t2 |
243 |
& + eosJMDCFw(4)*t3 |
244 |
& + eosJMDCFw(5)*t4 |
245 |
& + eosJMDCFw(6)*t4*t |
246 |
C density of sea water at the surface |
247 |
rsalt = |
248 |
& s*( |
249 |
& eosJMDCSw(1) |
250 |
& + eosJMDCSw(2)*t |
251 |
& + eosJMDCSw(3)*t2 |
252 |
& + eosJMDCSw(4)*t3 |
253 |
& + eosJMDCSw(5)*t4 |
254 |
& ) |
255 |
& + s3o2*( |
256 |
& eosJMDCSw(6) |
257 |
& + eosJMDCSw(7)*t |
258 |
& + eosJMDCSw(8)*t2 |
259 |
& ) |
260 |
& + eosJMDCSw(9)*s*s |
261 |
|
262 |
rhoP0(i,j) = rfresh + rsalt |
263 |
|
264 |
end do |
265 |
end do |
266 |
|
267 |
return |
268 |
end |
269 |
|
270 |
C !ROUTINE: FIND_BULKMOD |
271 |
C !INTERFACE: |
272 |
subroutine FIND_BULKMOD( |
273 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
274 |
I tFld, sFld, |
275 |
O bulkMod, |
276 |
I myThid ) |
277 |
|
278 |
C !DESCRIPTION: \bv |
279 |
C *==========================================================* |
280 |
C | o SUBROUTINE FIND_BULKMOD |
281 |
C | Calculates the secant bulk modulus K(S,T,p) of a slice |
282 |
C *==========================================================* |
283 |
C | |
284 |
C | k - is the surface level |
285 |
C | kRef - is the level to use to determine the pressure |
286 |
C | |
287 |
C *==========================================================* |
288 |
C \ev |
289 |
|
290 |
C !USES: |
291 |
implicit none |
292 |
C == Global variables == |
293 |
#include "SIZE.h" |
294 |
#include "EEPARAMS.h" |
295 |
#include "PARAMS.h" |
296 |
#include "EOS.h" |
297 |
|
298 |
C !INPUT/OUTPUT PARAMETERS: |
299 |
C == Routine arguments == |
300 |
integer bi,bj,iMin,iMax,jMin,jMax |
301 |
integer k ! Level of surface slice |
302 |
integer kRef ! Reference pressure level |
303 |
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
304 |
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
305 |
_RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
306 |
_RL bMfresh, bMsalt, bMpres |
307 |
integer myThid |
308 |
|
309 |
C !LOCAL VARIABLES: |
310 |
C == Local variables == |
311 |
integer i,j |
312 |
_RL t, t2, t3, t4, s, s3o2, p, p2 |
313 |
CEOP |
314 |
|
315 |
do j=jMin,jMax |
316 |
do i=iMin,iMax |
317 |
C abbreviations |
318 |
t = tFld(i,j,k,bi,bj) |
319 |
t2 = t*t |
320 |
t3 = t2*t |
321 |
t4 = t3*t |
322 |
|
323 |
s = sFld(i,j,k,bi,bj) |
324 |
if ( s .lt. 0. _d 0 ) then |
325 |
C issue a warning |
326 |
write(*,'(a,i3,a,i3,a,i3,a,i3,a,i3,a,e13.5)') |
327 |
& ' FIND_BULKMOD: WARNING, salinity(', |
328 |
& i,',',j,',',k,',',bi,',',bj,') = ', s |
329 |
s = 0. _d 0 |
330 |
end if |
331 |
s3o2 = s*sqrt(s) |
332 |
C |
333 |
p = pressure(i,j,kRef,bi,bj)*SItoBar |
334 |
p2 = p*p |
335 |
C secant bulk modulus of fresh water at the surface |
336 |
bMfresh = |
337 |
& eosJMDCKFw(1) |
338 |
& + eosJMDCKFw(2)*t |
339 |
& + eosJMDCKFw(3)*t2 |
340 |
& + eosJMDCKFw(4)*t3 |
341 |
& + eosJMDCKFw(5)*t4 |
342 |
C secant bulk modulus of sea water at the surface |
343 |
bMsalt = |
344 |
& s*( eosJMDCKSw(1) |
345 |
& + eosJMDCKSw(2)*t |
346 |
& + eosJMDCKSw(3)*t2 |
347 |
& + eosJMDCKSw(4)*t3 |
348 |
& ) |
349 |
& + s3o2*( eosJMDCKSw(5) |
350 |
& + eosJMDCKSw(6)*t |
351 |
& + eosJMDCKSw(7)*t2 |
352 |
& ) |
353 |
C secant bulk modulus of sea water at pressure p |
354 |
bMpres = |
355 |
& p*( eosJMDCKP(1) |
356 |
& + eosJMDCKP(2)*t |
357 |
& + eosJMDCKP(3)*t2 |
358 |
& + eosJMDCKP(4)*t3 |
359 |
& ) |
360 |
& + p*s*( eosJMDCKP(5) |
361 |
& + eosJMDCKP(6)*t |
362 |
& + eosJMDCKP(7)*t2 |
363 |
& ) |
364 |
& + p*s3o2*eosJMDCKP(8) |
365 |
& + p2*( eosJMDCKP(9) |
366 |
& + eosJMDCKP(10)*t |
367 |
& + eosJMDCKP(11)*t2 |
368 |
& ) |
369 |
& + p2*s*( eosJMDCKP(12) |
370 |
& + eosJMDCKP(13)*t |
371 |
& + eosJMDCKP(14)*t2 |
372 |
& ) |
373 |
|
374 |
bulkMod(i,j) = bMfresh + bMsalt + bMpres |
375 |
|
376 |
end do |
377 |
end do |
378 |
|
379 |
return |
380 |
end |
381 |
|
382 |
CBOP |
383 |
C !ROUTINE: FIND_RHONUM |
384 |
C !INTERFACE: |
385 |
subroutine FIND_RHONUM( |
386 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
387 |
I tFld, sFld, |
388 |
O rhoNum, |
389 |
I myThid ) |
390 |
|
391 |
C !DESCRIPTION: \bv |
392 |
C *==========================================================* |
393 |
C | o SUBROUTINE FIND_RHONUM |
394 |
C | Calculates the numerator of the McDougall et al. |
395 |
C | equation of state |
396 |
C | - the code is more or less a copy of MOM4 |
397 |
C *==========================================================* |
398 |
C | |
399 |
C | k - is the surface level |
400 |
C | kRef - is the level to use to determine the pressure |
401 |
C | |
402 |
C *==========================================================* |
403 |
C \ev |
404 |
|
405 |
C !USES: |
406 |
implicit none |
407 |
C == Global variables == |
408 |
#include "SIZE.h" |
409 |
#include "EEPARAMS.h" |
410 |
#include "PARAMS.h" |
411 |
#include "EOS.h" |
412 |
|
413 |
C !INPUT/OUTPUT PARAMETERS: |
414 |
C == Routine arguments == |
415 |
integer bi,bj,iMin,iMax,jMin,jMax |
416 |
integer k ! Level of surface slice |
417 |
integer kRef ! Reference pressure level |
418 |
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
419 |
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
420 |
_RL rhoNum(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
421 |
integer myThid |
422 |
|
423 |
C !LOCAL VARIABLES: |
424 |
C == Local variables == |
425 |
integer i,j |
426 |
_RL t1, t2, s1, p1 |
427 |
CEOP |
428 |
do j=jMin,jMax |
429 |
do i=iMin,iMax |
430 |
C abbreviations |
431 |
t1 = tFld(i,j,k,bi,bj) |
432 |
t2 = t1*t1 |
433 |
s1 = sFld(i,j,k,bi,bj) |
434 |
|
435 |
p1 = pressure(i,j,kRef,bi,bj)*SItodBar |
436 |
|
437 |
rhoNum(i,j) = eosMDJWFnum(0) |
438 |
& + t1*(eosMDJWFnum(1) |
439 |
& + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) ) |
440 |
& + s1*(eosMDJWFnum(4) |
441 |
& + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1) |
442 |
& + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2 |
443 |
& + eosMDJWFnum(9)*s1 |
444 |
& + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) ) |
445 |
|
446 |
end do |
447 |
end do |
448 |
|
449 |
return |
450 |
end |
451 |
|
452 |
CBOP |
453 |
C !ROUTINE: FIND_RHODEN |
454 |
C !INTERFACE: |
455 |
subroutine FIND_RHODEN( |
456 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
457 |
I tFld, sFld, |
458 |
O rhoDen, |
459 |
I myThid ) |
460 |
|
461 |
C !DESCRIPTION: \bv |
462 |
C *==========================================================* |
463 |
C | o SUBROUTINE FIND_RHODEN |
464 |
C | Calculates the denominator of the McDougall et al. |
465 |
C | equation of state |
466 |
C | - the code is more or less a copy of MOM4 |
467 |
C *==========================================================* |
468 |
C | |
469 |
C | k - is the surface level |
470 |
C | kRef - is the level to use to determine the pressure |
471 |
C | |
472 |
C *==========================================================* |
473 |
C \ev |
474 |
|
475 |
C !USES: |
476 |
implicit none |
477 |
C == Global variables == |
478 |
#include "SIZE.h" |
479 |
#include "EEPARAMS.h" |
480 |
#include "PARAMS.h" |
481 |
#include "EOS.h" |
482 |
|
483 |
C !INPUT/OUTPUT PARAMETERS: |
484 |
C == Routine arguments == |
485 |
integer bi,bj,iMin,iMax,jMin,jMax |
486 |
integer k ! Level of surface slice |
487 |
integer kRef ! Reference pressure level |
488 |
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
489 |
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
490 |
_RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
491 |
integer myThid |
492 |
|
493 |
C !LOCAL VARIABLES: |
494 |
C == Local variables == |
495 |
integer i,j |
496 |
_RL t1, t2, s1, sp5, p1, p1t1 |
497 |
_RL den, epsln |
498 |
parameter ( epsln = 0. _d 0 ) |
499 |
CEOP |
500 |
do j=jMin,jMax |
501 |
do i=iMin,iMax |
502 |
C abbreviations |
503 |
t1 = tFld(i,j,k,bi,bj) |
504 |
t2 = t1*t1 |
505 |
s1 = sFld(i,j,k,bi,bj) |
506 |
if ( s1 .lt. 0. _d 0 ) then |
507 |
C issue a warning |
508 |
write(*,'(a,i3,a,i3,a,i3,a,i3,a,i3,a,e13.5)') |
509 |
& ' FIND_RHODEN: WARNING, salinity(', |
510 |
& i,',',j,',',k,',',bi,',',bj,') = ', s1 |
511 |
s1 = 0. _d 0 |
512 |
end if |
513 |
sp5 = sqrt(s1) |
514 |
|
515 |
p1 = pressure(i,j,kRef,bi,bj)*SItodBar |
516 |
p1t1 = p1*t1 |
517 |
|
518 |
den = eosMDJWFden(0) |
519 |
& + t1*(eosMDJWFden(1) |
520 |
& + t1*(eosMDJWFden(2) |
521 |
& + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) ) |
522 |
& + s1*(eosMDJWFden(5) |
523 |
& + t1*(eosMDJWFden(6) |
524 |
& + eosMDJWFden(7)*t2) |
525 |
& + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) ) |
526 |
& + p1*(eosMDJWFden(10) |
527 |
& + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) ) |
528 |
|
529 |
rhoDen(i,j) = 1.0/(epsln+den) |
530 |
|
531 |
end do |
532 |
end do |
533 |
|
534 |
return |
535 |
end |
536 |
|
537 |
subroutine find_rho_scalar( |
538 |
I tLoc, sLoc, pLoc, |
539 |
O rhoLoc, |
540 |
I myThid ) |
541 |
|
542 |
C !DESCRIPTION: \bv |
543 |
C *==========================================================* |
544 |
C | o SUBROUTINE FIND_RHO_SCALAR |
545 |
C | Calculates [rho(S,T,p)-Rhonil] |
546 |
C *==========================================================* |
547 |
C \ev |
548 |
|
549 |
C !USES: |
550 |
implicit none |
551 |
C == Global variables == |
552 |
#include "SIZE.h" |
553 |
#include "EEPARAMS.h" |
554 |
#include "PARAMS.h" |
555 |
#include "EOS.h" |
556 |
|
557 |
C !INPUT/OUTPUT PARAMETERS: |
558 |
C == Routine arguments == |
559 |
_RL sLoc, tLoc, pLoc |
560 |
_RL rhoLoc |
561 |
integer myThid |
562 |
|
563 |
C !LOCAL VARIABLES: |
564 |
C == Local variables == |
565 |
|
566 |
_RL t1, t2, t3, t4, s1, s3o2, p1, p2, sp5, p1t1 |
567 |
_RL rfresh, rsalt, rhoP0 |
568 |
_RL bMfresh, bMsalt, bMpres, BulkMod |
569 |
_RL rhoNum, rhoDen, den, epsln |
570 |
parameter ( epsln = 0. _d 0 ) |
571 |
|
572 |
character*(max_len_mbuf) msgbuf |
573 |
CEOP |
574 |
|
575 |
rhoLoc = 0. _d 0 |
576 |
rhoP0 = 0. _d 0 |
577 |
bulkMod = 0. _d 0 |
578 |
rfresh = 0. _d 0 |
579 |
rsalt = 0. _d 0 |
580 |
bMfresh = 0. _d 0 |
581 |
bMsalt = 0. _d 0 |
582 |
bMpres = 0. _d 0 |
583 |
rhoNum = 0. _d 0 |
584 |
rhoDen = 0. _d 0 |
585 |
den = 0. _d 0 |
586 |
|
587 |
t1 = tLoc |
588 |
t2 = t1*t1 |
589 |
t3 = t2*t1 |
590 |
t4 = t3*t1 |
591 |
|
592 |
s1 = sLoc |
593 |
if ( s1 .lt. 0. _d 0 ) then |
594 |
C issue a warning |
595 |
write(*,'(a,i3,a,i3,a,i3,a,e13.5)') |
596 |
& ' FIND_RHO_SCALAR: WARNING, salinity = ', s1 |
597 |
s1 = 0. _d 0 |
598 |
end if |
599 |
|
600 |
if (equationOfState.eq.'LINEAR') then |
601 |
|
602 |
rhoLoc = 0. _d 0 |
603 |
|
604 |
elseif (equationOfState.eq.'POLY3') then |
605 |
|
606 |
C this is not correct, there is a field eosSig0 which should be use here |
607 |
C but I do not intent to include the reference level in this routine |
608 |
write(*,'(a)') |
609 |
& ' FIND_RHO_SCALAR: for POLY3, the density is not' |
610 |
write(*,'(a)') |
611 |
& ' computed correctly in this routine' |
612 |
rhoLoc = 0. _d 0 |
613 |
|
614 |
elseif ( equationOfState(1:5).eq.'JMD95' |
615 |
& .or. equationOfState.eq.'UNESCO' ) then |
616 |
C nonlinear equation of state in pressure coordinates |
617 |
|
618 |
s3o2 = s1*sqrt(s1) |
619 |
|
620 |
p1 = pLoc*SItoBar |
621 |
p2 = p1*p1 |
622 |
|
623 |
C density of freshwater at the surface |
624 |
rfresh = |
625 |
& eosJMDCFw(1) |
626 |
& + eosJMDCFw(2)*t1 |
627 |
& + eosJMDCFw(3)*t2 |
628 |
& + eosJMDCFw(4)*t3 |
629 |
& + eosJMDCFw(5)*t4 |
630 |
& + eosJMDCFw(6)*t4*t1 |
631 |
C density of sea water at the surface |
632 |
rsalt = |
633 |
& s1*( |
634 |
& eosJMDCSw(1) |
635 |
& + eosJMDCSw(2)*t1 |
636 |
& + eosJMDCSw(3)*t2 |
637 |
& + eosJMDCSw(4)*t3 |
638 |
& + eosJMDCSw(5)*t4 |
639 |
& ) |
640 |
& + s3o2*( |
641 |
& eosJMDCSw(6) |
642 |
& + eosJMDCSw(7)*t1 |
643 |
& + eosJMDCSw(8)*t2 |
644 |
& ) |
645 |
& + eosJMDCSw(9)*s1*s1 |
646 |
|
647 |
rhoP0 = rfresh + rsalt |
648 |
|
649 |
C secant bulk modulus of fresh water at the surface |
650 |
bMfresh = |
651 |
& eosJMDCKFw(1) |
652 |
& + eosJMDCKFw(2)*t1 |
653 |
& + eosJMDCKFw(3)*t2 |
654 |
& + eosJMDCKFw(4)*t3 |
655 |
& + eosJMDCKFw(5)*t4 |
656 |
C secant bulk modulus of sea water at the surface |
657 |
bMsalt = |
658 |
& s1*( eosJMDCKSw(1) |
659 |
& + eosJMDCKSw(2)*t1 |
660 |
& + eosJMDCKSw(3)*t2 |
661 |
& + eosJMDCKSw(4)*t3 |
662 |
& ) |
663 |
& + s3o2*( eosJMDCKSw(5) |
664 |
& + eosJMDCKSw(6)*t1 |
665 |
& + eosJMDCKSw(7)*t2 |
666 |
& ) |
667 |
C secant bulk modulus of sea water at pressure p |
668 |
bMpres = |
669 |
& p1*( eosJMDCKP(1) |
670 |
& + eosJMDCKP(2)*t1 |
671 |
& + eosJMDCKP(3)*t2 |
672 |
& + eosJMDCKP(4)*t3 |
673 |
& ) |
674 |
& + p1*s1*( eosJMDCKP(5) |
675 |
& + eosJMDCKP(6)*t1 |
676 |
& + eosJMDCKP(7)*t2 |
677 |
& ) |
678 |
& + p1*s3o2*eosJMDCKP(8) |
679 |
& + p2*( eosJMDCKP(9) |
680 |
& + eosJMDCKP(10)*t1 |
681 |
& + eosJMDCKP(11)*t2 |
682 |
& ) |
683 |
& + p2*s1*( eosJMDCKP(12) |
684 |
& + eosJMDCKP(13)*t1 |
685 |
& + eosJMDCKP(14)*t2 |
686 |
& ) |
687 |
|
688 |
bulkMod = bMfresh + bMsalt + bMpres |
689 |
|
690 |
C density of sea water at pressure p |
691 |
rhoLoc = rhoP0/(1. _d 0 - p1/bulkMod) - rhonil |
692 |
|
693 |
elseif ( equationOfState.eq.'MDJWF' ) then |
694 |
|
695 |
sp5 = sqrt(s1) |
696 |
|
697 |
p1 = pLoc*SItodBar |
698 |
p1t1 = p1*t1 |
699 |
|
700 |
rhoNum = eosMDJWFnum(0) |
701 |
& + t1*(eosMDJWFnum(1) |
702 |
& + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) ) |
703 |
& + s1*(eosMDJWFnum(4) |
704 |
& + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1) |
705 |
& + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2 |
706 |
& + eosMDJWFnum(9)*s1 |
707 |
& + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) ) |
708 |
|
709 |
|
710 |
den = eosMDJWFden(0) |
711 |
& + t1*(eosMDJWFden(1) |
712 |
& + t1*(eosMDJWFden(2) |
713 |
& + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) ) |
714 |
& + s1*(eosMDJWFden(5) |
715 |
& + t1*(eosMDJWFden(6) |
716 |
& + eosMDJWFden(7)*t2) |
717 |
& + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) ) |
718 |
& + p1*(eosMDJWFden(10) |
719 |
& + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) ) |
720 |
|
721 |
rhoDen = 1.0/(epsln+den) |
722 |
|
723 |
rhoLoc = rhoNum*rhoDen - rhoNil |
724 |
|
725 |
elseif( equationOfState .eq. 'IDEALG' ) then |
726 |
C |
727 |
else |
728 |
write(msgbuf,'(3a)') |
729 |
& ' FIND_RHO_SCALAR : equationOfState = "', |
730 |
& equationOfState,'"' |
731 |
call print_error( msgbuf, mythid ) |
732 |
stop 'ABNORMAL END: S/R FIND_RHO_SCALAR' |
733 |
endif |
734 |
|
735 |
return |
736 |
end |