source: trunk/SRC/ReadWrite/ncdf_getaxis.pro

Last change on this file was 495, checked in by pinsard, 10 years ago

fix thanks to coding rules; typo; dupe empty lines; trailing blanks

  • Property svn:keywords set to Id
File size: 13.9 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 XMINMESH {default=0L}{type=scalar}
54; Define common (cm_4mesh) variables ixminmesh used to define the localization
55; of the first point of the grid along the x direction in a zoom of the original grid
56;
57; @keyword YMINMESH {default=0L}{type=scalar}
58; Define common (cm_4mesh) variables iyminmesh used to define the localization
59; of the first point of the grid along the y direction in a zoom of the original grid
60;
61; @keyword XMAXMESH {default=jpiglo-1}{type=scalar}
62; Define common (cm_4mesh) variables ixmaxmesh used to define the localization
63; of the last point of the grid along the x direction in a zoom of the original grid
64; Note that if XMAXMESH < 0 then ixmaxmesh is defined as ixmaxmesh = jpiglo -1 + xmaxmesh
65;
66; @keyword YMAXMESH {default=jpjglo-1}{type=scalar}
67; Define common (cm_4mesh) variables iymaxmesh used to define the localization
68; of the last point of the grid along the y direction in a zoom of the original grid
69; Note that if YMAXMESH < 0 then iymaxmesh is defined as iymaxmesh = jpjglo -1 + ymaxmesh
70;
71; @keyword _EXTRA
72; Used to be able to call ncdf_getaxis with _extra
73;
74; @history
75; March 2007: Sebastien Masson (smasson\@locean-ipsl.upmc.fr)
76;
77; @version
78; $Id$
79;-
80PRO ncdf_getaxis, fileid, dimidx, dimidy, xaxis, yaxis $
81                  , XAXISNAME = xaxisname, YAXISNAME = yaxisname $
82                  , XDIMNAME = xdimname, YDIMNAME = ydimname $
83                  , XYINDEX = xyindex, START1 = start1 $
84                  , XMINMESH = xminmesh, XMAXMESH = xmaxmesh $
85                  , YMINMESH = yminmesh, YMAXMESH = ymaxmesh $
86                  , ROMSGRID = romsgrid, _EXTRA = ex
87;
88  compile_opt idl2, strictarrsubs
89;
90@cm_general                   ;   needed for iodir
91@cm_4mesh
92;
93  IF n_elements(xminmesh) NE 0 THEN ixminmesh = 0L > long(xminmesh[0]) ELSE ixminmesh  = 0l
94  IF n_elements(yminmesh) NE 0 THEN iyminmesh = 0L > long(yminmesh[0]) ELSE iyminmesh  = 0l
95  IF n_elements(zminmesh) NE 0 THEN izminmesh = 0L > long(zminmesh[0]) ELSE izminmesh  = 0l
96; should we open the file?
97  IF size(fileid, /type) EQ 7 THEN $
98     cdfid = ncdf_open(isafile(fileid, title = 'which file must be open by ncdf_getaxis?', IODIRECTORY = iodir, _extra = ex)) $
99  ELSE cdfid = fileid
100; what is inside the file
101  inside = ncdf_inquire(cdfid)
102;------------------------------------------------------------
103; name of all dimensions
104  namedim = strarr(inside.ndims)
105  for dimiq = 0, inside.ndims-1 do begin
106    ncdf_diminq, cdfid, dimiq, tmpname, value
107    namedim[dimiq] = strlowcase(tmpname)
108  ENDFOR
109;----------------------------------------------------------
110; name of the variables
111  namevar = strarr(inside.nvars)
112  for varid = 0, inside.nvars-1 do begin
113    invar = ncdf_varinq(cdfid, varid)
114    namevar[varid] = strlowcase(invar.name)
115  ENDFOR
116;----------------------------------------------------------
117; find the xaxis
118; try to get the variable that contains the xaxis
119  if keyword_set(xaxisname) then xaxisname = strlowcase(xaxisname) ELSE xaxisname = 'x'
120  xvarid = (where(namevar EQ xaxisname OR namevar EQ 'longitude' $
121                  OR namevar EQ 'nav_lon' OR namevar EQ 'lon' $
122                  OR namevar EQ 'lon_rho' OR namevar EQ 'nblongitudes'))[0]
123; no xaxis variable found, we will build a fake xaxis based on the size of the x dimension
124; -> we must find the x dimension
125  IF xvarid EQ -1 THEN BEGIN
126    dummy = report(['xaxis variable was not found within the default names:' $
127                    , '''x'', ''longitude'', ''nav_lon'', ''lon'', ''lon_rho'', ''nblongitudes''' $
128                    , ' we use a fake xaxis based on x dimension size (or use XAXISNAME keyword)'], /simple)
129    xaxisname = 'Not Found'
130; try to get the dimension corresponding to x
131; roms file?
132    dimidx = where(namedim EQ 'xi_rho' OR namedim EQ 'xi_u' OR namedim EQ 'xi_v' OR namedim EQ 'xi_psi')
133    IF dimidx[0] EQ -1 THEN BEGIN
134; we are looking for a x dimension with a name matching one of the following regular expression:
135      if keyword_set(xdimname) then testname = strlowcase(xdimname) $
136      ELSE testname = ['longitude', 'lon', 'x', 'longitude*', 'lon*', 'x*', '*longitude*', '*lon*', '*x*']
137      cnt = -1
138      ii = 0
139      WHILE cnt NE 1 AND ii LT n_elements(testname) DO BEGIN
140        dimidx = where(strmatch(namedim, testname[ii]) EQ 1, cnt)
141        ii = ii+1
142      ENDWHILE
143      CASE cnt OF
144        0:begin
145          dummy = report(['none of the dimensions name matches one of the following regular expression:' $
146                          , '''longitude'', ''lon'', ''x'', ''longitude*'', ''lon*'', ''x*'', ''*longitude*'', ''*lon*'', ''*x*''' $
147                          , ' => we cannot find the x dimension, use XDIMNAME keyword'], /simple)
148          stop
149        END
150        1:dimidx = dimidx[0]
151        ELSE:begin
152          dummy = report(['several (and not one unique) dimensions name matches the following regular expression:' $
153                          , '''longitude'', ''lon'', ''x'', ''longitude*'', ''lon*'', ''x*'', ''*longitude*'', ''*lon*'', ''*x*''' $
154                          , ' => we cannot find the x dimension, use XDIMNAME keyword'], /simple)
155          stop
156        ENDELSE
157      ENDCASE
158      romsgrid = 0b
159    ENDIF ELSE romsgrid = 1b
160  ENDIF ELSE BEGIN
161    romsgrid = strmid(namevar[xvarid], 0, 4) EQ 'lon_'
162    xinq = ncdf_varinq(cdfid, xvarid)
163    xaxisname = xinq.name
164    dimidx = xinq.dim[0]
165    xoffset = lonarr(xinq.ndims)
166    xcount = lonarr(xinq.ndims)
167    FOR i = 0, xinq.ndims -1 DO BEGIN
168      ncdf_diminq, cdfid, xinq.dim[i], blabla, tmpsz
169      xcount[i] = tmpsz
170    ENDFOR
171    jpiglo = xcount[0]
172    IF n_elements(xmaxmesh) NE 0 THEN BEGIN
173      IF xmaxmesh GE 0 THEN ixmaxmesh = long(xmaxmesh[0]) ELSE ixmaxmesh = xcount[0] -1 + long(xmaxmesh[0])
174    ENDIF ELSE ixmaxmesh = xcount[0] - 1L
175    ixmaxmesh  = 0 > ixmaxmesh < (xcount[0] - 1L) ; make sure ixmaxmesh is not too big
176    xoffset[0] = ixminmesh
177    xcount[0] = ixmaxmesh - ixminmesh + 1L
178    jpifromx = xcount[0]
179    jpi = jpifromx
180    IF xinq.ndims GE 2 THEN BEGIN
181      jpjglo = xcount[1]
182      IF n_elements(ymaxmesh) NE 0 THEN BEGIN
183        IF ymaxmesh GE 0 THEN iymaxmesh = long(ymaxmesh[0]) ELSE iymaxmesh = xcount[1] - 1L + long(ymaxmesh[0])
184      ENDIF ELSE iymaxmesh = xcount[1] - 1L
185      iymaxmesh  = 0 > iymaxmesh < (xcount[1] - 1L) ; make sure ixmaxmesh is not too big
186      xoffset[1] = iyminmesh
187      xcount[1] = iymaxmesh - iyminmesh + 1L
188      jpjfromx = xcount[1]
189      jpj = jpjfromx
190    ENDIF
191  ENDELSE
192  ncdf_diminq, cdfid, dimidx, xdimname, tmp
193  IF n_elements(jpifromx) EQ 0 THEN jpifromx = tmp
194;
195  IF arg_present(xaxis) THEN BEGIN
196; should we read or compute the xaxis?
197    IF keyword_set(xyindex) OR xvarid EQ -1 THEN BEGIN
198      xaxis = keyword_set(start1) + findgen(jpifromx)
199      xyindex = 1
200    ENDIF ELSE BEGIN
201; read the xaxis
202      ncdf_varget, cdfid, xvarid, xaxis, offset = xoffset, count = xcount
203      ncdf_getatt, cdfid, xvarid, ADD_OFFSET = add_offset, SCALE_FACTOR = scale_factor
204      IF scale_factor NE 1. THEN xaxis = temporary(xaxis) * scale_factor
205      IF add_offset NE 0. THEN xaxis = temporary(xaxis) + add_offset
206; make sure of the shape of xaxis
207      IF n_elements(jpjfromx) NE 0 THEN xaxis = reform(xaxis, jpifromx, jpjfromx, /over)
208    ENDELSE
209  ENDIF
210
211;----------------------------------------------------------
212; find the yaxis
213; try to get the variable that contains the yaxis
214  if keyword_set(yaxisname) then yaxisname = strlowcase(yaxisname) ELSE yaxisname = 'y'
215  yvarid = (where(namevar EQ yaxisname OR namevar EQ 'latitude' $
216                  OR namevar EQ 'nav_lat' OR namevar EQ 'lat' $
217                  OR namevar EQ 'lat_rho' OR namevar EQ 'nblatitudes'))[0]
218  yvarid = yvarid[0]
219; no yaxis variable found, we will build a fake yaxis based on the size of the y dimension
220; -> we must find the y dimension
221  if yvarid EQ -1 then begin
222    dummy = report(['yaxis variable was not found within the default names:' $
223                    , '''y'', ''latitude'', ''nav_lat'', ''lat'', ''lat_rho'', ''nblatitudes''' $
224                    , ' we use a fake yaxis based on y dimension size (or use YAXISNAME keyword)'], /simple)
225    yaxisname = 'Not Found'
226; try to get the dimension corresponding to y
227; roms file?
228    dimidy = where(namedim EQ 'eta_rho' OR namedim EQ 'eta_u' OR namedim EQ 'eta_v' OR namedim EQ 'eta_psi')
229    IF dimidy[0] EQ -1 THEN BEGIN
230; we are looking for a y dimension with a name matching one of the following regular expression:
231      if keyword_set(ydimname) then testname = strlowcase(ydimname) $
232      ELSE testname = ['latitude', 'lat', 'y', 'latitude*', 'lat*', 'y*', 'eta_*', '*latitude*', '*lat*', '*y*']
233      cnt = -1
234      ii = 0
235      WHILE cnt NE 1 AND ii LT n_elements(testname) DO BEGIN
236        dimidy = where(strmatch(namedim, testname[ii]) EQ 1, cnt)
237        ii = ii+1
238      ENDWHILE
239      CASE cnt OF
240        0:begin
241          dummy = report(['none of the dimensions name matches one of the following regular expression:' $
242                          , '''latitude'', ''lat'', ''y'', ''latitude*'', ''lat*'', ''y*'', ''eta_*'', ''*latitude*'', ''*lat*'', ''*y*''' $
243                          , ' => we cannot find the y dimension, use YDIMNAME keyword'], /simple)
244          stop
245        END
246        1:dimidy = dimidy[0]
247        ELSE:begin
248          dummy = report(['several (and not one unique) dimensions name matches the following regular expression:' $
249                          , '''latitude'', ''lat'', ''y'', ''latitude*'', ''lat*'', ''y*'', ''eta_*'', ''*latitude*'', ''*lat*'', ''*y*''' $
250                          , ' => we cannot find the y dimension, use YDIMNAME keyword'], /simple)
251          stop
252        ENDELSE
253      ENDCASE
254    ENDIF
255  ENDIF ELSE BEGIN
256    yinq = ncdf_varinq(cdfid, yvarid)
257    yaxisname = yinq.name
258    IF yinq.ndims GE 2 THEN dimidy = yinq.dim[1] ELSE dimidy = yinq.dim[0]
259    yoffset = lonarr(yinq.ndims)
260    ycount = lonarr(yinq.ndims)
261    FOR i = 0, yinq.ndims -1 DO BEGIN
262      ncdf_diminq, cdfid, yinq.dim[i], blabla, tmpsz
263      ycount[i] = tmpsz
264    ENDFOR
265    idy = yinq.ndims GE 2
266    jpjglo = ycount[idy]
267    IF n_elements(ymaxmesh) NE 0 THEN BEGIN
268      IF ymaxmesh GE 0 THEN iymaxmesh = long(ymaxmesh[0]) ELSE iymaxmesh = ycount[idy] - 1L + long(ymaxmesh[0])
269    ENDIF ELSE iymaxmesh = ycount[idy] - 1L
270    iymaxmesh  = iymaxmesh < (ycount[idy] - 1L) ; make sure ixmaxmesh is not too big
271    yoffset[idy] = iyminmesh
272    ycount[idy] = iymaxmesh - iyminmesh + 1L
273    jpjfromy = ycount[idy]
274    jpj = jpjfromy
275    IF yinq.ndims GE 2 THEN BEGIN
276      jpiglo = ycount[0]
277      IF n_elements(xmaxmesh) NE 0 THEN BEGIN
278        IF xmaxmesh GE 0 THEN ixmaxmesh = long(xmaxmesh[0]) ELSE ixmaxmesh = ycount[0] - 1L + long(xmaxmesh[0])
279      ENDIF ELSE ixmaxmesh = ycount[0] - 1L
280      ixmaxmesh  = ixmaxmesh < (ycount[0] - 1L) ; make sure ixmaxmesh is not too big
281      yoffset[0] = ixminmesh
282      ycount[0] = ixmaxmesh - ixminmesh + 1L
283      jpifromy = xcount[0]
284      jpi = jpifromy
285    ENDIF
286  ENDELSE
287  ncdf_diminq, cdfid, dimidy, ydimname, tmp
288  IF n_elements(jpjfromy) EQ 0 THEN jpjfromy = tmp
289;
290  IF arg_present(yaxis) THEN BEGIN
291    IF n_elements(jpifromy) NE 0 THEN BEGIN
292      IF n_elements(jpifromx) EQ 0 THEN ncdf_diminq, cdfid, dimidx, blabla, jpifromx
293      IF jpifromy NE jpifromx THEN BEGIN
294        dummy = report('x/y axes do not have the same x dimension...')
295        stop
296      ENDIF
297    ENDIF
298    IF n_elements(jpjfromx) NE 0 THEN BEGIN
299      IF n_elements(jpjfromy) EQ 0 THEN ncdf_diminq, cdfid, dimidy, blabla, jpjfromy
300      IF jpjfromy NE jpjfromx THEN BEGIN
301        dummy = report(' x/y axes do not have the same y dimension...')
302        stop
303      ENDIF
304    ENDIF
305; should we read or compute the xaxis?
306    IF keyword_set(xyindex) OR yvarid EQ -1 THEN BEGIN
307      yaxis = keyword_set(start1) + findgen(jpjfromy)
308    ENDIF ELSE BEGIN
309; read the yaxis
310      ncdf_varget, cdfid, yvarid, yaxis, offset = yoffset, count = ycount
311      ncdf_getatt, cdfid, yvarid, ADD_OFFSET = add_offset, SCALE_FACTOR = scale_factor
312      IF scale_factor NE 1. THEN yaxis = temporary(yaxis) * scale_factor
313      IF add_offset NE 0. THEN yaxis = temporary(yaxis) + add_offset
314; make sure of the shape of xaxis
315      IF n_elements(jpifromy) NE 0 THEN yaxis = reform(yaxis, jpifromy, jpjfromy, /over)
316    ENDELSE
317  ENDIF
318
319  IF size(fileid, /type) EQ 7 THEN ncdf_close, cdfid
320
321  return
322END
Note: See TracBrowser for help on using the repository browser.