1 |
C $Header: /u/gcmpack/MITgcm/model/src/seawater.F,v 1.7 2015/02/04 18:54:35 mmazloff Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
C-- File seawater.F: routines that compute quantities related to seawater. |
7 |
C-- Contents |
8 |
C-- o SW_PTMP: function to compute potential temperature |
9 |
C-- o SW_TEMP: function to compute in-situ temperature from pot. temp. |
10 |
C-- o SW_ADTG: function to compute adiabatic temperature gradient |
11 |
C-- (used by both SW_PTMP & SW_TEMP) |
12 |
|
13 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
14 |
|
15 |
CBOP |
16 |
C !ROUTINE: SW_PTMP |
17 |
C !INTERFACE: |
18 |
_RL FUNCTION SW_PTMP (S,T,P,PR) |
19 |
|
20 |
C !DESCRIPTION: \bv |
21 |
C *=============================================================* |
22 |
C | S/R SW_PTMP |
23 |
C | o compute potential temperature as per UNESCO 1983 report. |
24 |
C *=============================================================* |
25 |
C |
26 |
C started: |
27 |
C Armin Koehl akoehl@ucsd.edu |
28 |
C |
29 |
C ================================================================== |
30 |
C SUBROUTINE SW_PTMP |
31 |
C ================================================================== |
32 |
C S :: salinity [psu (PSS-78) ] |
33 |
C T :: temperature [degree C (IPTS-68)] |
34 |
C P :: pressure [db] |
35 |
C PR :: Reference pressure [db] |
36 |
|
37 |
C !USES: |
38 |
IMPLICIT NONE |
39 |
|
40 |
C !INPUT/OUTPUT PARAMETERS: |
41 |
_RL S,T,P,PR |
42 |
|
43 |
C !FUNCTIONS: |
44 |
_RL sw_adtg |
45 |
EXTERNAL sw_adtg |
46 |
|
47 |
C !LOCAL VARIABLES |
48 |
_RL del_P ,del_th, th, q |
49 |
_RL onehalf, two, three |
50 |
PARAMETER ( onehalf = 0.5 _d 0, two = 2. _d 0, three = 3. _d 0 ) |
51 |
CEOP |
52 |
|
53 |
C theta1 |
54 |
del_P = PR - P |
55 |
del_th = del_P*sw_adtg(S,T,P) |
56 |
th = T + onehalf*del_th |
57 |
q = del_th |
58 |
C theta2 |
59 |
del_th = del_P*sw_adtg(S,th,P+onehalf*del_P) |
60 |
|
61 |
th = th + (1 - 1/sqrt(two))*(del_th - q) |
62 |
q = (two-sqrt(two))*del_th + (-two+three/sqrt(two))*q |
63 |
|
64 |
C theta3 |
65 |
del_th = del_P*sw_adtg(S,th,P+onehalf*del_P) |
66 |
th = th + (1 + 1/sqrt(two))*(del_th - q) |
67 |
q = (two + sqrt(two))*del_th + (-two-three/sqrt(two))*q |
68 |
|
69 |
C theta4 |
70 |
del_th = del_P*sw_adtg(S,th,P+del_P) |
71 |
SW_PTMP = th + (del_th - two*q)/(two*three) |
72 |
|
73 |
RETURN |
74 |
END |
75 |
|
76 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
77 |
|
78 |
CBOP |
79 |
C !ROUTINE: SW_TEMP |
80 |
C !INTERFACE: |
81 |
_RL FUNCTION SW_TEMP( S, T, P, PR ) |
82 |
C !DESCRIPTION: \bv |
83 |
C *=============================================================* |
84 |
C | S/R SW_TEMP |
85 |
C | o compute in-situ temperature from potential temperature |
86 |
C *=============================================================* |
87 |
C |
88 |
C REFERENCES: |
89 |
C Fofonoff, P. and Millard, R.C. Jr |
90 |
C Unesco 1983. Algorithms for computation of fundamental properties of |
91 |
C seawater, 1983. _Unesco Tech. Pap. in Mar. Sci._, No. 44, 53 pp. |
92 |
C Eqn.(31) p.39 |
93 |
C |
94 |
C Bryden, H. 1973. |
95 |
C New Polynomials for thermal expansion, adiabatic temperature gradient |
96 |
C and potential temperature of sea water. |
97 |
C DEEP-SEA RES., 1973, Vol20,401-408. |
98 |
|
99 |
C !USES: |
100 |
IMPLICIT NONE |
101 |
C === Global variables === |
102 |
|
103 |
C !INPUT/OUTPUT PARAMETERS: |
104 |
C === Routine arguments === |
105 |
C S :: salinity |
106 |
C T :: potential temperature |
107 |
C P :: pressure |
108 |
C PR :: reference pressure |
109 |
_RL S, T, P, PR |
110 |
CEOP |
111 |
|
112 |
C !FUNCTIONS: |
113 |
_RL sw_adtg |
114 |
EXTERNAL sw_adtg |
115 |
|
116 |
C !LOCAL VARIABLES: |
117 |
_RL del_P ,del_th, th, q |
118 |
_RL onehalf, two, three |
119 |
PARAMETER ( onehalf = 0.5 _d 0, two = 2. _d 0, three = 3. _d 0 ) |
120 |
|
121 |
C theta1 |
122 |
C-- here we swap P and PR in order to get in-situ temperature |
123 |
C del_P = PR - P ! to get potential from in-situ temperature |
124 |
del_P = P - PR ! to get in-situ from potential temperature |
125 |
del_th = del_P*sw_adtg(S,T,PR) |
126 |
th = T + onehalf*del_th |
127 |
q = del_th |
128 |
C theta2 |
129 |
del_th = del_P*sw_adtg(S,th,PR+onehalf*del_P) |
130 |
|
131 |
th = th + (1 - 1/sqrt(two))*(del_th - q) |
132 |
q = (two-sqrt(two))*del_th + (-two+three/sqrt(two))*q |
133 |
|
134 |
C theta3 |
135 |
del_th = del_P*sw_adtg(S,th,PR+onehalf*del_P) |
136 |
th = th + (1 + 1/sqrt(two))*(del_th - q) |
137 |
q = (two + sqrt(two))*del_th + (-two-three/sqrt(two))*q |
138 |
|
139 |
C theta4 |
140 |
del_th = del_P*sw_adtg(S,th,PR+del_P) |
141 |
SW_temp= th + (del_th - two*q)/(two*three) |
142 |
|
143 |
RETURN |
144 |
END |
145 |
|
146 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
147 |
|
148 |
CBOP |
149 |
C !ROUTINE: SW_ADTG |
150 |
C !INTERFACE: |
151 |
_RL FUNCTION SW_ADTG (S,T,P) |
152 |
|
153 |
C !DESCRIPTION: \bv |
154 |
C *=============================================================* |
155 |
C | S/R SW_ADTG |
156 |
C | o compute adiabatic temperature gradient as per UNESCO 1983 routines. |
157 |
C *=============================================================* |
158 |
C |
159 |
C started: |
160 |
C Armin Koehl akoehl@ucsd.edu |
161 |
|
162 |
C !USES: |
163 |
IMPLICIT NONE |
164 |
|
165 |
C !INPUT/OUTPUT PARAMETERS: |
166 |
_RL S,T,P |
167 |
|
168 |
C !LOCAL VARIABLES: |
169 |
_RL a0,a1,a2,a3,b0,b1,c0,c1,c2,c3,d0,d1,e0,e1,e2 |
170 |
_RL sref |
171 |
CEOP |
172 |
|
173 |
sref = 35. _d 0 |
174 |
a0 = 3.5803 _d -5 |
175 |
a1 = +8.5258 _d -6 |
176 |
a2 = -6.836 _d -8 |
177 |
a3 = 6.6228 _d -10 |
178 |
|
179 |
b0 = +1.8932 _d -6 |
180 |
b1 = -4.2393 _d -8 |
181 |
|
182 |
c0 = +1.8741 _d -8 |
183 |
c1 = -6.7795 _d -10 |
184 |
c2 = +8.733 _d -12 |
185 |
c3 = -5.4481 _d -14 |
186 |
|
187 |
d0 = -1.1351 _d -10 |
188 |
d1 = 2.7759 _d -12 |
189 |
|
190 |
e0 = -4.6206 _d -13 |
191 |
e1 = +1.8676 _d -14 |
192 |
e2 = -2.1687 _d -16 |
193 |
|
194 |
SW_ADTG = a0 + (a1 + (a2 + a3*T)*T)*T |
195 |
& + (b0 + b1*T)*(S-sref) |
196 |
& + ( (c0 + (c1 + (c2 + c3*T)*T)*T) + (d0 + d1*T)*(S-sref) )*P |
197 |
& + ( e0 + (e1 + e2*T)*T )*P*P |
198 |
|
199 |
RETURN |
200 |
END |