-
Notifications
You must be signed in to change notification settings - Fork 0
/
gsw_CT_maxdensity.m
executable file
·124 lines (109 loc) · 4.71 KB
/
gsw_CT_maxdensity.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
function CT_maxdensity = gsw_CT_maxdensity(SA,p)
% gsw_CT_maxdensity Conservative Temperature of maximum
% density of seawater (75-term equation)
% =========================================================================
%
% USAGE:
% CT_maxdensity = gsw_CT_maxdensity(SA,p)
%
% DESCRIPTION:
% Calculates the Conservative Temperature of maximum density of seawater.
% This function returns the Conservative temperature at which the density
% of seawater is a maximum, at given Absolute Salinity, SA, and sea
% pressure, p (in dbar). This function uses the computationally-efficient
% expression for specific volume in terms of SA, CT and p
% (Roquet et al., 2015).
%
% Note that the 75-term equation has been fitted in a restricted range of
% parameter space, and is most accurate inside the "oceanographic funnel"
% described in McDougall et al. (2003). The GSW library function
% "gsw_infunnel(SA,CT,p)" is avaialble to be used if one wants to test if
% some of one's data lies outside this "funnel".
%
% INPUT:
% SA = Absolute Salinity [ g/kg ]
% p = sea pressure [ dbar ]
% ( i.e. absolute pressure - 10.1325 dbar )
%
% p may have dimensions 1x1 or Mx1 or 1xN or MxN, where SA is MxN.
%
% OUTPUT:
% CT_maxdensity = Conservative Temperature at which [ deg C ]
% the density of seawater is a maximum for
% given Absolute Salinity and pressure.
%
% AUTHOR:
% Trevor McDougall & Paul Barker [ [email protected] ]
%
% VERSION NUMBER: 3.05 (27th January 2015)
%
% REFERENCES:
% IOC, SCOR and IAPSO, 2010: The international thermodynamic equation of
% seawater - 2010: Calculation and use of thermodynamic properties.
% Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
% UNESCO (English), 196 pp. Available from http://www.TEOS-10.org
% See section 3.42 of this TEOS-10 Manual.
%
% McDougall T. J. and S. J. Wotherspoon, 2013: A simple modification of
% Newton's method to achieve convergence of order 1 + sqrt(2). Applied
% Mathematics Letters, 29, 20-25.
%
% Roquet, F., G. Madec, T.J. McDougall, P.M. Barker, 2015: Accurate
% polynomial expressions for the density and specifc volume of seawater
% using the TEOS-10 standard. Ocean Modelling.
%
% The software is available from http://www.TEOS-10.org
%
%==========================================================================
%--------------------------------------------------------------------------
% Check variables and resize if necessary
%--------------------------------------------------------------------------
if ~(nargin == 2)
error('gsw_CT_maxdensity: Requires two inputs')
end %if
[ms,ns] = size(SA);
[mp,np] = size(p);
if (mp == 1) & (np == 1) % p scalar - fill to size of SA
p = p*ones(size(SA));
elseif (ns == np) & (mp == 1) % p is row vector,
p = p(ones(1,ms), :); % copy down each column.
elseif (ms == mp) & (np == 1) % p is column vector,
p = p(:,ones(1,ns)); % copy across each row.
elseif (ns == mp) & (np == 1) % p is a transposed row vector,
p = p.'; % transposed then
p = p(ones(1,ms), :); % copy down each column.
elseif (ms == mp) & (ns == np)
% ok
else
error('gsw_CT_maxdensity: Inputs array dimensions arguments do not agree')
end %if
if ms == 1
SA = SA.';
p = p.';
transposed = 1;
else
transposed = 0;
end
%--------------------------------------------------------------------------
% Start of the calculation
%--------------------------------------------------------------------------
dCT = 0.001; % the Conservative Temperature increment.
CT = 3.978 - 0.22072*SA; % the initial guess of CT.
dalpha_dCT = 1.1e-5; % the initial guess for d(alpha)_dCT.
for Number_of_iterations = 1:3
CT_old = CT;
alpha = gsw_alpha(SA,CT_old,p);
CT = CT_old - alpha./dalpha_dCT; % this is half way through the modified method
CT_mean = 0.5*(CT + CT_old);
dalpha_dCT = (gsw_alpha(SA,CT_mean + dCT,p) ...
- gsw_alpha(SA,CT_mean - dCT,p))./(dCT + dCT);
CT = CT_old - alpha./dalpha_dCT;
end
% After three iterations of this modified Newton-Raphson (McDougall and
% Wotherspoon, 2013) iteration, the error in CT_maxdensity is typically no
% larger than 1x10^-15 degress C.
CT_maxdensity = CT;
if transposed
CT_maxdensity = CT_maxdensity.';
end
end