NCCOOS Trac Projects: Top | Web | Platforms | Processing | Viz | Sprints | Sandbox | (Wind)

root/gliderproc/trunk/MATLAB/seawater/sw_dens.m

Revision 495 (checked in by cbc, 12 years ago)

Initial import of Stark code.

Line 
1
2 function dens = sw_dens(S,T,P)
3
4 % SW_DENS    Density of sea water
5 %=========================================================================
6 % SW_DENS  $Id: sw_dens.m,v 1.1 2003/12/12 04:23:22 pen078 Exp $
7 %          Copyright (C) CSIRO, Phil Morgan 1992.
8 %
9 % USAGE:  dens = sw_dens(S,T,P)
10 %
11 % DESCRIPTION:
12 %    Density of Sea Water using UNESCO 1983 (EOS 80) polynomial.
13 %
14 % INPUT:  (all must have same dimensions)
15 %   S = salinity    [psu      (PSS-78)]
16 %   T = temperature [degree C (ITS-90)]
17 %   P = pressure    [db]
18 %       (P may have dims 1x1, mx1, 1xn or mxn for S(mxn) )
19 %
20 % OUTPUT:
21 %   dens = density  [kg/m^3]
22 %
23 % AUTHOR:  Phil Morgan 92-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
24 %
25 % DISCLAIMER:
26 %   This software is provided "as is" without warranty of any kind.
27 %   See the file sw_copy.m for conditions of use and licence.
28 %
29 % REFERENCES:
30 %    Fofonoff, P. and Millard, R.C. Jr
31 %    Unesco 1983. Algorithms for computation of fundamental properties of
32 %    seawater, 1983. _Unesco Tech. Pap. in Mar. Sci._, No. 44, 53 pp.
33 %
34 %    Millero, F.J., Chen, C.T., Bradshaw, A., and Schleicher, K.
35 %    " A new high pressure equation of state for seawater"
36 %    Deap-Sea Research., 1980, Vol27A, pp255-264.
37 %=========================================================================
38
39 % Modifications
40 % 99-06-25. Lindsay Pender, Fixed transpose of row vectors.
41 % 03-12-12. Lindsay Pender, Converted to ITS-90.
42
43 % CALLER: general purpose
44 % CALLEE: sw_dens0.m sw_seck.m
45
46 % UNESCO 1983. eqn.7  p.15
47
48 %----------------------
49 % CHECK INPUT ARGUMENTS
50 %----------------------
51 if nargin ~=3
52    error('sw_dens.m: Must pass 3 parameters')
53 end %if
54
55 % CHECK S,T,P dimensions and verify consistent
56 [ms,ns] = size(S);
57 [mt,nt] = size(T);
58 [mp,np] = size(P);
59
60
61 % CHECK THAT S & T HAVE SAME SHAPE
62 if (ms~=mt) | (ns~=nt)
63    error('check_stp: S & T must have same dimensions')
64 end %if
65
66 % CHECK OPTIONAL SHAPES FOR P
67 if     mp==1  & np==1      % P is a scalar.  Fill to size of S
68    P = P(1)*ones(ms,ns);
69 elseif np==ns & mp==1      % P is row vector with same cols as S
70    P = P( ones(1,ms), : ); %   Copy down each column.
71 elseif mp==ms & np==1      % P is column vector
72    P = P( :, ones(1,ns) ); %   Copy across each row
73 elseif mp==ms & np==ns     % PR is a matrix size(S)
74    % shape ok
75 else
76    error('check_stp: P has wrong dimensions')
77 end %if
78
79 %***check_stp
80
81 %------
82 % BEGIN
83 %------
84 densP0 = sw_dens0(S,T);
85 K      = sw_seck(S,T,P);
86 P      = P/10;  % convert from db to atm pressure units
87 dens   = densP0./(1-P./K);
88 return
89 %--------------------------------------------------------------------
90
91
Note: See TracBrowser for help on using the browser.