source: trunk/tools/density_binning/binning_neutral_and_co/bn2.pro @ 163

Last change on this file since 163 was 163, checked in by pinsard, 15 years ago

going on modify unformal header to idldoc 2. header syntax

File size: 3.9 KB
Line 
1;+
2;
3; Compute the local Brunt-Vaisala frequency
4;       
5;  Method :
6;  -------
7;       neos = 0  : UNESCO sea water properties
8;          The brunt-vaisala frequency is computed using the polynomial
9;       polynomial expression of McDougall (1987):
10;               N^2 = g * beta * ( alpha/beta*dk[ tn ] - dk[ sn ] )/e3w
11;
12;       neos = 1  : linear equation of state (temperature only)
13;               N^2 = g * ralpha * dk[ tn ]/e3w
14;
15;       neos = 2  : linear equation of state (temperature & salinity)
16;               N^2 = g * (ralpha * dk[ tn ] - rbeta * dk[ sn ] ) / e3w
17;
18;       The use of potential density to compute N^2 introduces e r r o r
19;      in the sign of N^2 at great depths. We recommand the use of
20;      neos = 0, except for academical studies.
21;
22;      N.B. N^2 is set to zero at the first level (JK=1) in inidtr
23;      and is never used at this level.
24;
25; @param TN {in}
26; potential temperature
27;
28; @param SN {in}
29; salinity sn
30;
31; @returns
32; brunt-vaisala frequency
33;
34;   References :
35;   -----------
36;       McDougall, T. J., J. Phys. Oceanogr., 17, 1950-1964, 1987.
37;
38; @history
39; Original  : 94-07 (G. Madec, M. Imbard)
40; Additions : 97-07 (G. Madec) introduction of statement functions
41; Idl version : 98-12 (M. Levy)
42;
43; @version
44; $Id$
45;
46;-
47FUNCTION bn2, tn,  sn, neos
48@common
49   case n_params() of
50      0:begin
51         print, 'rentrer nom des champs de temperature et/ou salinite'
52         goto, sortie
53      end
54      1:begin
55         print, 'pas de salinite : valeur de neos imposee a 1'
56         neos = 1
57      end
58      2:begin
59         print, 'calcul Unesco par default'
60         neos = 0
61      END
62      3:
63   endcase
64domdef
65;  coordonnees z : on s'assure qu'on a autant de points W que de points
66;  T sur le domaine defini par domdef.
67;  on va se servir dans les calculs de e3w (e3) et gdepw (gdep)
68   IF (nzt NE nzw ) THEN BEGIN
69      print, 'le domaine doit etre choisi avec autant de points'
70      print, 'sur  niveaux T que sur les niveaux W pour calculer N2'
71      z = -1
72      GOTO,  sortie
73   endif
74   IF (nzt LE 1 ) THEN BEGIN
75      print, 'il faut un minimum de deux points sur la verticale'
76      print, ' pour calculer N2'
77      z = -1
78      GOTO,  sortie
79   endif
80   cdgrid = 'W'
81   g = 9.81
82   grille, mask, glam, gphi, gdep, nx, ny,nz, premierx, premiery, premierz, dernierx,  derniery,  dernierz ,e1,e2,e3, TRI = tri
83   zgde3w = reform( (fltarr(nx*ny)+1)#(g/e3),  nx, ny, nz)
84;
85varname='Brunt-Vaisala frequency'
86varunit='s-1'
87;vargrid='W'
88;
89; mask tn and sn
90        tn=tn*tmask
91        sn=sn*tmask
92;
93;  Interior points ( 2=< jk =< jpkm1 )
94;  ----------------
95;
96   case 1 of
97      (neos EQ 0): begin
98; ... UNESCO seawater equation of state
99;   ... temperature, salinity anomalie
100         zzt = 0.5*( tn + shift( tn, 0, 0,  1 ))
101         zzs = 0.5*( sn + shift( sn, 0, 0,  1) ) - 35.0
102         zzp =  reform( (fltarr(nx*ny)+1)#gdep, nx , ny, nz)
103;   ... ratio alpha/beta and beta
104         zalbet = fsalbt( zzt, zzs, zzp )
105         zbeta  = fsbeta( zzt, zzs, zzp )
106;   ... N^2 = g/e3w * beta * ( alpha/beta * dk[tn] - dk[sn])
107         z = zgde3w * zbeta * mask $
108          * ( zalbet * ( shift(tn, 0, 0, 1) - tn ) $
109              - ( shift(sn, 0, 0, 1) - sn ) )
110      end
111      ( neos EQ 1 ): begin
112; ... First Linear density formulation (function of temperature only)
113;
114        ralpha  = 1.176e-4
115        z = zgde3w * ralpha * mask* (shift(tn, 0, 0, 1) - tn)
116         
117      end
118      ( neos EQ 2 ): begin
119; ... Second linear density formulation (function of temp. and salinity)
120;
121         ralpha = 2.e-4
122         rbeta = 0.001     
123         z = zgde3w * mask $
124          * (  ralpha * (shift(tn, 0, 0, 1) - tn) $
125               - rbeta  * (shift(sn, 0, 0, 1) - sn)  )
126      end   
127      (neos NE 0) AND (neos NE 1) AND (neos NE 2): BEGIN
128         print,  'neos doit etre egal a 0, 1 ou 2'
129         GOTO,  sortie
130      end
131   endcase
132;
133;
134;  first and last levels
135; ------------------------
136; bn^2=0. at first vertical point and jk=jpk
137   z(*, *, 0) = 0
138   IF (gdep(nz-1) eq gdepw(jpk-1)) THEN z(*, *, nz-1) = 0
139sortie:
140   return,  z
141END
Note: See TracBrowser for help on using the repository browser.