[221] | 1 | ;+ |
---|
[232] | 2 | ; |
---|
[221] | 3 | ; @file_comments |
---|
| 4 | ; get the x/y dimension Id and x/y axes from a netcdf file |
---|
| 5 | ; |
---|
| 6 | ; @categories |
---|
[271] | 7 | ; Read NetCDF file |
---|
[221] | 8 | ; |
---|
[325] | 9 | ; @param fileid {in}{required}{type=scalar} |
---|
[221] | 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} |
---|
[224] | 28 | ; Index the axis from 1 instead of 0 when using /xyindex |
---|
[221] | 29 | ; |
---|
[224] | 30 | ; @keyword XDIMNAME {default='longitude', 'lon', 'x', 'longitude*', 'lon*', 'x*', '*longitude*', '*lon*' or '*x*'}{type=scalar string} |
---|
[271] | 31 | ; A string giving the name of the x dimension or/and a named variable |
---|
| 32 | ; in which x dimension name is returned. |
---|
[224] | 33 | ; |
---|
| 34 | ; @keyword YDIMNAME {default='latitude', 'lat', 'y', 'latitude*', 'lat*', 'y*', 'eta_*', '*latitude*', '*lat*', '*y*'}{type=scalar string} |
---|
[271] | 35 | ; A string giving the name of the y dimension or/and a named variable |
---|
| 36 | ; in which y dimension name is returned. |
---|
[224] | 37 | ; |
---|
[221] | 38 | ; @keyword XAXISNAME {default='x', 'longitude', 'nav_lon', 'lon', 'lon_rho' or 'NbLongitudes'}{type=scalar string} |
---|
[224] | 39 | ; A string giving the name of the variable in the file |
---|
[271] | 40 | ; that contains the x axis or/and a named variable |
---|
| 41 | ; in which this variable name is returned. |
---|
[224] | 42 | ; |
---|
[221] | 43 | ; @keyword YAXISNAME {default='y', 'latitude', 'nav_lat','lat', 'lat_rho' or 'NbLatitudes'}{type=scalar string} |
---|
[224] | 44 | ; A string giving the name of the variable in the file |
---|
[271] | 45 | ; that contains the y axis or/and a named variable |
---|
| 46 | ; in which this variable name is returned. |
---|
[221] | 47 | ; |
---|
| 48 | ; @keyword XYINDEX {default=0}{type=scalar: 0 or 1} |
---|
| 49 | ; To define the x/y axis with index instead of using |
---|
[224] | 50 | ; the values contained in X/YAXISNAME. |
---|
| 51 | ; x/yaxis = keyword_set(start1) + findgen(jpi/jpj) |
---|
[221] | 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 | ;- |
---|
[271] | 63 | PRO ncdf_getaxis, fileid, dimidx, dimidy, xaxis, yaxis $ |
---|
[327] | 64 | , XAXISNAME=xaxisname, YAXISNAME=yaxisname $ |
---|
| 65 | , XDIMNAME=xdimname, YDIMNAME=ydimname $ |
---|
| 66 | , XYINDEX=xyindex, START1=start1 $ |
---|
| 67 | , ROMSGRID=romsgrid, _EXTRA=ex |
---|
[221] | 68 | ; |
---|
| 69 | compile_opt idl2, strictarrsubs |
---|
| 70 | ; |
---|
[271] | 71 | ; should we open the file? |
---|
| 72 | IF size(fileid, /type) EQ 7 THEN cdfid = ncdf_open(fileid) ELSE cdfid = fileid |
---|
[221] | 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 |
---|
[224] | 79 | ncdf_diminq, cdfid, dimiq, tmpname, value |
---|
[221] | 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] |
---|
[224] | 96 | ; no xaxis variable found, we will build a fake xaxis based on the size of the x dimension |
---|
[221] | 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) |
---|
[271] | 102 | xaxisname = 'Not Found' |
---|
[221] | 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 |
---|
[271] | 131 | romsgrid = 0b |
---|
| 132 | ENDIF ELSE romsgrid = 1b |
---|
[224] | 133 | ENDIF ELSE BEGIN |
---|
[221] | 134 | romsgrid = strmid(namevar[xvarid], 0, 4) EQ 'lon_' |
---|
| 135 | xinq = ncdf_varinq(cdfid, xvarid) |
---|
[271] | 136 | xaxisname = xinq.name |
---|
[221] | 137 | dimidx = xinq.dim[0] |
---|
| 138 | IF xinq.ndims GE 2 THEN ncdf_diminq, cdfid, xinq.dim[1], blabla, jpjfromx |
---|
| 139 | ENDELSE |
---|
[271] | 140 | IF arg_present(xdimname) THEN ncdf_diminq, cdfid, dimidx, xdimname, jpifromx |
---|
[221] | 141 | ; |
---|
| 142 | IF arg_present(xaxis) THEN BEGIN |
---|
| 143 | ncdf_diminq, cdfid, dimidx, blabla, jpifromx |
---|
| 144 | ; should we read or compute the xaxis? |
---|
[224] | 145 | IF keyword_set(xyindex) OR xvarid EQ -1 THEN BEGIN |
---|
[221] | 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 | ; make sure of the shape of xaxis |
---|
| 152 | IF n_elements(jpjfromx) NE 0 THEN xaxis = reform(xaxis, jpifromx, jpjfromx, /over) |
---|
[224] | 153 | ENDELSE |
---|
[221] | 154 | ENDIF |
---|
| 155 | |
---|
| 156 | ;---------------------------------------------------------- |
---|
| 157 | ; find the yaxis |
---|
| 158 | ; try to get the variable that contains the yaxis |
---|
| 159 | if keyword_set(yaxisname) then yaxisname = strlowcase(yaxisname) ELSE yaxisname = 'y' |
---|
| 160 | yvarid = (where(namevar EQ yaxisname OR namevar EQ 'latitude' $ |
---|
| 161 | OR namevar EQ 'nav_lat' OR namevar EQ 'lat' $ |
---|
| 162 | OR namevar EQ 'lat_rho' OR namevar EQ 'nblatitudes'))[0] |
---|
| 163 | yvarid = yvarid[0] |
---|
[224] | 164 | ; no yaxis variable found, we will build a fake yaxis based on the size of the y dimension |
---|
[221] | 165 | ; -> we must find the y dimension |
---|
| 166 | if yvarid EQ -1 then begin |
---|
| 167 | dummy = report(['yaxis variable was not found within the default names:' $ |
---|
| 168 | , '''latitude'', ''nav_lat'', ''lat'', ''lat_rho'', ''nblatitudes''' $ |
---|
| 169 | , ' we use a fake yaxis based on y dimension size (or use YAXISNAME keyword)'], /simple) |
---|
[271] | 170 | yaxisname = 'Not Found' |
---|
[221] | 171 | ; try to get the dimension corresponding to y |
---|
| 172 | ; roms file? |
---|
| 173 | dimidy = where(namedim EQ 'eta_rho' OR namedim EQ 'eta_u' OR namedim EQ 'eta_v' OR namedim EQ 'eta_psi') |
---|
| 174 | IF dimidy[0] EQ -1 THEN BEGIN |
---|
| 175 | ; we are looking for a y dimension with a name matching one of the following regular expression: |
---|
| 176 | if keyword_set(ydimname) then testname = strlowcase(ydimname) $ |
---|
| 177 | ELSE testname = ['latitude', 'lat', 'y', 'latitude*', 'lat*', 'y*', 'eta_*', '*latitude*', '*lat*', '*y*'] |
---|
| 178 | cnt = -1 |
---|
| 179 | ii = 0 |
---|
| 180 | WHILE cnt NE 1 AND ii LT n_elements(testname) DO BEGIN |
---|
| 181 | dimidy = where(strmatch(namedim, testname[ii]) EQ 1, cnt) |
---|
| 182 | ii = ii+1 |
---|
| 183 | ENDWHILE |
---|
| 184 | CASE cnt OF |
---|
| 185 | 0:begin |
---|
| 186 | dummy = report(['none of the dimensions name matches one of the following regular expression:' $ |
---|
| 187 | , '''latitude'', ''lat'', ''y'', ''latitude*'', ''lat*'', ''y*'', ''eta_*'', ''*latitude*'', ''*lat*'', ''*y*''' $ |
---|
| 188 | , ' => we cannot find the y dimension, use YDIMNAME keyword'], /simple) |
---|
| 189 | stop |
---|
| 190 | END |
---|
| 191 | 1:dimidy = dimidy[0] |
---|
| 192 | ELSE:begin |
---|
| 193 | dummy = report(['several (and not one unique) dimensions name matches 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 | ENDELSE |
---|
| 198 | ENDCASE |
---|
| 199 | ENDIF |
---|
[224] | 200 | ENDIF ELSE BEGIN |
---|
[221] | 201 | yinq = ncdf_varinq(cdfid, yvarid) |
---|
[271] | 202 | yaxisname = yinq.name |
---|
[221] | 203 | IF yinq.ndims GE 2 THEN BEGIN |
---|
| 204 | ncdf_diminq, cdfid, yinq.dim[0], blabla, jpifromy |
---|
| 205 | dimidy = yinq.dim[1] |
---|
| 206 | ENDIF ELSE dimidy = yinq.dim[0] |
---|
| 207 | ENDELSE |
---|
[271] | 208 | IF arg_present(ydimname) THEN ncdf_diminq, cdfid, dimidy, ydimname, jpjfromy |
---|
[221] | 209 | ; |
---|
| 210 | IF arg_present(yaxis) THEN BEGIN |
---|
[225] | 211 | IF n_elements(jpifromy) NE 0 THEN BEGIN |
---|
| 212 | IF jpifromy NE jpifromx THEN BEGIN |
---|
| 213 | dummy = report('x/y axes do not have the same x dimension...') |
---|
| 214 | stop |
---|
| 215 | ENDIF |
---|
| 216 | ENDIF |
---|
[221] | 217 | ncdf_diminq, cdfid, dimidy, blabla, jpjfromy |
---|
| 218 | IF n_elements(jpjfromx) NE 0 THEN BEGIN |
---|
[224] | 219 | IF jpjfromy NE jpjfromx THEN BEGIN |
---|
[221] | 220 | dummy = report(' x/y axes do not have the same y dimension...') |
---|
| 221 | stop |
---|
[224] | 222 | ENDIF |
---|
[221] | 223 | ENDIF |
---|
| 224 | ; should we read or compute the xaxis? |
---|
[224] | 225 | IF keyword_set(xyindex) OR yvarid EQ -1 THEN BEGIN |
---|
[221] | 226 | yaxis = keyword_set(start1) + findgen(jpjfromy) |
---|
| 227 | ENDIF ELSE BEGIN |
---|
| 228 | ; read the yaxis |
---|
| 229 | ncdf_varget, cdfid, yvarid, yaxis |
---|
| 230 | ; make sure of the shape of xaxis |
---|
| 231 | IF n_elements(jpifromy) NE 0 THEN yaxis = reform(yaxis, jpifromy, jpjfromy, /over) |
---|
| 232 | ENDELSE |
---|
| 233 | ENDIF |
---|
| 234 | |
---|
[271] | 235 | IF size(fileid, /type) EQ 7 THEN ncdf_close, cdfid |
---|
| 236 | |
---|
[221] | 237 | return |
---|
| 238 | END |
---|