source: trunk/SRC/ReadWrite/ncdf_getaxis.pro @ 224

Last change on this file since 224 was 224, checked in by pinsard, 17 years ago

improvements/corrections of some *.pro headers

  • Property svn:keywords set to Id
File size: 9.0 KB
Line 
1;+
2; @file_comments
3; get the x/y dimension Id and x/y axes from a netcdf file
4;
5; @categories
6; Reading
7;
8; @param cdfid {in}{required}{type=scalar}
9; the id of the netcdf file
10;
11; @param dimidx {out}{type=scalar (long)}
12; id of the x dimension
13;
14; @param dimidy {out}{type=scalar (long)}
15; id of the y dimension
16;
17; @param xaxis {out}{type=1D or 2D array}
18; the x axis
19;
20; @param yaxis {out}{type=1D or 2D array}
21; the y axis
22;
23; @keyword ROMSGRID {out}{type=scalar: 0 or 1}
24; gives back if we are using a ROMS grid (1) or not (0)
25;
26; @keyword START1 {default=0}{type=scalar: 0 or 1}
27; Index the axis from 1 instead of 0 when using /xyindex
28;
29; @keyword XDIMNAME {default='longitude', 'lon', 'x', 'longitude*', 'lon*', 'x*', '*longitude*', '*lon*' or '*x*'}{type=scalar string}
30; A string giving the name of the x dimension
31;
32; @keyword YDIMNAME {default='latitude', 'lat', 'y', 'latitude*', 'lat*', 'y*', 'eta_*', '*latitude*', '*lat*', '*y*'}{type=scalar string}
33; A string giving the name of the y dimension
34;
35; @keyword XAXISNAME {default='x', 'longitude', 'nav_lon', 'lon', 'lon_rho' or 'NbLongitudes'}{type=scalar string}
36; A string giving the name of the variable in the file
37; that contains the [xyz]axis.
38;
39; @keyword YAXISNAME {default='y', 'latitude', 'nav_lat','lat', 'lat_rho' or 'NbLatitudes'}{type=scalar string}
40; A string giving the name of the variable in the file
41; that contains the [xyz]axis.
42;
43; @keyword XYINDEX {default=0}{type=scalar: 0 or 1}
44; To define the x/y axis with index instead of using
45; the values contained in X/YAXISNAME.
46; x/yaxis = keyword_set(start1) + findgen(jpi/jpj)
47;
48; @keyword XYINDEX {default=0}{type=scalar: 0 or 1}
49;
50; @keyword _EXTRA
51; Used to be able to call ncdf_getaxis with _extra
52;
53; @history
54; March 2007: Sebastien Masson (smasson\@locean-ipsl.upmc.fr)
55;
56;
57; @version
58; $Id$
59;-
60PRO ncdf_getaxis, cdfid, dimidx, dimidy, xaxis, yaxis $
61                  , XAXISNAME = xaxisname, YAXISNAME = yaxisname $
62                  , XDIMNAME = xdimname, YDIMNAME = ydimname $
63                  , XYINDEX = xyindex, START1 = start1, ROMSGRID = romsgrid, _extra = ex
64;
65  compile_opt idl2, strictarrsubs
66;
67
68; what is inside the file
69  inside = ncdf_inquire(cdfid)
70;------------------------------------------------------------
71; name of all dimensions
72  namedim = strarr(inside.ndims)
73  for dimiq = 0, inside.ndims-1 do begin
74    ncdf_diminq, cdfid, dimiq, tmpname, value
75    namedim[dimiq] = strlowcase(tmpname)
76  ENDFOR
77;----------------------------------------------------------
78; name of the variables
79  namevar = strarr(inside.nvars)
80  for varid = 0, inside.nvars-1 do begin
81    invar = ncdf_varinq(cdfid, varid)
82    namevar[varid] = strlowcase(invar.name)
83  ENDFOR
84;----------------------------------------------------------
85; find the xaxis
86; try to get the variable that contains the xaxis
87  if keyword_set(xaxisname) then xaxisname = strlowcase(xaxisname) ELSE xaxisname = 'x'
88  xvarid = (where(namevar EQ xaxisname OR namevar EQ 'longitude' $
89                  OR namevar EQ 'nav_lon' OR namevar EQ 'lon' $
90                  OR namevar EQ 'lon_rho' OR namevar EQ 'nblongitudes'))[0]
91; no xaxis variable found, we will build a fake xaxis based on the size of the x dimension
92; -> we must find the x dimension
93  IF xvarid EQ -1 THEN BEGIN
94    dummy = report(['xaxis variable was not found within the default names:' $
95                    , '''longitude'', ''nav_lon'', ''lon'', ''lon_rho'', ''nblongitudes''' $
96                    , ' we use a fake xaxis based on x dimension size (or use XAXISNAME keyword)'], /simple)
97; try to get the dimension corresponding to x
98; roms file?
99    dimidx = where(namedim EQ 'xi_rho' OR namedim EQ 'xi_u' OR namedim EQ 'xi_v' OR namedim EQ 'xi_psi')
100    IF dimidx[0] EQ -1 THEN BEGIN
101; we are looking for a x dimension with a name matching one of the following regular expression:
102      if keyword_set(xdimname) then testname = strlowcase(xdimname) $
103      ELSE testname = ['longitude', 'lon', 'x', 'longitude*', 'lon*', 'x*', '*longitude*', '*lon*', '*x*']
104      cnt = -1
105      ii = 0
106      WHILE cnt NE 1 AND ii LT n_elements(testname) DO BEGIN
107        dimidx = where(strmatch(namedim, testname[ii]) EQ 1, cnt)
108        ii = ii+1
109      ENDWHILE
110      CASE cnt OF
111        0:begin
112          dummy = report(['none of the dimensions name matches one of the following regular expression:' $
113                          , '''longitude'', ''lon'', ''x'', ''longitude*'', ''lon*'', ''x*'', ''*longitude*'', ''*lon*'', ''*x*''' $
114                          , ' => we cannot find the x dimension, use XDIMNAME keyword'], /simple)
115          stop
116        END
117        1:dimidx = dimidx[0]
118        ELSE:begin
119          dummy = report(['several (and not one unique) dimensions name matches the following regular expression:' $
120                          , '''longitude'', ''lon'', ''x'', ''longitude*'', ''lon*'', ''x*'', ''*longitude*'', ''*lon*'', ''*x*''' $
121                          , ' => we cannot find the x dimension, use XDIMNAME keyword'], /simple)
122          stop
123        ENDELSE
124      ENDCASE
125    ENDIF
126    romsgrid = 0b
127  ENDIF ELSE BEGIN
128    romsgrid = strmid(namevar[xvarid], 0, 4) EQ 'lon_'
129    xinq = ncdf_varinq(cdfid, xvarid)
130    dimidx = xinq.dim[0]
131    IF xinq.ndims GE 2 THEN ncdf_diminq, cdfid, xinq.dim[1], blabla, jpjfromx
132  ENDELSE
133;
134  IF arg_present(xaxis) THEN BEGIN
135    ncdf_diminq, cdfid, dimidx, blabla, jpifromx
136; should we read or compute the xaxis?
137    IF keyword_set(xyindex) OR xvarid EQ -1 THEN BEGIN
138      xaxis = keyword_set(start1) + findgen(jpifromx)
139      xyindex = 1
140    ENDIF ELSE BEGIN
141; read the xaxis
142      ncdf_varget, cdfid, xvarid, xaxis
143; make sure of the shape of xaxis
144      IF n_elements(jpjfromx) NE 0 THEN xaxis = reform(xaxis, jpifromx, jpjfromx, /over)
145    ENDELSE
146  ENDIF
147
148;----------------------------------------------------------
149; find the yaxis
150; try to get the variable that contains the yaxis
151  if keyword_set(yaxisname) then yaxisname = strlowcase(yaxisname) ELSE yaxisname = 'y'
152  yvarid = (where(namevar EQ yaxisname OR namevar EQ 'latitude' $
153                  OR namevar EQ 'nav_lat' OR namevar EQ 'lat' $
154                  OR namevar EQ 'lat_rho' OR namevar EQ 'nblatitudes'))[0]
155  yvarid = yvarid[0]
156; no yaxis variable found, we will build a fake yaxis based on the size of the y dimension
157; -> we must find the y dimension
158  if yvarid EQ -1 then begin
159    dummy = report(['yaxis variable was not found within the default names:' $
160                    , '''latitude'', ''nav_lat'', ''lat'', ''lat_rho'', ''nblatitudes''' $
161                    , ' we use a fake yaxis based on y dimension size (or use YAXISNAME keyword)'], /simple)
162; try to get the dimension corresponding to y
163; roms file?
164    dimidy = where(namedim EQ 'eta_rho' OR namedim EQ 'eta_u' OR namedim EQ 'eta_v' OR namedim EQ 'eta_psi')
165    IF dimidy[0] EQ -1 THEN BEGIN
166; we are looking for a y dimension with a name matching one of the following regular expression:
167      if keyword_set(ydimname) then testname = strlowcase(ydimname) $
168      ELSE testname = ['latitude', 'lat', 'y', 'latitude*', 'lat*', 'y*', 'eta_*', '*latitude*', '*lat*', '*y*']
169      cnt = -1
170      ii = 0
171      WHILE cnt NE 1 AND ii LT n_elements(testname) DO BEGIN
172        dimidy = where(strmatch(namedim, testname[ii]) EQ 1, cnt)
173        ii = ii+1
174      ENDWHILE
175      CASE cnt OF
176        0:begin
177          dummy = report(['none of the dimensions name matches one of the following regular expression:' $
178                          , '''latitude'', ''lat'', ''y'', ''latitude*'', ''lat*'', ''y*'', ''eta_*'', ''*latitude*'', ''*lat*'', ''*y*''' $
179                          , ' => we cannot find the y dimension, use YDIMNAME keyword'], /simple)
180          stop
181        END
182        1:dimidy = dimidy[0]
183        ELSE:begin
184          dummy = report(['several (and not one unique) dimensions name matches the following regular expression:' $
185                          , '''latitude'', ''lat'', ''y'', ''latitude*'', ''lat*'', ''y*'', ''eta_*'', ''*latitude*'', ''*lat*'', ''*y*''' $
186                          , ' => we cannot find the y dimension, use YDIMNAME keyword'], /simple)
187          stop
188        ENDELSE
189      ENDCASE
190    ENDIF
191  ENDIF ELSE BEGIN
192    yinq = ncdf_varinq(cdfid, yvarid)
193    IF yinq.ndims GE 2 THEN BEGIN
194      ncdf_diminq, cdfid, yinq.dim[0], blabla, jpifromy
195      IF jpifromy NE jpifromx THEN BEGIN
196        dummy = report('x/y axes do not have the same x dimension...')
197        stop
198      ENDIF
199      dimidy = yinq.dim[1]
200    ENDIF ELSE dimidy = yinq.dim[0]
201  ENDELSE
202;
203  IF arg_present(yaxis) THEN BEGIN
204    ncdf_diminq, cdfid, dimidy, blabla, jpjfromy
205    IF n_elements(jpjfromx) NE 0 THEN BEGIN
206      IF jpjfromy NE jpjfromx THEN BEGIN
207        dummy = report(' x/y axes do not have the same y dimension...')
208        stop
209      ENDIF
210    ENDIF
211; should we read or compute the xaxis?
212    IF keyword_set(xyindex) OR yvarid EQ -1 THEN BEGIN
213      yaxis = keyword_set(start1) + findgen(jpjfromy)
214    ENDIF ELSE BEGIN
215; read the yaxis
216      ncdf_varget, cdfid, yvarid, yaxis
217; make sure of the shape of xaxis
218      IF n_elements(jpifromy) NE 0 THEN yaxis = reform(yaxis, jpifromy, jpjfromy, /over)
219    ENDELSE
220  ENDIF
221
222  return
223END
Note: See TracBrowser for help on using the repository browser.