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

Last change on this file since 385 was 385, checked in by smasson, 16 years ago

take into account add_offset and scale_factor in ncdf_getaxis

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