REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_NCInterface.cpp
Go to the documentation of this file.
1#include <cstdio>
2
4#include <AMReX.H>
5#include <AMReX_Print.H>
6
7#define abort_func amrex::Abort
8
9namespace ncutils {
10
11namespace {
12
13char recname[NC_MAX_NAME + 1];
14
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");
19 }
20}
21} // namespace
22
23std::string NCDim::name() const {
24 check_ncmpi_error(ncmpi_inq_dimname(ncid, dimid, recname));
25 return std::string(recname);
26}
27
28MPI_Offset NCDim::len() const {
29 MPI_Offset dlen;
30 check_ncmpi_error(ncmpi_inq_dimlen(ncid, dimid, &dlen));
31 return dlen;
32}
33
34std::string NCVar::name() const {
35 check_ncmpi_error(ncmpi_inq_varname(ncid, varid, recname));
36 return std::string(recname);
37}
38
39int NCVar::ndim() const {
40 int ndims;
41 check_ncmpi_error(ncmpi_inq_varndims(ncid, varid, &ndims));
42 return ndims;
43}
44
45std::vector<MPI_Offset> NCVar::shape() const {
46 int ndims = ndim();
47 std::vector<int> dimids(ndims);
48 std::vector<MPI_Offset> vshape(ndims);
49
50 check_ncmpi_error(ncmpi_inq_vardimid(ncid, varid, dimids.data()));
51
52 for (int i = 0; i < ndims; ++i) {
53 check_ncmpi_error(ncmpi_inq_dimlen(ncid, dimids[i], &vshape[i]));
54 }
55
56 return vshape;
57}
58
59void NCVar::put(const double *ptr) const {
60 check_ncmpi_error(ncmpi_put_var_double(ncid, varid, ptr));
61}
62
63void NCVar::put(const float *ptr) const {
64 check_ncmpi_error(ncmpi_put_var_float(ncid, varid, ptr));
65}
66
67void NCVar::put(const int *ptr) const {
68 check_ncmpi_error(ncmpi_put_var_int(ncid, varid, ptr));
69}
70
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));
73}
74
75//! Write out a slice of data, collective
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));
78}
79
80//! Write out a slice of data, non-blocking
81void NCVar::iput(const double *dptr, const std::vector<MPI_Offset> &start, const std::vector<MPI_Offset> &count,
82 int *request) const {
83 check_ncmpi_error(ncmpi_iput_vara_double(ncid, varid, start.data(), count.data(), dptr, request));
84}
85
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));
89}
90
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));
94}
95
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));
98}
99
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));
102}
103
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));
107}
108
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));
112}
113
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));
116}
117
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));
120}
121
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));
125}
126
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));
130}
131
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));
134}
135
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));
139}
140
141void NCVar::get(double *ptr) const {
142 check_ncmpi_error(ncmpi_get_var_double(ncid, varid, ptr));
143}
144
145void NCVar::get(float *ptr) const {
146 check_ncmpi_error(ncmpi_get_var_float(ncid, varid, ptr));
147}
148
149void NCVar::get(int *ptr) const {
150 check_ncmpi_error(ncmpi_get_var_int(ncid, varid, ptr));
151}
152
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));
155}
156
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));
159}
160
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));
164}
165
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));
169}
170
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));
173}
174
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));
177}
178
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));
182}
183
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));
187}
188
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));
191}
192
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));
195}
196
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));
200}
201
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));
205}
206
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));
209}
210
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));
214}
215
216bool NCVar::has_attr(const std::string &name) const {
217 int ierr;
218 MPI_Offset lenp;
219 ierr = ncmpi_inq_att(ncid, varid, name.data(), NULL, &lenp);
220 return (ierr == NC_NOERR);
221}
222
223void NCVar::put_attr(const std::string &name, const std::string &value) const {
224 check_ncmpi_error(ncmpi_put_att_text(ncid, varid, name.data(), value.size(), value.data()));
225}
226
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()));
229}
230
231void NCVar::put_attr(const std::string &name, const std::vector<float> &value) const {
232 check_ncmpi_error(ncmpi_put_att_float(ncid, varid, name.data(), NC_FLOAT, value.size(), value.data()));
233}
234
235void NCVar::put_attr(const std::string &name, const std::vector<int> &value) const {
236 check_ncmpi_error(ncmpi_put_att_int(ncid, varid, name.data(), NC_INT, value.size(), value.data()));
237}
238
239std::string NCVar::get_attr(const std::string &name) const {
240 MPI_Offset lenp;
241 std::vector<char> aval;
242 check_ncmpi_error(ncmpi_inq_attlen(ncid, varid, name.data(), &lenp));
243 aval.resize(lenp);
244 check_ncmpi_error(ncmpi_get_att_text(ncid, varid, name.data(), aval.data()));
245 return std::string { aval.begin(), aval.end() };
246}
247
248void NCVar::get_attr(const std::string &name, std::vector<double> &values) const {
249 MPI_Offset lenp;
250 check_ncmpi_error(ncmpi_inq_attlen(ncid, varid, name.data(), &lenp));
251 values.resize(lenp);
252 check_ncmpi_error(ncmpi_get_att_double(ncid, varid, name.data(), values.data()));
253}
254
255void NCVar::get_attr(const std::string &name, std::vector<float> &values) const {
256 MPI_Offset lenp;
257 check_ncmpi_error(ncmpi_inq_attlen(ncid, varid, name.data(), &lenp));
258 values.resize(lenp);
259 check_ncmpi_error(ncmpi_get_att_float(ncid, varid, name.data(), values.data()));
260}
261
262void NCVar::get_attr(const std::string &name, std::vector<int> &values) const {
263 MPI_Offset lenp;
264 check_ncmpi_error(ncmpi_inq_attlen(ncid, varid, name.data(), &lenp));
265 values.resize(lenp);
266 check_ncmpi_error(ncmpi_get_att_int(ncid, varid, name.data(), values.data()));
267}
268
269NCDim NCFile::dim(const std::string &name) const {
270 int newid;
271 check_ncmpi_error(ncmpi_inq_dimid(ncid, name.data(), &newid));
272 return NCDim { ncid, newid };
273}
274
275NCDim NCFile::def_dim(const std::string &name, const size_t len) const {
276 int newid;
277 check_ncmpi_error(ncmpi_def_dim(ncid, name.data(), (MPI_Offset) len, &newid));
278 return NCDim { ncid, newid };
279}
280
281NCVar NCFile::def_scalar(const std::string &name, const nc_type dtype) const {
282 int newid;
283 check_ncmpi_error(ncmpi_def_var(ncid, name.data(), dtype, 0, NULL, &newid));
284 return NCVar { ncid, newid };
285}
286
287NCVar NCFile::def_array(const std::string &name, const nc_type dtype, const std::vector<std::string> &dnames) const {
288 int newid;
289 int ndims = dnames.size();
290 std::vector<int> dimids(ndims);
291 for (int i = 0; i < ndims; ++i) {
292 dimids[i] = dim(dnames[i]).dimid;
293 }
294
295 check_ncmpi_error(ncmpi_def_var(ncid, name.data(), dtype, ndims, dimids.data(), &newid));
296 return NCVar { ncid, newid };
297}
298
299NCVar NCFile::def_array_fill(const std::string &name, const nc_type dtype, const std::vector<std::string> &dnames,
300 const void *fill_val) const {
301 int newid;
302 int ndims = dnames.size();
303 std::vector<int> dimids(ndims);
304 for (int i = 0; i < ndims; ++i) {
305 dimids[i] = dim(dnames[i]).dimid;
306 }
307
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));
310 return NCVar { ncid, newid };
311}
312
313NCVar NCFile::var(const std::string &name) const {
314 int varid;
315 check_ncmpi_error(ncmpi_inq_varid(ncid, name.data(), &varid));
316 return NCVar { ncid, varid };
317}
318
320 int ndims;
321 check_ncmpi_error(ncmpi_inq(ncid, &ndims, NULL, NULL, NULL));
322 return ndims;
323}
324
326 int nattrs;
327 check_ncmpi_error(ncmpi_inq(ncid, NULL, NULL, &nattrs, NULL));
328 return nattrs;
329}
330
332 int nvars;
333 check_ncmpi_error(ncmpi_inq(ncid, NULL, &nvars, NULL, NULL));
334 return nvars;
335}
336
337bool NCFile::has_dim(const std::string &name) const {
338 int ierr = ncmpi_inq_dimid(ncid, name.data(), NULL);
339 return (ierr == NC_NOERR);
340}
341
342bool NCFile::has_var(const std::string &name) const {
343 int ierr = ncmpi_inq_varid(ncid, name.data(), NULL);
344 return (ierr == NC_NOERR);
345}
346
347bool NCFile::has_attr(const std::string &name) const {
348 int ierr;
349 MPI_Offset lenp;
350 ierr = ncmpi_inq_att(ncid, NC_GLOBAL, name.data(), NULL, &lenp);
351 return (ierr == NC_NOERR);
352}
353
354void NCFile::put_attr(const std::string &name, const std::string &value) const {
355 check_ncmpi_error(ncmpi_put_att_text(ncid, NC_GLOBAL, name.data(), value.size(), value.data()));
356}
357
358void NCFile::put_attr(const std::string &name, const std::vector<double> &value) const {
359 check_ncmpi_error(ncmpi_put_att_double(ncid, NC_GLOBAL, name.data(), NC_DOUBLE, value.size(), value.data()));
360}
361
362void NCFile::put_attr(const std::string &name, const std::vector<float> &value) const {
363 check_ncmpi_error(ncmpi_put_att_float(ncid, NC_GLOBAL, name.data(), NC_FLOAT, value.size(), value.data()));
364}
365
366void NCFile::put_attr(const std::string &name, const std::vector<int> &value) const {
367 check_ncmpi_error(ncmpi_put_att_int(ncid, NC_GLOBAL, name.data(), NC_INT, value.size(), value.data()));
368}
369
370std::string NCFile::get_attr(const std::string &name) const {
371 MPI_Offset lenp;
372 std::vector<char> aval;
373 check_ncmpi_error(ncmpi_inq_attlen(ncid, NC_GLOBAL, name.data(), &lenp));
374 aval.resize(lenp);
375 check_ncmpi_error(ncmpi_get_att_text(ncid, NC_GLOBAL, name.data(), aval.data()));
376 return std::string { aval.begin(), aval.end() };
377}
378
379void NCFile::get_attr(const std::string &name, std::vector<double> &values) const {
380 MPI_Offset lenp;
381 check_ncmpi_error(ncmpi_inq_attlen(ncid, NC_GLOBAL, name.data(), &lenp));
382 values.resize(lenp);
383 check_ncmpi_error(ncmpi_get_att_double(ncid, NC_GLOBAL, name.data(), values.data()));
384}
385
386void NCFile::get_attr(const std::string &name, std::vector<float> &values) const {
387 MPI_Offset lenp;
388 check_ncmpi_error(ncmpi_inq_attlen(ncid, NC_GLOBAL, name.data(), &lenp));
389 values.resize(lenp);
390 check_ncmpi_error(ncmpi_get_att_float(ncid, NC_GLOBAL, name.data(), values.data()));
391}
392
393void NCFile::get_attr(const std::string &name, std::vector<int> &values) const {
394 MPI_Offset lenp;
395 check_ncmpi_error(ncmpi_inq_attlen(ncid, NC_GLOBAL, name.data(), &lenp));
396 values.resize(lenp);
397 check_ncmpi_error(ncmpi_get_att_int(ncid, NC_GLOBAL, name.data(), values.data()));
398}
399
400std::vector<NCDim> NCFile::all_dims() const {
401 std::vector<NCDim> adims;
402 int ndims = num_dimensions();
403 adims.reserve(ndims);
404 for (int i = 0; i < ndims; ++i) {
405 adims.emplace_back(NCDim { ncid, i });
406 }
407 return adims;
408}
409
410std::vector<NCVar> NCFile::all_vars() const {
411 std::vector<NCVar> avars;
412 int nvars = num_variables();
413 avars.reserve(nvars);
414 for (int i = 0; i < nvars; ++i) {
415 avars.emplace_back(NCVar { ncid, i });
416 }
417 return avars;
418}
419
421 int ierr;
422 ierr = ncmpi_redef(ncid);
423
424 // Ignore already in define mode error
425 if (ierr == NC_EINDEFINE)
426 return;
427 // Handle all other errors
428 check_ncmpi_error(ierr);
429}
430
432 check_ncmpi_error(ncmpi_enddef(ncid));
433}
434
435NCFile NCFile::create(const std::string &name, const int cmode, MPI_Comm comm, MPI_Info info) {
436 int ncid;
437 check_ncmpi_error(ncmpi_create(comm, name.data(), cmode, info, &ncid));
438 return NCFile(ncid);
439}
440
441NCFile NCFile::open(const std::string &name, const int cmode, MPI_Comm comm, MPI_Info info) {
442 int ncid;
443 check_ncmpi_error(ncmpi_open(comm, name.data(), cmode, info, &ncid));
444 return NCFile(ncid);
445}
446
447void NCFile::wait_all(int num_requests, int *requests) {
448 std::vector<int> statuses(num_requests);
449 ncmpi_wait_all(ncid, num_requests, requests, &statuses[0]);
450}
451
453 if (is_open)
454 check_ncmpi_error(ncmpi_close(ncid));
455}
456
458 is_open = false;
459 check_ncmpi_error(ncmpi_close(ncid));
460}
461} // namespace ncutils
#define abort_func
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.