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 | ;- |
---|
63 | PRO 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 | @cm_general ; needed for iodir |
---|
72 | ; |
---|
73 | ; should we open the file? |
---|
74 | IF size(fileid, /type) EQ 7 THEN $ |
---|
75 | cdfid = ncdf_open(isafile(fileid, title = 'which file must be open by ncdf_getaxis?', IODIRECTORY = iodir, _extra = ex)) $ |
---|
76 | ELSE cdfid = fileid |
---|
77 | ; what is inside the file |
---|
78 | inside = ncdf_inquire(cdfid) |
---|
79 | ;------------------------------------------------------------ |
---|
80 | ; name of all dimensions |
---|
81 | namedim = strarr(inside.ndims) |
---|
82 | for dimiq = 0, inside.ndims-1 do begin |
---|
83 | ncdf_diminq, cdfid, dimiq, tmpname, value |
---|
84 | namedim[dimiq] = strlowcase(tmpname) |
---|
85 | ENDFOR |
---|
86 | ;---------------------------------------------------------- |
---|
87 | ; name of the variables |
---|
88 | namevar = strarr(inside.nvars) |
---|
89 | for varid = 0, inside.nvars-1 do begin |
---|
90 | invar = ncdf_varinq(cdfid, varid) |
---|
91 | namevar[varid] = strlowcase(invar.name) |
---|
92 | ENDFOR |
---|
93 | ;---------------------------------------------------------- |
---|
94 | ; find the xaxis |
---|
95 | ; try to get the variable that contains the xaxis |
---|
96 | if keyword_set(xaxisname) then xaxisname = strlowcase(xaxisname) ELSE xaxisname = 'x' |
---|
97 | xvarid = (where(namevar EQ xaxisname OR namevar EQ 'longitude' $ |
---|
98 | OR namevar EQ 'nav_lon' OR namevar EQ 'lon' $ |
---|
99 | OR namevar EQ 'lon_rho' OR namevar EQ 'nblongitudes'))[0] |
---|
100 | ; no xaxis variable found, we will build a fake xaxis based on the size of the x dimension |
---|
101 | ; -> we must find the x dimension |
---|
102 | IF xvarid EQ -1 THEN BEGIN |
---|
103 | dummy = report(['xaxis variable was not found within the default names:' $ |
---|
104 | , '''longitude'', ''nav_lon'', ''lon'', ''lon_rho'', ''nblongitudes''' $ |
---|
105 | , ' we use a fake xaxis based on x dimension size (or use XAXISNAME keyword)'], /simple) |
---|
106 | xaxisname = 'Not Found' |
---|
107 | ; try to get the dimension corresponding to x |
---|
108 | ; roms file? |
---|
109 | dimidx = where(namedim EQ 'xi_rho' OR namedim EQ 'xi_u' OR namedim EQ 'xi_v' OR namedim EQ 'xi_psi') |
---|
110 | IF dimidx[0] EQ -1 THEN BEGIN |
---|
111 | ; we are looking for a x dimension with a name matching one of the following regular expression: |
---|
112 | if keyword_set(xdimname) then testname = strlowcase(xdimname) $ |
---|
113 | ELSE testname = ['longitude', 'lon', 'x', 'longitude*', 'lon*', 'x*', '*longitude*', '*lon*', '*x*'] |
---|
114 | cnt = -1 |
---|
115 | ii = 0 |
---|
116 | WHILE cnt NE 1 AND ii LT n_elements(testname) DO BEGIN |
---|
117 | dimidx = where(strmatch(namedim, testname[ii]) EQ 1, cnt) |
---|
118 | ii = ii+1 |
---|
119 | ENDWHILE |
---|
120 | CASE cnt OF |
---|
121 | 0:begin |
---|
122 | dummy = report(['none of the dimensions name matches one of the following regular expression:' $ |
---|
123 | , '''longitude'', ''lon'', ''x'', ''longitude*'', ''lon*'', ''x*'', ''*longitude*'', ''*lon*'', ''*x*''' $ |
---|
124 | , ' => we cannot find the x dimension, use XDIMNAME keyword'], /simple) |
---|
125 | stop |
---|
126 | END |
---|
127 | 1:dimidx = dimidx[0] |
---|
128 | ELSE:begin |
---|
129 | dummy = report(['several (and not one unique) dimensions name matches the following regular expression:' $ |
---|
130 | , '''longitude'', ''lon'', ''x'', ''longitude*'', ''lon*'', ''x*'', ''*longitude*'', ''*lon*'', ''*x*''' $ |
---|
131 | , ' => we cannot find the x dimension, use XDIMNAME keyword'], /simple) |
---|
132 | stop |
---|
133 | ENDELSE |
---|
134 | ENDCASE |
---|
135 | romsgrid = 0b |
---|
136 | ENDIF ELSE romsgrid = 1b |
---|
137 | ENDIF ELSE BEGIN |
---|
138 | romsgrid = strmid(namevar[xvarid], 0, 4) EQ 'lon_' |
---|
139 | xinq = ncdf_varinq(cdfid, xvarid) |
---|
140 | xaxisname = xinq.name |
---|
141 | dimidx = xinq.dim[0] |
---|
142 | IF xinq.ndims GE 2 THEN ncdf_diminq, cdfid, xinq.dim[1], blabla, jpjfromx |
---|
143 | ENDELSE |
---|
144 | IF arg_present(xdimname) THEN ncdf_diminq, cdfid, dimidx, xdimname, jpifromx |
---|
145 | ; |
---|
146 | IF arg_present(xaxis) THEN BEGIN |
---|
147 | ncdf_diminq, cdfid, dimidx, blabla, jpifromx |
---|
148 | ; should we read or compute the xaxis? |
---|
149 | IF keyword_set(xyindex) OR xvarid EQ -1 THEN BEGIN |
---|
150 | xaxis = keyword_set(start1) + findgen(jpifromx) |
---|
151 | xyindex = 1 |
---|
152 | ENDIF ELSE BEGIN |
---|
153 | ; read the xaxis |
---|
154 | ncdf_varget, cdfid, xvarid, xaxis |
---|
155 | ncdf_getatt, cdfid, xvarid, ADD_OFFSET = add_offset, SCALE_FACTOR = scale_factor |
---|
156 | IF scale_factor NE 1. THEN xaxis = temporary(xaxis) * scale_factor |
---|
157 | IF add_offset NE 0. THEN xaxis = temporary(xaxis) + add_offset |
---|
158 | ; make sure of the shape of xaxis |
---|
159 | IF n_elements(jpjfromx) NE 0 THEN xaxis = reform(xaxis, jpifromx, jpjfromx, /over) |
---|
160 | ENDELSE |
---|
161 | ENDIF |
---|
162 | |
---|
163 | ;---------------------------------------------------------- |
---|
164 | ; find the yaxis |
---|
165 | ; try to get the variable that contains the yaxis |
---|
166 | if keyword_set(yaxisname) then yaxisname = strlowcase(yaxisname) ELSE yaxisname = 'y' |
---|
167 | yvarid = (where(namevar EQ yaxisname OR namevar EQ 'latitude' $ |
---|
168 | OR namevar EQ 'nav_lat' OR namevar EQ 'lat' $ |
---|
169 | OR namevar EQ 'lat_rho' OR namevar EQ 'nblatitudes'))[0] |
---|
170 | yvarid = yvarid[0] |
---|
171 | ; no yaxis variable found, we will build a fake yaxis based on the size of the y dimension |
---|
172 | ; -> we must find the y dimension |
---|
173 | if yvarid EQ -1 then begin |
---|
174 | dummy = report(['yaxis variable was not found within the default names:' $ |
---|
175 | , '''latitude'', ''nav_lat'', ''lat'', ''lat_rho'', ''nblatitudes''' $ |
---|
176 | , ' we use a fake yaxis based on y dimension size (or use YAXISNAME keyword)'], /simple) |
---|
177 | yaxisname = 'Not Found' |
---|
178 | ; try to get the dimension corresponding to y |
---|
179 | ; roms file? |
---|
180 | dimidy = where(namedim EQ 'eta_rho' OR namedim EQ 'eta_u' OR namedim EQ 'eta_v' OR namedim EQ 'eta_psi') |
---|
181 | IF dimidy[0] EQ -1 THEN BEGIN |
---|
182 | ; we are looking for a y dimension with a name matching one of the following regular expression: |
---|
183 | if keyword_set(ydimname) then testname = strlowcase(ydimname) $ |
---|
184 | ELSE testname = ['latitude', 'lat', 'y', 'latitude*', 'lat*', 'y*', 'eta_*', '*latitude*', '*lat*', '*y*'] |
---|
185 | cnt = -1 |
---|
186 | ii = 0 |
---|
187 | WHILE cnt NE 1 AND ii LT n_elements(testname) DO BEGIN |
---|
188 | dimidy = where(strmatch(namedim, testname[ii]) EQ 1, cnt) |
---|
189 | ii = ii+1 |
---|
190 | ENDWHILE |
---|
191 | CASE cnt OF |
---|
192 | 0:begin |
---|
193 | dummy = report(['none of the dimensions name matches one of the following regular expression:' $ |
---|
194 | , '''latitude'', ''lat'', ''y'', ''latitude*'', ''lat*'', ''y*'', ''eta_*'', ''*latitude*'', ''*lat*'', ''*y*''' $ |
---|
195 | , ' => we cannot find the y dimension, use YDIMNAME keyword'], /simple) |
---|
196 | stop |
---|
197 | END |
---|
198 | 1:dimidy = dimidy[0] |
---|
199 | ELSE:begin |
---|
200 | dummy = report(['several (and not one unique) dimensions name matches the following regular expression:' $ |
---|
201 | , '''latitude'', ''lat'', ''y'', ''latitude*'', ''lat*'', ''y*'', ''eta_*'', ''*latitude*'', ''*lat*'', ''*y*''' $ |
---|
202 | , ' => we cannot find the y dimension, use YDIMNAME keyword'], /simple) |
---|
203 | stop |
---|
204 | ENDELSE |
---|
205 | ENDCASE |
---|
206 | ENDIF |
---|
207 | ENDIF ELSE BEGIN |
---|
208 | yinq = ncdf_varinq(cdfid, yvarid) |
---|
209 | yaxisname = yinq.name |
---|
210 | IF yinq.ndims GE 2 THEN BEGIN |
---|
211 | ncdf_diminq, cdfid, yinq.dim[0], blabla, jpifromy |
---|
212 | dimidy = yinq.dim[1] |
---|
213 | ENDIF ELSE dimidy = yinq.dim[0] |
---|
214 | ENDELSE |
---|
215 | IF arg_present(ydimname) THEN ncdf_diminq, cdfid, dimidy, ydimname, jpjfromy |
---|
216 | ; |
---|
217 | IF arg_present(yaxis) THEN BEGIN |
---|
218 | IF n_elements(jpifromy) NE 0 THEN BEGIN |
---|
219 | IF jpifromy NE jpifromx THEN BEGIN |
---|
220 | dummy = report('x/y axes do not have the same x dimension...') |
---|
221 | stop |
---|
222 | ENDIF |
---|
223 | ENDIF |
---|
224 | ncdf_diminq, cdfid, dimidy, blabla, jpjfromy |
---|
225 | IF n_elements(jpjfromx) NE 0 THEN BEGIN |
---|
226 | IF jpjfromy NE jpjfromx THEN BEGIN |
---|
227 | dummy = report(' x/y axes do not have the same y dimension...') |
---|
228 | stop |
---|
229 | ENDIF |
---|
230 | ENDIF |
---|
231 | ; should we read or compute the xaxis? |
---|
232 | IF keyword_set(xyindex) OR yvarid EQ -1 THEN BEGIN |
---|
233 | yaxis = keyword_set(start1) + findgen(jpjfromy) |
---|
234 | ENDIF ELSE BEGIN |
---|
235 | ; read the yaxis |
---|
236 | ncdf_varget, cdfid, yvarid, yaxis |
---|
237 | ncdf_getatt, cdfid, yvarid, ADD_OFFSET = add_offset, SCALE_FACTOR = scale_factor |
---|
238 | IF scale_factor NE 1. THEN yaxis = temporary(yaxis) * scale_factor |
---|
239 | IF add_offset NE 0. THEN yaxis = temporary(yaxis) + add_offset |
---|
240 | ; make sure of the shape of xaxis |
---|
241 | IF n_elements(jpifromy) NE 0 THEN yaxis = reform(yaxis, jpifromy, jpjfromy, /over) |
---|
242 | ENDELSE |
---|
243 | ENDIF |
---|
244 | |
---|
245 | IF size(fileid, /type) EQ 7 THEN ncdf_close, cdfid |
---|
246 | |
---|
247 | return |
---|
248 | END |
---|