read_netcdf_var_3d_2d.c File Reference

Read a 2D field from a 3D NetCDF variable. More...

#include <io.h>
Include dependency graph for read_netcdf_var_3d_2d.c:

Go to the source code of this file.

Functions

int read_netcdf_var_3d_2d (double **buf, info_field_struct *info_field, proj_struct *proj, char *filename, char *varname, char *dimxname, char *dimyname, char *timename, int t, int *nlon, int *nlat, int *ntime, int outinfo)
 Read a 2D field from a 3D variable in a NetCDF file, and return information in info_field_struct structure and proj_struct.

Detailed Description

Read a 2D field from a 3D NetCDF variable.

Definition in file read_netcdf_var_3d_2d.c.


Function Documentation

int read_netcdf_var_3d_2d ( double **  buf,
info_field_struct info_field,
proj_struct proj,
char *  filename,
char *  varname,
char *  dimxname,
char *  dimyname,
char *  timename,
int  t,
int *  nlon,
int *  nlat,
int *  ntime,
int  outinfo 
)

Read a 2D field from a 3D variable in a NetCDF file, and return information in info_field_struct structure and proj_struct.

Parameters:
[out] buf 2D variable
[out] info_field Information about the output variable
[out] proj Information about the horizontal projection of the output variable
[in] filename NetCDF input filename
[in] varname NetCDF variable name
[in] dimxname Longitude dimension name
[in] dimyname Latitude dimension name
[in] timename Time dimension name
[in] t Time index to retrieve
[out] nlon Longitude dimension length
[out] nlat Latitude dimension length
[out] ntime Time dimension length
[in] outinfo TRUE if we want information output, FALSE if not
Returns:
Status.

Read data variable

Definition at line 68 of file read_netcdf_var_3d_2d.c.

References alloc_error(), info_field_struct::coordinates, proj_struct::false_easting, proj_struct::false_northing, info_field_struct::fillvalue, get_attribute_str(), info_field_struct::grid_mapping, proj_struct::grid_mapping_name, handle_netcdf_error(), info_field_struct::height, proj_struct::lat0, proj_struct::latin1, proj_struct::latin2, proj_struct::latpole, proj_struct::lonc, info_field_struct::long_name, proj_struct::lonpole, MAXPATH, proj_struct::name, TRUE, and info_field_struct::units.

Referenced by output_downscaled_analog(), read_field_subdomain_period(), and read_obs_period().

00069                                                                                                                             {
00088   int istat; /* Diagnostic status */
00089 
00090   size_t dimval; /* Variable used to retrieve dimension length */
00091 
00092   int ncinid; /* NetCDF input file handle ID */
00093   int varinid; /* NetCDF variable ID */
00094   int projinid; /* Projection variable ID */
00095   nc_type vartype_main; /* Type of the variable (NC_FLOAT, NC_DOUBLE, etc.) */
00096   int varndims; /* Number of dimensions of variable */
00097   int vardimids[NC_MAX_VAR_DIMS]; /* Variable dimension ids */
00098   int timediminid; /* Time dimension ID */
00099   int londiminid; /* Longitude dimension ID */
00100   int latdiminid; /* Latitude dimension ID */
00101 
00102   size_t start[3]; /* Start position to read */
00103   size_t count[3]; /* Number of elements to read */
00104 
00105   float valf; /* Variable used to retrieve fillvalue */
00106   int vali; /* Variable used to retrieve integer values */
00107   char *tmpstr = NULL; /* Temporary string */
00108   size_t t_len; /* Length of string attribute */
00109 
00110   float *proj_latin = NULL; /* Parallel latitudes of projection */
00111   int npts = 0; /* Number of points for 1D variables */
00112   char *grid_mapping = NULL;
00113 
00114   /* Allocate memory */
00115   tmpstr = (char *) malloc(MAXPATH * sizeof(char));
00116   if (tmpstr == NULL) alloc_error(__FILE__, __LINE__);
00117 
00118   /* Read data in NetCDF file */
00119 
00120   /* Open NetCDF file for reading */
00121   if (outinfo == TRUE)
00122     printf("%s: Opening for reading NetCDF input file %s\n", __FILE__, filename);
00123   istat = nc_open(filename, NC_NOWRITE, &ncinid);  /* open for reading */
00124   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00125 
00126   /* Get dimensions length */
00127   istat = nc_inq_dimid(ncinid, timename, &timediminid);  /* get ID for time dimension */
00128   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00129   istat = nc_inq_dimlen(ncinid, timediminid, &dimval); /* get time length */
00130   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00131   *ntime = (int) dimval;
00132   /* Verify timestep provided */
00133   if (t < 0 || t > ((*ntime)-1)) {
00134     (void) free(tmpstr);
00135     istat = ncclose(ncinid);
00136     (void) fprintf(stderr, "%s: Invalid timestep provided: %d. Maximum value is %d\n", __FILE__, t, *ntime);
00137     return -1;
00138   }
00139 
00140   if (outinfo == TRUE)
00141     printf("%s: READ %s %s time=%d %d.\n", __FILE__, varname, filename, t, *ntime);
00142 
00143   istat = nc_inq_dimid(ncinid, dimyname, &latdiminid);  /* get ID for lat dimension */
00144   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00145   istat = nc_inq_dimlen(ncinid, latdiminid, &dimval); /* get lat length */
00146   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00147   *nlat = (int) dimval;
00148 
00149   istat = nc_inq_dimid(ncinid, dimxname, &londiminid);  /* get ID for lon dimension */
00150   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00151   istat = nc_inq_dimlen(ncinid, londiminid, &dimval); /* get lon length */
00152   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00153   *nlon = (int) dimval;
00154 
00155   /* Get main variable ID */
00156   istat = nc_inq_varid(ncinid, varname, &varinid);
00157   if (istat != NC_NOERR) {
00158     (void) fprintf(stderr, "%s: Error with variable %s in file %s\n", __FILE__, varname, filename);
00159     handle_netcdf_error(istat, __FILE__, __LINE__);
00160   }
00161 
00164   /* Get variable information */
00165   istat = nc_inq_var(ncinid, varinid, (char *) NULL, &vartype_main, &varndims, vardimids, (int *) NULL);
00166   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00167 
00168   /* Verify that variable is really 3D or 2D */
00169   if (varndims != 3 && varndims != 2) {
00170     (void) fprintf(stderr, "%s: Error NetCDF type and/or dimensions nlon %d nlat %d.\n", __FILE__, *nlon, *nlat);
00171     (void) free(tmpstr);
00172     istat = ncclose(ncinid);
00173     if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00174     return -1;
00175   }
00176 
00177   if (varndims == 2) {
00178     if ((*nlat) != (*nlon)) {
00179       (void) fprintf(stderr, "%s: Error NetCDF type and/or dimensions nlon %d nlat %d.\n", __FILE__, *nlon, *nlat);
00180       (void) free(tmpstr);
00181       istat = ncclose(ncinid);
00182       if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00183       return -1;
00184     }
00185     npts = *nlon;
00186     *nlat = 0;
00187   }
00188 
00189   /* If info_field si not NULL, get some information about the read variable */
00190   if (info_field != NULL) {
00191     /* Get missing value */
00192     if (vartype_main == NC_FLOAT) {
00193       istat = nc_get_att_float(ncinid, varinid, "missing_value", &valf);
00194       if (istat != NC_NOERR)
00195         info_field->fillvalue = -9999.0;
00196       else
00197         info_field->fillvalue = (double) valf;
00198     }
00199     else if (vartype_main == NC_DOUBLE) {
00200       istat = nc_get_att_double(ncinid, varinid, "missing_value", &(info_field->fillvalue));
00201       if (istat != NC_NOERR)
00202         info_field->fillvalue = -9999.0;
00203     }
00204 
00205     /* Get coordinates */
00206     istat = nc_inq_attlen(ncinid, varinid, "coordinates", &t_len);
00207     if (istat == NC_NOERR) {
00208       istat = nc_get_att_text(ncinid, varinid, "coordinates", tmpstr);
00209       if (istat == NC_NOERR) {
00210         if (tmpstr[t_len-1] != '\0')
00211           tmpstr[t_len] = '\0';
00212         info_field->coordinates = strdup(tmpstr);
00213       }
00214       else
00215         info_field->coordinates = strdup("lon lat");
00216     }
00217     else
00218       info_field->coordinates = strdup("lon lat");
00219 
00220     /* Get grid projection */
00221     istat = nc_inq_attlen(ncinid, varinid, "grid_mapping", &t_len);
00222     if (istat == NC_NOERR) {
00223       handle_netcdf_error(istat, __FILE__, __LINE__);
00224       istat = nc_get_att_text(ncinid, varinid, "grid_mapping", tmpstr);
00225       if (istat == NC_NOERR) {
00226         if (tmpstr[t_len-1] != '\0')
00227           tmpstr[t_len] = '\0';
00228         info_field->grid_mapping = strdup(tmpstr);
00229       }
00230       else
00231         info_field->grid_mapping = strdup("unknown");
00232     }
00233     else
00234       info_field->grid_mapping = strdup("unknown");
00235 
00236     /* Get units */
00237     istat = nc_inq_attlen(ncinid, varinid, "units", &t_len);
00238     if (istat == NC_NOERR) {
00239       handle_netcdf_error(istat, __FILE__, __LINE__);
00240       istat = nc_get_att_text(ncinid, varinid, "units", tmpstr);
00241       if (istat == NC_NOERR) {
00242         if (tmpstr[t_len-1] != '\0')
00243           tmpstr[t_len] = '\0';
00244         info_field->units = strdup(tmpstr);
00245       }
00246       else
00247         info_field->units = strdup("unknown");
00248     }
00249     else
00250       info_field->units = strdup("unknown");
00251 
00252     /* Get height */
00253     istat = nc_inq_attlen(ncinid, varinid, "height", &t_len);
00254     if (istat == NC_NOERR) {
00255       handle_netcdf_error(istat, __FILE__, __LINE__);
00256       istat = nc_get_att_text(ncinid, varinid, "height", tmpstr);
00257       if (istat == NC_NOERR) {
00258         if (tmpstr[t_len-1] != '\0')
00259           tmpstr[t_len] = '\0';
00260         info_field->height = strdup(tmpstr);
00261       }
00262       else
00263         info_field->height = strdup("unknown");
00264     }
00265     else
00266       info_field->height = strdup("unknown");
00267 
00268     /* Get long name */
00269     istat = nc_inq_attlen(ncinid, varinid, "long_name", &t_len);
00270     if (istat == NC_NOERR) {
00271       handle_netcdf_error(istat, __FILE__, __LINE__);
00272       istat = nc_get_att_text(ncinid, varinid, "long_name", tmpstr);
00273       if (istat == NC_NOERR) {
00274         if (tmpstr[t_len-1] != '\0')
00275           tmpstr[t_len] = '\0';
00276         info_field->long_name = strdup(tmpstr);
00277       }
00278       else
00279         info_field->long_name = strdup(varname);
00280     }
00281     else
00282       info_field->long_name = strdup(varname);
00283   }
00284 
00285   /* if proj is not NULL, retrieve informations about the horizontal projection parameters */
00286   if (proj != NULL) {
00287     if (info_field == NULL)
00288       grid_mapping = strdup("unknown");
00289     else
00290       grid_mapping = strdup(info_field->grid_mapping);
00291     /* Get projection variable ID */
00292     if ( !strcmp(grid_mapping, "Lambert_Conformal") ) {
00293       istat = nc_inq_varid(ncinid, grid_mapping, &projinid); /* get projection variable ID */
00294       if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00295     }
00296     else if ( !strcmp(grid_mapping, "rotated_latitude_longitude") ) {
00297       istat = nc_inq_varid(ncinid, grid_mapping, &projinid); /* get projection variable ID */
00298       if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00299     }
00300     if (proj->name != NULL) {
00301       if ( strcmp(grid_mapping, "Lambert_Conformal") && strcmp(grid_mapping, "Latitude_Longitude") &&
00302            ( !strcmp(proj->name, "Latitude_Longitude") || !strcmp(proj->name, "Lambert_Conformal") ) ) {
00303         (void) free(grid_mapping);
00304         grid_mapping = strdup(proj->name);
00305       }
00306       else if ( strcmp(grid_mapping, "rotated_latitude_longitude") && !strcmp(proj->name, "rotated_pole") ) {
00307         (void) free(grid_mapping);
00308         grid_mapping = strdup("rotated_latitude_longitude");
00309       }
00310     }
00311     if ( !strcmp(grid_mapping, "Lambert_Conformal") ) {
00312       /*              int Lambert_Conformal ;
00313                       Lambert_Conformal:grid_mapping_name = "lambert_conformal_conic" ;
00314                       Lambert_Conformal:standard_parallel = 45.89892f, 47.69601f ;
00315                       Lambert_Conformal:longitude_of_central_meridian = 2.337229f ;
00316                       Lambert_Conformal:latitude_of_projection_origin = 46.8f ;
00317                       Lambert_Conformal:false_easting = 600000.f ;
00318                       Lambert_Conformal:false_northing = 2200000.f ;
00319       */
00320       istat = nc_get_var1_int(ncinid, projinid, 0, &vali);
00321       if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);  
00322       proj->name = strdup(grid_mapping);
00323       istat = get_attribute_str(&(proj->grid_mapping_name), ncinid, projinid, "grid_mapping_name");
00324     
00325       proj_latin = (float *) malloc(2 * sizeof(float));
00326       if (proj_latin == NULL) alloc_error(__FILE__, __LINE__);
00327       istat = nc_get_att_float(ncinid, projinid, "standard_parallel", proj_latin);
00328       proj->latin1 = (double) proj_latin[0];
00329       proj->latin2 = (double) proj_latin[1];
00330       (void) free(proj_latin);
00331     
00332       istat = nc_get_att_double(ncinid, projinid, "longitude_of_central_meridian", &(proj->lonc));
00333       istat = nc_get_att_double(ncinid, projinid, "latitude_of_projection_origin", &(proj->lat0));
00334       istat = nc_get_att_double(ncinid, projinid, "false_easting", &(proj->false_easting));
00335       istat = nc_get_att_double(ncinid, projinid, "false_northing", &(proj->false_northing));
00336     
00337     }
00338     else if ( !strcmp(grid_mapping, "rotated_latitude_longitude") ) {
00339       /*        char rotated_pole ;
00340                 rotated_pole:grid_mapping_name = "rotated_latitude_longitude" ;
00341                 rotated_pole:grid_north_pole_latitude = 39.25 ;
00342                 rotated_pole:grid_north_pole_longitude = -162. ;
00343       */
00344       istat = nc_get_var1_int(ncinid, projinid, 0, &vali);
00345       if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);  
00346       proj->name = strdup("rotated_pole");
00347       istat = get_attribute_str(&(proj->grid_mapping_name), ncinid, projinid, "grid_mapping_name");
00348     
00349       istat = nc_get_att_double(ncinid, projinid, "grid_north_pole_latitude", &(proj->latpole));
00350       istat = nc_get_att_double(ncinid, projinid, "grid_north_pole_longitude", &(proj->lonpole));
00351     
00352     }
00353     else if ( !strcmp(grid_mapping, "Latitude_Longitude") ) {
00354       proj->name = strdup(grid_mapping);
00355       proj->grid_mapping_name = strdup("Latitude_Longitude");
00356       proj->latin1 = 0.0;
00357       proj->latin2 = 0.0;
00358       proj->lonc = 0.0;
00359       proj->lat0 = 0.0;
00360       proj->false_easting = 0.0;
00361       proj->false_northing = 0.0;
00362     }
00363     else if ( !strcmp(grid_mapping, "list") ) {
00364       proj->name = strdup(grid_mapping);
00365       proj->grid_mapping_name = strdup("list");
00366       proj->latin1 = 0.0;
00367       proj->latin2 = 0.0;
00368       proj->lonc = 0.0;
00369       proj->lat0 = 0.0;
00370       proj->false_easting = 0.0;
00371       proj->false_northing = 0.0;
00372     }
00373     else {
00374       (void) fprintf(stderr, "%s: WARNING: No projection parameter available for %s.\n", __FILE__, grid_mapping);
00375       (void) fprintf(stderr, "%s: WARNING: Assuming list of longitude and latitude points.\n", __FILE__);
00376       proj->name = strdup("list");
00377       proj->grid_mapping_name = strdup("list");
00378       proj->latin1 = 0.0;
00379       proj->latin2 = 0.0;
00380       proj->lonc = 0.0;
00381       proj->lat0 = 0.0;
00382       proj->false_easting = 0.0;
00383       proj->false_northing = 0.0;
00384     }
00385     (void) free(grid_mapping);
00386   }
00387 
00388   if (varndims == 3) {
00389     /* Allocate memory and set start and count */
00390     start[0] = t;
00391     start[1] = 0; /* Timestep to retrieve */
00392     start[2] = 0;
00393     count[0] = 1; /* Only get one timestep */
00394     count[1] = (size_t) *nlat;
00395     count[2] = (size_t) *nlon;
00396     /* Allocate memory */
00397     (*buf) = (double *) malloc((*nlat)*(*nlon) * sizeof(double));
00398     if ((*buf) == NULL) alloc_error(__FILE__, __LINE__);
00399   }
00400   else if (varndims == 2) {
00401     /* Allocate memory and set start and count */
00402     start[0] = t;
00403     start[1] = 0; /* Timestep to retrieve */
00404     start[2] = 0;
00405     count[0] = 1; /* Only get one timestep */
00406     count[1] = (size_t) npts; /* List of latitude+longitude points only */
00407     count[2] = 0;
00408     /* Allocate memory */
00409     (*buf) = (double *) malloc(npts * sizeof(double));
00410     if ((*buf) == NULL) alloc_error(__FILE__, __LINE__);
00411   }
00412 
00413   /* Read values from netCDF variable */
00414   istat = nc_get_vara_double(ncinid, varinid, start, count, *buf);
00415   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00416 
00417   /* Close the input netCDF file. */
00418   istat = ncclose(ncinid);
00419   if (istat != NC_NOERR) handle_netcdf_error(istat, __FILE__, __LINE__);
00420 
00421   /* Free memory */
00422   (void) free(tmpstr);
00423 
00424   /* Success status */
00425   return 0;
00426 }


Generated on 12 May 2016 for DSCLIM by  doxygen 1.6.1