source: trunk/tools/density_binning/density_bin_IDL_gm/bining2.pro @ 161

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

conversion from IDL help header syntax to IDLDOC 2. header syntax

File size: 10.0 KB
Line 
1PRO bining2, $
2   density, $
3   x1, $
4   sobwlmax, $
5   sig_bowl, $
6   depth_bin, $
7   thick_bin, $
8   x1_bin, $
9   bowl_bin, $
10   SIGMA = sigma, $
11   DEPTH_T = depth_t, $
12   DEPTH_W = depth_w,  $
13   E3T = e3t,          $
14   E3W = e3w,          $
15   TMASK = tmask
16;;
17;;  PRINCIPE DU BINING
18;;
19;;  - etant donnee une liste de densite (sigma) on cherche les  profondeurs
20;;    correspondants a ces valeurs. Pour cela on utilise une interpolation lineaire.
21;;
22;;  - puis, connaissant les profondeurs des isopycnes, on interpole les valeurs du contenu
23;;    vertical de x a ces profondeurs. Il est preferable d' interpoler le contenu de x plutot
24;;    que x afin d'assurer la conservation de x.
25;;
26;;  - on deduit la valeur moyenne de x dans chaque couche (i.e. entre deux isopycnes). Il ne
27;;    s'agit donc pas de la projection de x selon les isopycnes mais de la valeur moyenne entre
28;;    deux isopycnes. Les proprietes sont meilleures (conservation de la quantite totale de x).
29;;
30;; REMARQUES
31;;
32;;  - il n'y a aucun appel a un fichier de common (d'ou la liste des mots clef un peu longue)
33;;
34;;  - le calcul n'est rigoureux que pour les points T. Le bining de la vitesse doit etre adapte...
35;;
36;;  - si l'utilisateur donne N interfaces, cela defini N+1 couches, d'ou les decalages dans les indices
37;;    et l'apparition d'un +1 dans la declaration des tableaux de sortie
38;;    la couche 0   contient tous les points dont la densite est inferieure a la densite minimale du bining
39;;    la couche N+1  "                   "                     " superieure      "       maximale
40;;
41;;
42; @param DENSITY {in}
43; density at T points (3D array)
44;
45; @param X1 {in}
46; field at T points (e.g. Salinity or Temperature)
47;
48; @param SOBWLMAX {in}
49; bowl depth array
50;
51; @param SIG_BOWL {in}
52; switch for bowl overlay
53;
54; @param DEPTH_BIN {out}
55; depth of layers (3D array)
56;
57; @param THICK_BIN {out}
58; thickness of layers (3D array)
59;
60; @param X1_BIN {out}
61; X averaged for each sigma-layer (3D array)
62;
63; @param BOWL_BIN {out}
64; bowl depth binned on density
65;
66; @keyword SIGMA {in}{optional}
67; bining values
68;
69; @keyword DEPTH_T {in}{required}
70; depth of T level
71;
72; @keyword DEPTH_W {in}{required}
73; depth of W level
74;
75; @keyword E3T
76;
77; @keyword E3W
78;
79; @keyword TMASK{in}{required}
80; tmask (3D array)
81;
82; @history
83; 19/10/99 G. Roullet
84;-
85
86   size3d = size(density)
87   jpi = size3d(1)
88   jpj = size3d(2)
89   jpk = size3d(3)
90;print,  'size3d=', size3d
91;print,  'size(tmask)=',  size(tmask)
92
93   indm = where(tmask eq 0)
94
95   x1(indm) = !VALUES.F_NAN
96   density(indm) = !VALUES.F_NAN
97
98   N_z = jpk
99   
100; adjust profiles to avoid static instabilities (there are almost none IF using pac_40)
101IF jpi EQ 180 THEN BEGIN
102   density =  npc(density)
103ENDIF
104
105;  N_s, defini ci dessous est le nombre d'interface
106
107   IF keyword_set(sigma) THEN BEGIN
108      s_s = sigma
109      N_s = n_elements(s_s)
110   ENDIF ELSE BEGIN
111;
112;  Definition d'un bining arbitraire (si density est un sigma potentiel, cela repere
113;   les eaux intermediaires)
114;
115      N_s = 20
116      s_s = 27.6+0.2*findgen(N_s)/(N_s-1)
117      N_s = 3
118      s_s = [26.8, 27.6, 27.8]
119      s_s = [22.8, 27.6, 27.8]
120   ENDELSE
121
122;   print, 'Density bining :', s_s
123;;
124;; Definition des variables
125;;
126;   Profils selon les niveaux du modele (suffixe _z)
127;
128   c1_z = fltarr(N_z)    ; profil du contenu vertical de x
129   s_z = fltarr(N_z)    ; profil de la densite
130   z_zt = depth_t  ; profondeur au point T  (k=0 -> 5m)
131   z_zw = depth_w  ;                     W  (k=0 -> 10m)
132
133;   Profils selon les couches de densite (suffixe _s)
134;
135   z_s  = fltarr(N_s+1)
136   c1_s = fltarr(N_s+1)
137   x1_s = fltarr(N_s+1)
138   bowl_s = 0.
139
140; Tableaux de sorties
141;
142   depth_bin  = fltarr(jpi, jpj, N_s+1)
143   thick_bin  = fltarr(jpi, jpj, N_s+1)
144   x1_bin = fltarr(jpi, jpj, N_s+1)
145   bowl_bin = fltarr(jpi, jpj)
146
147   x1_bin(*, *, *) = !VALUES.F_NAN
148   depth_bin(*, *, *) = !VALUES.F_NAN
149   thick_bin(*, *, *) = !VALUES.F_NAN
150   bowl_bin(*, *) = !VALUES.F_NAN
151
152
153;  Calcul du contenu de x sur la verticale (permet d assurer la conservation integrale)
154;
155   x1_content = fltarr(jpi, jpj, jpk)
156
157
158;   x1_content(*, *, 0) = e3t(0) * x1(*, *, 0) * tmask (*, *, 0)
159;   FOR k = 1, jpk-1 DO x1_content(*, *, k) = x1_content(*, *, k-1) $
160;                        + e3t(k) * x1(*, *, k)* tmask (*, *, k)
161   x1_content =  x1
162
163;
164;  Boucle sur le domaine 2D
165;
166
167   FOR i = 0, (jpi-1) DO BEGIN
168;  FOR i = 20, 20 DO BEGIN
169;   FOR i = 162, 162 DO BEGIN
170      FOR j = 0, (jpj-1) DO BEGIN
171;     FOR j = 20, 20 DO BEGIN
172;      FOR j = 116, 116 DO BEGIN
173;
174;  Indices des points T dans l ocean
175;print,  '    '
176;print,  ' ==============================='
177;print, '    POINT I,J =', i, j
178;
179;stop
180         i_ocean = where(tmask(i, j, *) EQ 1)
181;print, 'tmask(i, j, *)=', tmask(i, j, *)
182;print, 'i_ocean=', i_ocean
183         z_s =  z_s*0.
184         c1_s(*) = !VALUES.F_NAN
185         bowl_s = !VALUES.F_NAN
186;         x1_s(*) = !VALUES.F_NAN
187
188
189         IF i_ocean[0] NE -1 THEN BEGIN ; on n entre que si il y a des points ocean
190;            print,  'ocean point'
191            i_bottom = i_ocean(n_elements(i_ocean)-1)
192;print,  'i_bottom=',  i_bottom
193            z_s(N_s) = z_zw(i_bottom)
194            c1_s(N_s) = x1_content(i, j, jpk-1)
195           
196            s_z(*) = density(i, j, *)
197            c1_z(*) = x1_content(i, j, *)
198;            print, 'density profile s_z', s_z
199;            print, 'field profile c1_z', c1_z
200
201;            print, 'bins s_s', s_s
202
203; extraction d'un sous profil strictement croissant
204;            print,  's_z(i_ocean)=',  s_z(i_ocean)
205            mini = min(s_z(i_ocean))
206            maxi = max(s_z(i_ocean))
207            i_min = where(s_z(i_ocean) EQ mini)
208            i_max = where(s_z(i_ocean) EQ maxi)
209;            print, 'mini, maxi', mini, maxi
210;            print, 'i_min, i_max', i_min, i_max
211;   on prend le plus grand des indices min
212;         et le plus petit des indices max
213            i_min = i_min(n_elements(i_min)-1)
214            i_max = i_max[0]
215;            print, 'i_min, i_max', i_min, i_max
216;IF i_min GT i_max THEN print, 'WARNING: i_min, i_max=', i_min, i_max,  ' at i,j=,', i, j
217;   on prend le plus grand des indices min
218
219;            IF i_max GE jpk-1 THEN print, i, j, i_max
220;            IF i_min LE 1 THEN print, i, j, i_min
221           
222; Si la valeur du niveau (s_s) est plus faible que la densite de surface,
223; l isopycne est mise en surface (z_s=0)
224;
225;            print, 's_z(i_min)', s_z(i_min)
226;print,  's_s=',  s_s
227            ind = where(s_s LT s_z(i_min))
228            IF ind[0] NE -1 THEN BEGIN
229;               IF i_min GT i_max THEN print, 'min reached at sigma indices', ind
230               z_s(ind) = 0
231               c1_s(ind) = !VALUES.F_NAN
232            ENDIF
233
234; Si la valeur du niveau (s_s) est plus elevee que la densite du fond,
235; l isopycne est mise au fond (z_s=z_zw(i_bottom))
236;
237;            print, 's_z(i_max)', s_z(i_max)
238            ind = where(s_s GT s_z(i_max))
239
240            IF ind[0] NE -1 THEN BEGIN
241;               IF i_min GT i_max THEN print, 'max reached at sigma indices', ind
242               z_s(ind) = z_s(N_s)
243               c1_s(ind) = c1_s(N_s)
244               c1_s(ind) = !VALUES.F_NAN
245            ENDIF
246; cas general
247            ind = where( (s_s GE s_z(i_min)) AND (s_s LE s_z(i_max)) )
248;            IF i_min GT i_max THEN print, 's_s indices inside min/max', ind
249            IF ind[0] NE -1 THEN BEGIN
250               i_profil = i_ocean(i_min:i_max) ; problem line
251
252               z_s(ind) = interpol(z_zt(i_profil), s_z(i_profil), s_s(ind)) ; original
253;               z_s(ind) = interpol(z_zw(i_profil), s_z(i_profil), s_s(ind)) ; changed to z_zw 1/7/04
254
255;
256; j'utilise la fonction spline pour interpoler le contenu
257; l'interpolation lineaire marche aussi. Elle donne des resultats differents. Elle
258; me semble moins propre. Je la donne en remarque.
259;
260               IF n_elements(i_profil) GT 2 THEN BEGIN
261;                  c1_s(ind) = spline(z_zw(i_profil), c1_z(i_profil), z_s(ind), 1)     ; cubic spline
262                  c1_s(ind) = interpol(c1_z(i_profil), z_zt(i_profil), z_s(ind))      ; linear interpolation
263;                   SPLINE_P, z_zw(i_profil), c1_z(i_profil), z_s(ind),   c1_s(ind)     ; generalization of spline(*)
264;print,  z_zw(i_profil), c1_z(i_profil),  z_s(ind), c1_s(ind)
265               ENDIF ELSE BEGIN
266                  c1_s(ind) = interpol(c1_z(i_profil), z_zt(i_profil), z_s(ind))
267               ENDELSE
268;
269;     bowl depth binning
270
271               IF sig_bowl EQ 1 THEN BEGIN
272                  bowl_s = interpol(s_z(i_profil), z_zt(i_profil),  sobwlmax(i, j))
273               ENDIF ELSE bowl_s = 0
274                 
275;stop
276;
277;
278;               x1_s(ind) = (c1_s(ind+1)-c1_s(ind))/(z_s(ind+1)-z_s(ind))
279;                     x1_s = c1_s
280
281            ENDIF 
282;ELSE print, 'ind[0] = -1',  ind,  i_bottom
283         ENDIF 
284;ELSE print, ' land ', tmask(i, j, 1)
285         depth_bin (i, j, *) = z_s
286         thick_bin (i, j, 0) = z_s(0)
287         thick_bin (i, j, 1:N_s) = z_s(1:N_s)-z_s(0:N_s-1)
288         IF total(thick_bin (i, j, 1:N_s) LT 0) GT 0 THEN BEGIN
289            ;print, 'WARNING: negative thick_bin (i, j, 1:N_s) at i,j= ',  i, j
290            ;print, 'depth_bin (i, j, 1:N_s)= ',  depth_bin (i, j, 1:N_s)
291            ;print, 'thick_bin (i, j, 1:N_s)= ',  thick_bin (i, j, 1:N_s)
292            ;print, 's_z= ',  s_z
293            ;stop
294         ENDIF
295         x1_bin    (i, j, *) = c1_s
296         bowl_bin  (i, j) = bowl_s
297;         print, 'c1_s', c1_s
298;
299      ENDFOR
300   ENDFOR
301; mask depth_bin with undefined values of x_bin
302   depth_bin(where(finite(x1_bin, /nan) EQ 1)) = !VALUES.F_NAN
303   thick_bin(where(finite(x1_bin, /nan) EQ 1)) = !VALUES.F_NAN
304
305; OPTIONAL: compute spiciness by removing isopycnal mean for each sigma (T or S)
306   IF 0 THEN BEGIN
307      print, '    Computing spiciness...'
308      FOR i = 0, (size(x1_bin))[3]-1 DO BEGIN
309         x1_bin(*,*,i) = x1_bin(*,*,i) - mean(x1_bin(*,*,i),/NaN)
310      ENDFOR
311   ENDIF
312
313; OPTIONAL: remove domain mean (S)
314   IF 0 THEN BEGIN
315      print, '    Removing domain mean...'
316      x1_bin(*,*,*) = x1_bin(*,*,*) - mean(x1_bin(*,*,*),/NaN)
317   ENDIF
318
319; OPTIONAL: remove 34.6psu (S)
320   IF 1 THEN BEGIN
321      print, '    Removing 34.6psu...'
322      x1_bin(*,*,*) = x1_bin(*,*,*) - 34.6
323   ENDIF
324
325END
Note: See TracBrowser for help on using the repository browser.