1 |
mlosch |
1.1 |
C$Header: $ |
2 |
|
|
C$Name: $ |
3 |
|
|
|
4 |
|
|
#include "CPP_OPTIONS.h" |
5 |
|
|
C |
6 |
|
|
C This file contains routines that compute quantities related to |
7 |
|
|
C seawater: |
8 |
|
|
C find_rho_scalar: in-situ density for individual points |
9 |
|
|
C sw_ptmp: function to compute potential temperature |
10 |
|
|
C sw_adtg: function to compute adiabatic tmperature gradient |
11 |
|
|
C used by sw_ptmp |
12 |
|
|
C |
13 |
|
|
SUBROUTINE FIND_RHO_SCALAR( |
14 |
|
|
I tLoc, sLoc, pLoc, |
15 |
|
|
O rhoLoc, |
16 |
|
|
I myThid ) |
17 |
|
|
|
18 |
|
|
C !DESCRIPTION: \bv |
19 |
|
|
C *==========================================================* |
20 |
|
|
C | o SUBROUTINE FIND_RHO_SCALAR |
21 |
|
|
C | Calculates [rho(S,T,p)-rhoConst] |
22 |
|
|
C *==========================================================* |
23 |
|
|
C \ev |
24 |
|
|
|
25 |
|
|
C !USES: |
26 |
|
|
IMPLICIT NONE |
27 |
|
|
C == Global variables == |
28 |
|
|
#include "SIZE.h" |
29 |
|
|
#include "EEPARAMS.h" |
30 |
|
|
#include "PARAMS.h" |
31 |
|
|
#include "EOS.h" |
32 |
|
|
|
33 |
|
|
C !INPUT/OUTPUT PARAMETERS: |
34 |
|
|
C == Routine arguments == |
35 |
|
|
_RL sLoc, tLoc, pLoc |
36 |
|
|
_RL rhoLoc |
37 |
|
|
INTEGER myThid |
38 |
|
|
|
39 |
|
|
C !LOCAL VARIABLES: |
40 |
|
|
C == Local variables == |
41 |
|
|
|
42 |
|
|
_RL t1, t2, t3, t4, s1, s3o2, p1, p2, sp5, p1t1 |
43 |
|
|
_RL rfresh, rsalt, rhoP0 |
44 |
|
|
_RL bMfresh, bMsalt, bMpres, BulkMod |
45 |
|
|
_RL rhoNum, rhoDen, den, epsln |
46 |
|
|
parameter ( epsln = 0. _d 0 ) |
47 |
|
|
|
48 |
|
|
CHARACTER*(MAX_LEN_MBUF) msgBuf |
49 |
|
|
CEOP |
50 |
|
|
|
51 |
|
|
rhoLoc = 0. _d 0 |
52 |
|
|
rhoP0 = 0. _d 0 |
53 |
|
|
bulkMod = 0. _d 0 |
54 |
|
|
rfresh = 0. _d 0 |
55 |
|
|
rsalt = 0. _d 0 |
56 |
|
|
bMfresh = 0. _d 0 |
57 |
|
|
bMsalt = 0. _d 0 |
58 |
|
|
bMpres = 0. _d 0 |
59 |
|
|
rhoNum = 0. _d 0 |
60 |
|
|
rhoDen = 0. _d 0 |
61 |
|
|
den = 0. _d 0 |
62 |
|
|
|
63 |
|
|
t1 = tLoc |
64 |
|
|
t2 = t1*t1 |
65 |
|
|
t3 = t2*t1 |
66 |
|
|
t4 = t3*t1 |
67 |
|
|
|
68 |
|
|
s1 = sLoc |
69 |
|
|
IF ( s1 .LT. 0. _d 0 ) THEN |
70 |
|
|
C issue a warning |
71 |
|
|
WRITE(msgBuf,'(A,E13.5)') |
72 |
|
|
& ' FIND_RHO_SCALAR: WARNING, salinity = ', s1 |
73 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
74 |
|
|
& SQUEEZE_RIGHT , myThid ) |
75 |
|
|
s1 = 0. _d 0 |
76 |
|
|
ENDIF |
77 |
|
|
|
78 |
|
|
IF (equationOfState.EQ.'LINEAR') THEN |
79 |
|
|
|
80 |
|
|
rholoc = rhoNil*( |
81 |
|
|
& sBeta *(sLoc-sRef(1)) |
82 |
|
|
& -tAlpha*(tLoc-tRef(1)) |
83 |
|
|
& ) + (rhoNil-rhoConst) |
84 |
|
|
c rhoLoc = 0. _d 0 |
85 |
|
|
|
86 |
|
|
ELSEIF (equationOfState.EQ.'POLY3') THEN |
87 |
|
|
|
88 |
|
|
C this is not correct, there is a field eosSig0 which should be use here |
89 |
|
|
C but I DO not intent to include the reference level in this routine |
90 |
|
|
WRITE(msgBuf,'(A)') |
91 |
|
|
& ' FIND_RHO_SCALAR: for POLY3, the density is not' |
92 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
93 |
|
|
& SQUEEZE_RIGHT , myThid ) |
94 |
|
|
WRITE(msgBuf,'(A)') |
95 |
|
|
& ' computed correctly in this routine' |
96 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
97 |
|
|
& SQUEEZE_RIGHT , myThid ) |
98 |
|
|
rhoLoc = 0. _d 0 |
99 |
|
|
|
100 |
|
|
ELSEIF ( equationOfState(1:5).EQ.'JMD95' |
101 |
|
|
& .OR. equationOfState.EQ.'UNESCO' ) THEN |
102 |
|
|
C nonlinear equation of state in pressure coordinates |
103 |
|
|
|
104 |
|
|
s3o2 = s1*SQRT(s1) |
105 |
|
|
|
106 |
|
|
p1 = pLoc*SItoBar |
107 |
|
|
p2 = p1*p1 |
108 |
|
|
|
109 |
|
|
C density of freshwater at the surface |
110 |
|
|
rfresh = |
111 |
|
|
& eosJMDCFw(1) |
112 |
|
|
& + eosJMDCFw(2)*t1 |
113 |
|
|
& + eosJMDCFw(3)*t2 |
114 |
|
|
& + eosJMDCFw(4)*t3 |
115 |
|
|
& + eosJMDCFw(5)*t4 |
116 |
|
|
& + eosJMDCFw(6)*t4*t1 |
117 |
|
|
C density of sea water at the surface |
118 |
|
|
rsalt = |
119 |
|
|
& s1*( |
120 |
|
|
& eosJMDCSw(1) |
121 |
|
|
& + eosJMDCSw(2)*t1 |
122 |
|
|
& + eosJMDCSw(3)*t2 |
123 |
|
|
& + eosJMDCSw(4)*t3 |
124 |
|
|
& + eosJMDCSw(5)*t4 |
125 |
|
|
& ) |
126 |
|
|
& + s3o2*( |
127 |
|
|
& eosJMDCSw(6) |
128 |
|
|
& + eosJMDCSw(7)*t1 |
129 |
|
|
& + eosJMDCSw(8)*t2 |
130 |
|
|
& ) |
131 |
|
|
& + eosJMDCSw(9)*s1*s1 |
132 |
|
|
|
133 |
|
|
rhoP0 = rfresh + rsalt |
134 |
|
|
|
135 |
|
|
C secant bulk modulus of fresh water at the surface |
136 |
|
|
bMfresh = |
137 |
|
|
& eosJMDCKFw(1) |
138 |
|
|
& + eosJMDCKFw(2)*t1 |
139 |
|
|
& + eosJMDCKFw(3)*t2 |
140 |
|
|
& + eosJMDCKFw(4)*t3 |
141 |
|
|
& + eosJMDCKFw(5)*t4 |
142 |
|
|
C secant bulk modulus of sea water at the surface |
143 |
|
|
bMsalt = |
144 |
|
|
& s1*( eosJMDCKSw(1) |
145 |
|
|
& + eosJMDCKSw(2)*t1 |
146 |
|
|
& + eosJMDCKSw(3)*t2 |
147 |
|
|
& + eosJMDCKSw(4)*t3 |
148 |
|
|
& ) |
149 |
|
|
& + s3o2*( eosJMDCKSw(5) |
150 |
|
|
& + eosJMDCKSw(6)*t1 |
151 |
|
|
& + eosJMDCKSw(7)*t2 |
152 |
|
|
& ) |
153 |
|
|
C secant bulk modulus of sea water at pressure p |
154 |
|
|
bMpres = |
155 |
|
|
& p1*( eosJMDCKP(1) |
156 |
|
|
& + eosJMDCKP(2)*t1 |
157 |
|
|
& + eosJMDCKP(3)*t2 |
158 |
|
|
& + eosJMDCKP(4)*t3 |
159 |
|
|
& ) |
160 |
|
|
& + p1*s1*( eosJMDCKP(5) |
161 |
|
|
& + eosJMDCKP(6)*t1 |
162 |
|
|
& + eosJMDCKP(7)*t2 |
163 |
|
|
& ) |
164 |
|
|
& + p1*s3o2*eosJMDCKP(8) |
165 |
|
|
& + p2*( eosJMDCKP(9) |
166 |
|
|
& + eosJMDCKP(10)*t1 |
167 |
|
|
& + eosJMDCKP(11)*t2 |
168 |
|
|
& ) |
169 |
|
|
& + p2*s1*( eosJMDCKP(12) |
170 |
|
|
& + eosJMDCKP(13)*t1 |
171 |
|
|
& + eosJMDCKP(14)*t2 |
172 |
|
|
& ) |
173 |
|
|
|
174 |
|
|
bulkMod = bMfresh + bMsalt + bMpres |
175 |
|
|
|
176 |
|
|
C density of sea water at pressure p |
177 |
|
|
rhoLoc = rhoP0/(1. _d 0 - p1/bulkMod) - rhoConst |
178 |
|
|
|
179 |
|
|
ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN |
180 |
|
|
|
181 |
|
|
sp5 = SQRT(s1) |
182 |
|
|
|
183 |
|
|
p1 = pLoc*SItodBar |
184 |
|
|
p1t1 = p1*t1 |
185 |
|
|
|
186 |
|
|
rhoNum = eosMDJWFnum(0) |
187 |
|
|
& + t1*(eosMDJWFnum(1) |
188 |
|
|
& + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) ) |
189 |
|
|
& + s1*(eosMDJWFnum(4) |
190 |
|
|
& + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1) |
191 |
|
|
& + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2 |
192 |
|
|
& + eosMDJWFnum(9)*s1 |
193 |
|
|
& + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) ) |
194 |
|
|
|
195 |
|
|
|
196 |
|
|
den = eosMDJWFden(0) |
197 |
|
|
& + t1*(eosMDJWFden(1) |
198 |
|
|
& + t1*(eosMDJWFden(2) |
199 |
|
|
& + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) ) |
200 |
|
|
& + s1*(eosMDJWFden(5) |
201 |
|
|
& + t1*(eosMDJWFden(6) |
202 |
|
|
& + eosMDJWFden(7)*t2) |
203 |
|
|
& + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) ) |
204 |
|
|
& + p1*(eosMDJWFden(10) |
205 |
|
|
& + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) ) |
206 |
|
|
|
207 |
|
|
rhoDen = 1.0/(epsln+den) |
208 |
|
|
|
209 |
|
|
rhoLoc = rhoNum*rhoDen - rhoConst |
210 |
|
|
|
211 |
|
|
ELSEIF( equationOfState .EQ. 'IDEALG' ) THEN |
212 |
|
|
C |
213 |
|
|
ELSE |
214 |
|
|
WRITE(msgBuf,'(3A)') |
215 |
|
|
& ' FIND_RHO_SCALAR : equationOfState = "', |
216 |
|
|
& equationOfState,'"' |
217 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
218 |
|
|
STOP 'ABNORMAL END: S/R FIND_RHO_SCALAR' |
219 |
|
|
ENDIF |
220 |
|
|
|
221 |
|
|
RETURN |
222 |
|
|
END |
223 |
|
|
|
224 |
|
|
C================================================================= |
225 |
|
|
|
226 |
|
|
_RL function SW_PTMP (S,T,P,PR) |
227 |
|
|
|
228 |
|
|
c ================================================================== |
229 |
|
|
c SUBROUTINE SW_PTMP |
230 |
|
|
c ================================================================== |
231 |
|
|
c |
232 |
|
|
c o Calculates potential temperature as per UNESCO 1983 report. |
233 |
|
|
c |
234 |
|
|
c started: |
235 |
|
|
c |
236 |
|
|
c Armin Koehl akoehl@ucsd.edu |
237 |
|
|
c |
238 |
|
|
c ================================================================== |
239 |
|
|
c SUBROUTINE SW_PTMP |
240 |
|
|
c ================================================================== |
241 |
|
|
C S = salinity [psu (PSS-78) ] |
242 |
|
|
C T = temperature [degree C (IPTS-68)] |
243 |
|
|
C P = pressure [db] |
244 |
|
|
C PR = Reference pressure [db] |
245 |
|
|
|
246 |
|
|
implicit none |
247 |
|
|
|
248 |
|
|
c routine arguments |
249 |
|
|
_RL S,T,P,PR |
250 |
|
|
|
251 |
|
|
c local arguments |
252 |
|
|
_RL del_P ,del_th, th, q |
253 |
|
|
_RL onehalf, two, three |
254 |
|
|
parameter ( onehalf = 0.5 _d 0, two = 2. _d 0, three = 3. _d 0 ) |
255 |
|
|
|
256 |
|
|
c externals |
257 |
|
|
_RL sw_adtg |
258 |
|
|
external sw_adtg |
259 |
|
|
c theta1 |
260 |
|
|
del_P = PR - P |
261 |
|
|
del_th = del_P*sw_adtg(S,T,P) |
262 |
|
|
th = T + onehalf*del_th |
263 |
|
|
q = del_th |
264 |
|
|
c theta2 |
265 |
|
|
del_th = del_P*sw_adtg(S,th,P+onehalf*del_P) |
266 |
|
|
|
267 |
|
|
th = th + (1 - 1/sqrt(two))*(del_th - q) |
268 |
|
|
q = (two-sqrt(two))*del_th + (-two+three/sqrt(two))*q |
269 |
|
|
|
270 |
|
|
c theta3 |
271 |
|
|
del_th = del_P*sw_adtg(S,th,P+onehalf*del_P) |
272 |
|
|
th = th + (1 + 1/sqrt(two))*(del_th - q) |
273 |
|
|
q = (two + sqrt(two))*del_th + (-two-three/sqrt(two))*q |
274 |
|
|
|
275 |
|
|
c theta4 |
276 |
|
|
del_th = del_P*sw_adtg(S,th,P+del_P) |
277 |
|
|
SW_PTMP = th + (del_th - two*q)/(two*three) |
278 |
|
|
return |
279 |
|
|
end |
280 |
|
|
|
281 |
|
|
C====================================================================== |
282 |
|
|
|
283 |
|
|
CBOP |
284 |
|
|
C !ROUTINE: SW_TEMP |
285 |
|
|
C !INTERFACE: |
286 |
|
|
_RL FUNCTION SW_TEMP( s, t, p, pr ) |
287 |
|
|
C !DESCRIPTION: \bv |
288 |
|
|
C *=============================================================* |
289 |
|
|
C | S/R SW_TEMP |
290 |
|
|
C | o compute in-situ temperature from potential temperature |
291 |
|
|
C *=============================================================* |
292 |
|
|
C |
293 |
|
|
C REFERENCES: |
294 |
|
|
C Fofonoff, P. and Millard, R.C. Jr |
295 |
|
|
C Unesco 1983. Algorithms for computation of fundamental properties of |
296 |
|
|
C seawater, 1983. _Unesco Tech. Pap. in Mar. Sci._, No. 44, 53 pp. |
297 |
|
|
C Eqn.(31) p.39 |
298 |
|
|
C |
299 |
|
|
C Bryden, H. 1973. |
300 |
|
|
C "New Polynomials for thermal expansion, adiabatic temperature gradient |
301 |
|
|
C and potential temperature of sea water." |
302 |
|
|
C DEEP-SEA RES., 1973, Vol20,401-408. |
303 |
|
|
C |
304 |
|
|
|
305 |
|
|
C !USES: |
306 |
|
|
IMPLICIT NONE |
307 |
|
|
|
308 |
|
|
C === Global variables === |
309 |
|
|
CML#include "SIZE.h" |
310 |
|
|
CML#include "EEPARAMS.h" |
311 |
|
|
CML#include "PARAMS.h" |
312 |
|
|
CML#include "GRID.h" |
313 |
|
|
CML#include "DYNVARS.h" |
314 |
|
|
CML#include "FFIELDS.h" |
315 |
|
|
CML#include "SHELFICE.h" |
316 |
|
|
|
317 |
|
|
C !INPUT/OUTPUT PARAMETERS: |
318 |
|
|
C === Routine arguments === |
319 |
|
|
C s :: salinity |
320 |
|
|
C t :: potential temperature |
321 |
|
|
C p :: pressure |
322 |
|
|
c pr :: reference pressure |
323 |
|
|
C myIter :: iteration counter for this thread |
324 |
|
|
C myTime :: time counter for this thread |
325 |
|
|
C myThid :: thread number for this instance of the routine. |
326 |
|
|
_RL s, t, p, pr |
327 |
|
|
_RL myTime |
328 |
|
|
INTEGER myIter |
329 |
|
|
INTEGER myThid |
330 |
|
|
CEOP |
331 |
|
|
|
332 |
|
|
C !LOCAL VARIABLES |
333 |
|
|
C === Local variables === |
334 |
|
|
_RL del_P ,del_th, th, q |
335 |
|
|
_RL onehalf, two, three |
336 |
|
|
PARAMETER ( onehalf = 0.5 _d 0, two = 2. _d 0, three = 3. _d 0 ) |
337 |
|
|
|
338 |
|
|
c externals |
339 |
|
|
_RL sw_adtg |
340 |
|
|
EXTERNAL sw_adtg |
341 |
|
|
c theta1 |
342 |
|
|
C-- here we swap P and PR in order to get in-situ temperature |
343 |
|
|
C del_P = PR - P ! to get potential from in-situ temperature |
344 |
|
|
del_P = P - PR ! to get in-situ from potential temperature |
345 |
|
|
del_th = del_P*sw_adtg(S,T,P) |
346 |
|
|
th = T + onehalf*del_th |
347 |
|
|
q = del_th |
348 |
|
|
c theta2 |
349 |
|
|
del_th = del_P*sw_adtg(S,th,P+onehalf*del_P) |
350 |
|
|
|
351 |
|
|
th = th + (1 - 1/sqrt(two))*(del_th - q) |
352 |
|
|
q = (two-sqrt(two))*del_th + (-two+three/sqrt(two))*q |
353 |
|
|
|
354 |
|
|
c theta3 |
355 |
|
|
del_th = del_P*sw_adtg(S,th,P+onehalf*del_P) |
356 |
|
|
th = th + (1 + 1/sqrt(two))*(del_th - q) |
357 |
|
|
q = (two + sqrt(two))*del_th + (-two-three/sqrt(two))*q |
358 |
|
|
|
359 |
|
|
c theta4 |
360 |
|
|
del_th = del_P*sw_adtg(S,th,P+del_P) |
361 |
|
|
SW_temp= th + (del_th - two*q)/(two*three) |
362 |
|
|
|
363 |
|
|
RETURN |
364 |
|
|
END |
365 |
|
|
|
366 |
|
|
C====================================================================== |
367 |
|
|
|
368 |
|
|
_RL function SW_ADTG (S,T,P) |
369 |
|
|
|
370 |
|
|
c ================================================================== |
371 |
|
|
c SUBROUTINE SW_ADTG |
372 |
|
|
c ================================================================== |
373 |
|
|
c |
374 |
|
|
c o Calculates adiabatic temperature gradient as per UNESCO 1983 routines. |
375 |
|
|
c |
376 |
|
|
c started: |
377 |
|
|
c |
378 |
|
|
c Armin Koehl akoehl@ucsd.edu |
379 |
|
|
c |
380 |
|
|
c ================================================================== |
381 |
|
|
c SUBROUTINE SW_ADTG |
382 |
|
|
c ================================================================== |
383 |
|
|
|
384 |
|
|
implicit none |
385 |
|
|
_RL a0,a1,a2,a3,b0,b1,c0,c1,c2,c3,d0,d1,e0,e1,e2 |
386 |
|
|
_RL S,T,P |
387 |
|
|
_RL sref |
388 |
|
|
|
389 |
|
|
sref = 35. _d 0 |
390 |
|
|
a0 = 3.5803 _d -5 |
391 |
|
|
a1 = +8.5258 _d -6 |
392 |
|
|
a2 = -6.836 _d -8 |
393 |
|
|
a3 = 6.6228 _d -10 |
394 |
|
|
|
395 |
|
|
b0 = +1.8932 _d -6 |
396 |
|
|
b1 = -4.2393 _d -8 |
397 |
|
|
|
398 |
|
|
c0 = +1.8741 _d -8 |
399 |
|
|
c1 = -6.7795 _d -10 |
400 |
|
|
c2 = +8.733 _d -12 |
401 |
|
|
c3 = -5.4481 _d -14 |
402 |
|
|
|
403 |
|
|
d0 = -1.1351 _d -10 |
404 |
|
|
d1 = 2.7759 _d -12 |
405 |
|
|
|
406 |
|
|
e0 = -4.6206 _d -13 |
407 |
|
|
e1 = +1.8676 _d -14 |
408 |
|
|
e2 = -2.1687 _d -16 |
409 |
|
|
|
410 |
|
|
SW_ADTG = a0 + (a1 + (a2 + a3*T)*T)*T |
411 |
|
|
& + (b0 + b1*T)*(S-sref) |
412 |
|
|
& + ( (c0 + (c1 + (c2 + c3*T)*T)*T) + (d0 + d1*T)*(S-sref) )*P |
413 |
|
|
& + ( e0 + (e1 + e2*T)*T )*P*P |
414 |
|
|
end |