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