5#include <AMReX_Print.H>
7#define abort_func amrex::Abort
13char recname[NC_MAX_NAME + 1];
15void check_ncmpi_error(
int ierr) {
16 if (ierr != NC_NOERR) {
17 printf(
"\n%s\n\n", ncmpi_strerror(ierr));
18 abort_func(
"Encountered NetCDF error; aborting");
24 check_ncmpi_error(ncmpi_inq_dimname(
ncid,
dimid, recname));
25 return std::string(recname);
30 check_ncmpi_error(ncmpi_inq_dimlen(
ncid,
dimid, &dlen));
35 check_ncmpi_error(ncmpi_inq_varname(
ncid,
varid, recname));
36 return std::string(recname);
41 check_ncmpi_error(ncmpi_inq_varndims(
ncid,
varid, &ndims));
47 std::vector<int> dimids(ndims);
48 std::vector<MPI_Offset> vshape(ndims);
50 check_ncmpi_error(ncmpi_inq_vardimid(
ncid,
varid, dimids.data()));
52 for (
int i = 0; i < ndims; ++i) {
53 check_ncmpi_error(ncmpi_inq_dimlen(
ncid, dimids[i], &vshape[i]));
60 check_ncmpi_error(ncmpi_put_var_double(
ncid,
varid, ptr));
64 check_ncmpi_error(ncmpi_put_var_float(
ncid,
varid, ptr));
68 check_ncmpi_error(ncmpi_put_var_int(
ncid,
varid, ptr));
71void NCVar::put(
const double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
72 check_ncmpi_error(ncmpi_put_vara_double(
ncid,
varid, start.data(), count.data(), dptr));
76void NCVar::put_all(
const double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
77 check_ncmpi_error(ncmpi_put_vara_double_all(
ncid,
varid, start.data(), count.data(), dptr));
81void NCVar::iput(
const double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
83 check_ncmpi_error(ncmpi_iput_vara_double(
ncid,
varid, start.data(), count.data(), dptr, request));
86void NCVar::put(
const double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
87 const std::vector<MPI_Offset> &stride)
const {
88 check_ncmpi_error(ncmpi_put_vars_double(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
91void NCVar::put_all(
const double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
92 const std::vector<MPI_Offset> &stride)
const {
93 check_ncmpi_error(ncmpi_put_vars_double_all(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
96void NCVar::put(
const float *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
97 check_ncmpi_error(ncmpi_put_vara_float(
ncid,
varid, start.data(), count.data(), dptr));
100void NCVar::put_all(
const float *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
101 check_ncmpi_error(ncmpi_put_vara_float_all(
ncid,
varid, start.data(), count.data(), dptr));
104void NCVar::put(
const float *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
105 const std::vector<MPI_Offset> &stride)
const {
106 check_ncmpi_error(ncmpi_put_vars_float(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
109void NCVar::put_all(
const float *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
110 const std::vector<MPI_Offset> &stride)
const {
111 check_ncmpi_error(ncmpi_put_vars_float_all(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
114void NCVar::put(
const int *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
115 check_ncmpi_error(ncmpi_put_vara_int(
ncid,
varid, start.data(), count.data(), dptr));
118void NCVar::put_all(
const int *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
119 check_ncmpi_error(ncmpi_put_vara_int_all(
ncid,
varid, start.data(), count.data(), dptr));
122void NCVar::put(
const int *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
123 const std::vector<MPI_Offset> &stride)
const {
124 check_ncmpi_error(ncmpi_put_vars_int(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
127void NCVar::put_all(
const int *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
128 const std::vector<MPI_Offset> &stride)
const {
129 check_ncmpi_error(ncmpi_put_vars_int_all(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
132void NCVar::put(
const char **dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
133 check_ncmpi_error(ncmpi_put_vara_text(
ncid,
varid, start.data(), count.data(), *dptr));
136void NCVar::put(
const char **dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
137 const std::vector<MPI_Offset> &stride)
const {
138 check_ncmpi_error(ncmpi_put_vars_text(
ncid,
varid, start.data(), count.data(), stride.data(), *dptr));
142 check_ncmpi_error(ncmpi_get_var_double(
ncid,
varid, ptr));
146 check_ncmpi_error(ncmpi_get_var_float(
ncid,
varid, ptr));
150 check_ncmpi_error(ncmpi_get_var_int(
ncid,
varid, ptr));
153void NCVar::get(
double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
154 check_ncmpi_error(ncmpi_get_vara_double(
ncid,
varid, start.data(), count.data(), dptr));
157void NCVar::get_all(
double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
158 check_ncmpi_error(ncmpi_get_vara_double_all(
ncid,
varid, start.data(), count.data(), dptr));
161void NCVar::get(
double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
162 const std::vector<MPI_Offset> &stride)
const {
163 check_ncmpi_error(ncmpi_get_vars_double(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
166void NCVar::get_all(
double *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
167 const std::vector<MPI_Offset> &stride)
const {
168 check_ncmpi_error(ncmpi_get_vars_double_all(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
171void NCVar::get(
float *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
172 check_ncmpi_error(ncmpi_get_vara_float(
ncid,
varid, start.data(), count.data(), dptr));
175void NCVar::get_all(
float *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
176 check_ncmpi_error(ncmpi_get_vara_float_all(
ncid,
varid, start.data(), count.data(), dptr));
179void NCVar::get(
float *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
180 const std::vector<MPI_Offset> &stride)
const {
181 check_ncmpi_error(ncmpi_get_vars_float(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
184void NCVar::get_all(
float *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
185 const std::vector<MPI_Offset> &stride)
const {
186 check_ncmpi_error(ncmpi_get_vars_float_all(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
189void NCVar::get(
int *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
190 check_ncmpi_error(ncmpi_get_vara_int(
ncid,
varid, start.data(), count.data(), dptr));
193void NCVar::get_all(
int *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
194 check_ncmpi_error(ncmpi_get_vara_int(
ncid,
varid, start.data(), count.data(), dptr));
197void NCVar::get(
int *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
198 const std::vector<MPI_Offset> &stride)
const {
199 check_ncmpi_error(ncmpi_get_vars_int(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
202void NCVar::get_all(
int *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
203 const std::vector<MPI_Offset> &stride)
const {
204 check_ncmpi_error(ncmpi_get_vars_int_all(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
207void NCVar::get(
char *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count)
const {
208 check_ncmpi_error(ncmpi_get_vara_text(
ncid,
varid, start.data(), count.data(), dptr));
211void NCVar::get(
char *dptr,
const std::vector<MPI_Offset> &start,
const std::vector<MPI_Offset> &count,
212 const std::vector<MPI_Offset> &stride)
const {
213 check_ncmpi_error(ncmpi_get_vars_text(
ncid,
varid, start.data(), count.data(), stride.data(), dptr));
220 return (ierr == NC_NOERR);
224 check_ncmpi_error(ncmpi_put_att_text(
ncid,
varid,
name.data(), value.size(), value.data()));
227void NCVar::put_attr(
const std::string &name,
const std::vector<double> &value)
const {
228 check_ncmpi_error(ncmpi_put_att_double(
ncid,
varid,
name.data(), NC_DOUBLE, value.size(), value.data()));
232 check_ncmpi_error(ncmpi_put_att_float(
ncid,
varid,
name.data(), NC_FLOAT, value.size(), value.data()));
236 check_ncmpi_error(ncmpi_put_att_int(
ncid,
varid,
name.data(), NC_INT, value.size(), value.data()));
241 std::vector<char> aval;
242 check_ncmpi_error(ncmpi_inq_attlen(
ncid,
varid,
name.data(), &lenp));
244 check_ncmpi_error(ncmpi_get_att_text(
ncid,
varid,
name.data(), aval.data()));
245 return std::string { aval.begin(), aval.end() };
250 check_ncmpi_error(ncmpi_inq_attlen(
ncid,
varid,
name.data(), &lenp));
252 check_ncmpi_error(ncmpi_get_att_double(
ncid,
varid,
name.data(), values.data()));
257 check_ncmpi_error(ncmpi_inq_attlen(
ncid,
varid,
name.data(), &lenp));
259 check_ncmpi_error(ncmpi_get_att_float(
ncid,
varid,
name.data(), values.data()));
264 check_ncmpi_error(ncmpi_inq_attlen(
ncid,
varid,
name.data(), &lenp));
266 check_ncmpi_error(ncmpi_get_att_int(
ncid,
varid,
name.data(), values.data()));
271 check_ncmpi_error(ncmpi_inq_dimid(
ncid, name.data(), &newid));
277 check_ncmpi_error(ncmpi_def_dim(
ncid, name.data(), (MPI_Offset) len, &newid));
283 check_ncmpi_error(ncmpi_def_var(
ncid, name.data(), dtype, 0, NULL, &newid));
289 int ndims = dnames.size();
290 std::vector<int> dimids(ndims);
291 for (
int i = 0; i < ndims; ++i) {
295 check_ncmpi_error(ncmpi_def_var(
ncid, name.data(), dtype, ndims, dimids.data(), &newid));
300 const void *fill_val)
const {
302 int ndims = dnames.size();
303 std::vector<int> dimids(ndims);
304 for (
int i = 0; i < ndims; ++i) {
308 check_ncmpi_error(ncmpi_def_var(
ncid, name.data(), dtype, ndims, dimids.data(), &newid));
309 check_ncmpi_error(ncmpi_def_var_fill(
ncid, newid, NC_FILL, fill_val));
315 check_ncmpi_error(ncmpi_inq_varid(
ncid, name.data(), &varid));
321 check_ncmpi_error(ncmpi_inq(
ncid, &ndims, NULL, NULL, NULL));
327 check_ncmpi_error(ncmpi_inq(
ncid, NULL, NULL, &nattrs, NULL));
333 check_ncmpi_error(ncmpi_inq(
ncid, NULL, &nvars, NULL, NULL));
338 int ierr = ncmpi_inq_dimid(
ncid, name.data(), NULL);
339 return (ierr == NC_NOERR);
343 int ierr = ncmpi_inq_varid(
ncid, name.data(), NULL);
344 return (ierr == NC_NOERR);
350 ierr = ncmpi_inq_att(
ncid, NC_GLOBAL, name.data(), NULL, &lenp);
351 return (ierr == NC_NOERR);
355 check_ncmpi_error(ncmpi_put_att_text(
ncid, NC_GLOBAL, name.data(), value.size(), value.data()));
359 check_ncmpi_error(ncmpi_put_att_double(
ncid, NC_GLOBAL, name.data(), NC_DOUBLE, value.size(), value.data()));
363 check_ncmpi_error(ncmpi_put_att_float(
ncid, NC_GLOBAL, name.data(), NC_FLOAT, value.size(), value.data()));
367 check_ncmpi_error(ncmpi_put_att_int(
ncid, NC_GLOBAL, name.data(), NC_INT, value.size(), value.data()));
372 std::vector<char> aval;
373 check_ncmpi_error(ncmpi_inq_attlen(
ncid, NC_GLOBAL, name.data(), &lenp));
375 check_ncmpi_error(ncmpi_get_att_text(
ncid, NC_GLOBAL, name.data(), aval.data()));
376 return std::string { aval.begin(), aval.end() };
381 check_ncmpi_error(ncmpi_inq_attlen(
ncid, NC_GLOBAL, name.data(), &lenp));
383 check_ncmpi_error(ncmpi_get_att_double(
ncid, NC_GLOBAL, name.data(), values.data()));
388 check_ncmpi_error(ncmpi_inq_attlen(
ncid, NC_GLOBAL, name.data(), &lenp));
390 check_ncmpi_error(ncmpi_get_att_float(
ncid, NC_GLOBAL, name.data(), values.data()));
395 check_ncmpi_error(ncmpi_inq_attlen(
ncid, NC_GLOBAL, name.data(), &lenp));
397 check_ncmpi_error(ncmpi_get_att_int(
ncid, NC_GLOBAL, name.data(), values.data()));
401 std::vector<NCDim> adims;
403 adims.reserve(ndims);
404 for (
int i = 0; i < ndims; ++i) {
411 std::vector<NCVar> avars;
413 avars.reserve(nvars);
414 for (
int i = 0; i < nvars; ++i) {
422 ierr = ncmpi_redef(
ncid);
425 if (ierr == NC_EINDEFINE)
428 check_ncmpi_error(ierr);
432 check_ncmpi_error(ncmpi_enddef(
ncid));
437 check_ncmpi_error(ncmpi_create(comm, name.data(), cmode, info, &
ncid));
443 check_ncmpi_error(ncmpi_open(comm, name.data(), cmode, info, &
ncid));
448 std::vector<int> statuses(num_requests);
449 ncmpi_wait_all(
ncid, num_requests, requests, &statuses[0]);
454 check_ncmpi_error(ncmpi_close(
ncid));
459 check_ncmpi_error(ncmpi_close(
ncid));
void put_attr(const std::string &name, const std::string &value) const
Set file attribute to value.
std::vector< NCDim > all_dims() const
Return a list of all dimensions defined in this group.
NCVar def_scalar(const std::string &name, const nc_type dtype) const
Define a scalar variable, i.e., 0-dimensional array.
NCVar def_array(const std::string &name, const nc_type dtype, const std::vector< std::string > &) const
Define an array.
std::vector< NCVar > all_vars() const
Return a list of all variables defined in this group.
NCDim dim(const std::string &) const
Get the dimension instance by name.
void enter_def_mode() const
Enter definition mode (not needed for NetCDF4 format)
bool has_var(const std::string &) const
Check if a variable exists by name.
bool has_attr(const std::string &) const
Check if an attribute exists.
NCDim def_dim(const std::string &, const size_t len) const
Define new dimension.
static NCFile create(const std::string &name, const int cmode=NC_CLOBBER|NC_MPIIO, MPI_Comm comm=MPI_COMM_WORLD, MPI_Info info=MPI_INFO_NULL)
Create a file.
void exit_def_mode() const
Exit definition mode.
std::string get_attr(const std::string &name) const
Read file attribute from file.
int num_dimensions() const
Number of dimensions.
void wait_all(int num_requests, int *requests)
wait for non-blocking calls to finish
bool has_dim(const std::string &) const
Check if a dimension exists by name.
int num_attributes() const
Number of attributes.
int num_variables() const
Number of variables.
static NCFile open(const std::string &name, const int cmode=NC_NOWRITE, MPI_Comm comm=MPI_COMM_WORLD, MPI_Info info=MPI_INFO_NULL)
Open an existing file.
void close()
Close file object.
NCVar var(const std::string &) const
Get the variable instance by name.
NCVar def_array_fill(const std::string &name, const nc_type dtype, const std::vector< std::string > &dnames, const void *fill_val) const
Define an array with a fill value.
Representation of NetCDF dimension.
const int dimid
Dimension ID used with NetCDF API.
const int ncid
File/Group Identifier.
MPI_Offset len() const
Length of this dimension.
std::string name() const
Name of this dimension.
Representation of a NetCDF variable.
bool has_attr(const std::string &name) const
Whether a variable has an attribute with name.
const int ncid
File/Group identifier.
const int varid
Variable ID used with NetCDF API.
void iput(const double *dptr, const std::vector< MPI_Offset > &start, const std::vector< MPI_Offset > &count, int *request) const
Write out a slice of data with with strides (see hyperslab definition in NetCDF)
std::string name() const
Name of this variable.
void put(const double *ptr) const
Write out the entire double variable.
void get(double *ptr) const
Read the entire variable from file.
std::string get_attr(const std::string &name) const
Read attribute from file.
void put_all(const double *dptr, const std::vector< MPI_Offset > &start, const std::vector< MPI_Offset > &count) const
Write out a slice of data, collective.
void put_attr(const std::string &name, const std::string &value) const
Set attribute "name" to "value".
std::vector< MPI_Offset > shape() const
Shape of the array (size in each array dimension)
void get_all(double *dptr, const std::vector< MPI_Offset > &start, const std::vector< MPI_Offset > &count) const
Read a chunk of data from the file, collective.
int ndim() const
Number of array dimensions for this variable.