21#include "pism/util/io/NC_Serial.hh"
32#include "pism/util/pism_utilities.hh"
33#include "pism/util/error_handling.hh"
35#include "pism/util/io/pism_type_conversion.hh"
42 if (return_code != NC_NOERR) {
49 if (return_code != NC_NOERR) {
50 fprintf(stderr,
"%s:%d: %s\n", where.
filename, where.
line_number, nc_strerror(return_code));
64 fprintf(stderr,
"NC_Serial::~NC_Serial: NetCDF file %s is still open\n",
83 stat = nc_open(fname.c_str(), open_mode, &
m_file_id);
88 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
98 stat = nc_create(fname.c_str(), NC_CLOBBER | NC_64BIT_OFFSET, &
m_file_id);
103 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
119 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
133 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
143 int header_size = 200 * 1024;
146 stat = nc__enddef(
m_file_id, header_size, 4, 0, 4);
150 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
164 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
176 stat = nc_def_dim(
m_file_id, name.c_str(), length, &dimid);
180 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
189 stat = nc_inq_dimid(
m_file_id, dimension_name.c_str(), &flag);
191 flag = (stat == NC_NOERR) ? 1 : 0;
194 MPI_Bcast(&flag, 1, MPI_INT, 0,
m_com);
196 exists = (flag == 1);
208 stat = nc_inq_dimid(
m_file_id, dimension_name.c_str(), &dimid);
210 if (stat == NC_NOERR) {
211 stat = nc_inq_dimlen(
m_file_id, dimid, &length);
212 result =
static_cast<unsigned int>(length);
217 MPI_Bcast(&result, 1, MPI_UNSIGNED, 0,
m_com);
218 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
226 std::vector<char> dimname(NC_MAX_NAME + 1, 0);
230 stat = nc_inq_unlimdim(
m_file_id, &dimid);
233 if (stat == NC_NOERR and dimid != -1) {
234 stat = nc_inq_dimname(
m_file_id, dimid, dimname.data());
243 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
244 MPI_Bcast(dimname.data(), NC_MAX_NAME, MPI_CHAR, 0,
m_com);
248 result = dimname.data();
253 const std::vector<std::string> &dims)
const {
257 std::vector<int> dimids;
260 for (
auto d : dims) {
262 stat = nc_inq_dimid(
m_file_id, d.c_str(), &dimid);
263 if (stat != NC_NOERR) {
266 dimids.push_back(dimid);
269 if (stat == NC_NOERR) {
271 static_cast<int>(dims.size()), dimids.data(), &varid);
276 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
282 const std::vector<unsigned int> &start,
283 const std::vector<unsigned int> &
count,
double *op)
const {
289 const std::vector<unsigned int> &start_input,
290 const std::vector<unsigned int> &count_input,
292 std::vector<unsigned int> start = start_input;
293 std::vector<unsigned int>
count = count_input;
295 const int start_tag = 1, count_tag = 2, data_tag = 3, chunk_size_tag = 4;
296 int stat = NC_NOERR, com_size, ndims =
static_cast<int>(start.size());
298 unsigned int local_chunk_size = 1, processor_0_chunk_size = 0;
301 MPI_Comm_size(
m_com, &com_size);
304 for (
int k = 0;
k < ndims; ++
k) {
305 local_chunk_size *=
count[
k];
310 MPI_Reduce(&local_chunk_size, &processor_0_chunk_size, 1, MPI_UNSIGNED, MPI_MAX, 0,
m_com);
314 std::vector<double> processor_0_buffer;
318 processor_0_buffer.resize(processor_0_chunk_size);
323 std::vector<size_t> nc_start(ndims), nc_count(ndims);
326 stat = nc_inq_varid(
m_file_id, variable_name.c_str(), &varid);
329 for (
int r = 0; r < com_size; ++r) {
334 MPI_Recv(start.data(), ndims, MPI_UNSIGNED, r, start_tag,
m_com, &mpi_stat);
335 MPI_Recv(
count.data(), ndims, MPI_UNSIGNED, r, count_tag,
m_com, &mpi_stat);
336 MPI_Recv(&local_chunk_size, 1, MPI_UNSIGNED, r, chunk_size_tag,
m_com, &mpi_stat);
341 for (
int k = 0;
k < ndims; ++
k) {
342 nc_start[
k] = start[
k];
346 stat = nc_get_vara_double(
m_file_id, varid, nc_start.data(), nc_count.data(),
347 processor_0_buffer.data());
351 MPI_Send(processor_0_buffer.data(), (
int)local_chunk_size, MPI_DOUBLE, r, data_tag,
m_com);
353 for (
unsigned int k = 0;
k < local_chunk_size; ++
k) {
354 ip[
k] = processor_0_buffer[
k];
360 MPI_Send(start.data(), ndims, MPI_UNSIGNED, 0, start_tag,
m_com);
361 MPI_Send(
count.data(), ndims, MPI_UNSIGNED, 0, count_tag,
m_com);
362 MPI_Send(&local_chunk_size, 1, MPI_UNSIGNED, 0, chunk_size_tag,
m_com);
364 MPI_Recv(ip, (
int)local_chunk_size, MPI_DOUBLE, 0, data_tag,
m_com, &mpi_stat);
369 const std::vector<unsigned int> &start_input,
370 const std::vector<unsigned int> &count_input,
371 const double *op)
const {
373 std::vector<unsigned int> start = start_input;
374 std::vector<unsigned int>
count = count_input;
375 const int start_tag = 1, count_tag = 2, data_tag = 3, chunk_size_tag = 4;
376 int stat = NC_NOERR, com_size = 0, ndims =
static_cast<int>(start.size());
378 unsigned int local_chunk_size = 1, processor_0_chunk_size = 0;
381 MPI_Comm_size(
m_com, &com_size);
384 for (
int k = 0;
k < ndims; ++
k) {
385 local_chunk_size *=
count[
k];
390 MPI_Reduce(&local_chunk_size, &processor_0_chunk_size, 1, MPI_UNSIGNED, MPI_MAX, 0,
m_com);
394 std::vector<double> processor_0_buffer;
395 processor_0_buffer.resize(processor_0_chunk_size);
400 std::vector<size_t> nc_start(ndims), nc_count(ndims);
403 stat = nc_inq_varid(
m_file_id, variable_name.c_str(), &varid);
406 for (
int r = 0; r < com_size; ++r) {
411 MPI_Recv(start.data(), ndims, MPI_UNSIGNED, r, start_tag,
m_com, &mpi_stat);
412 MPI_Recv(
count.data(), ndims, MPI_UNSIGNED, r, count_tag,
m_com, &mpi_stat);
413 MPI_Recv(&local_chunk_size, 1, MPI_UNSIGNED, r, chunk_size_tag,
m_com, &mpi_stat);
415 MPI_Recv(processor_0_buffer.data(), local_chunk_size, MPI_DOUBLE, r, data_tag,
m_com,
418 for (
unsigned int k = 0;
k < local_chunk_size; ++
k) {
419 processor_0_buffer[
k] = op[
k];
425 for (
int k = 0;
k < ndims; ++
k) {
426 nc_start[
k] = start[
k];
430 stat = nc_put_vara_double(
m_file_id, varid, nc_start.data(), nc_count.data(),
431 processor_0_buffer.data());
435 MPI_Send(start.data(), ndims, MPI_UNSIGNED, 0, start_tag,
m_com);
436 MPI_Send(
count.data(), ndims, MPI_UNSIGNED, 0, count_tag,
m_com);
437 MPI_Send(&local_chunk_size, 1, MPI_UNSIGNED, 0, chunk_size_tag,
m_com);
439 MPI_Send(
const_cast<double *
>(op), (
int)local_chunk_size, MPI_DOUBLE, 0, data_tag,
m_com);
452 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
455 MPI_Bcast(&result, 1, MPI_INT, 0,
m_com);
460 std::vector<std::string> &result)
const {
461 int stat, ndims, varid = -1;
462 std::vector<int> dimids;
465 stat = nc_inq_varid(
m_file_id, variable_name.c_str(), &varid);
467 if (stat == NC_NOERR) {
468 stat = nc_inq_varndims(
m_file_id, varid, &ndims);
472 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
475 MPI_Bcast(&ndims, 1, MPI_INT, 0,
m_com);
482 result.resize(ndims);
483 dimids.resize(ndims);
486 stat = nc_inq_vardimid(
m_file_id, varid, dimids.data());
489 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
494 for (
int k = 0;
k < ndims; ++
k) {
495 std::vector<char> name(NC_MAX_NAME + 1, 0);
498 stat = nc_inq_dimname(
m_file_id, dimids[
k], name.data());
501 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
505 MPI_Bcast(name.data(), name.size(), MPI_CHAR, 0,
m_com);
507 result[
k] = name.data();
521 if (varid >= NC_GLOBAL) {
522 stat = nc_inq_varnatts(
m_file_id, varid, &result);
529 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
532 MPI_Bcast(&result, 1, MPI_INT, 0,
m_com);
540 stat = nc_inq_varid(
m_file_id, variable_name.c_str(), &flag);
542 flag = (stat == NC_NOERR) ? 1 : 0;
545 MPI_Bcast(&flag, 1, MPI_INT, 0,
m_com);
547 exists = (flag == 1);
552 std::vector<char> varname(NC_MAX_NAME + 1, 0);
555 stat = nc_inq_varname(
m_file_id, j, varname.data());
560 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
561 MPI_Bcast(varname.data(), NC_MAX_NAME, MPI_CHAR, 0,
m_com);
565 result = varname.data();
573 std::vector<double> &result)
const {
574 int stat = NC_NOERR, len = 0;
582 if (varid >= NC_GLOBAL) {
583 stat = nc_inq_attlen(
m_file_id, varid, att_name.c_str(), &attlen);
588 if (stat == NC_NOERR) {
589 len =
static_cast<int>(attlen);
594 MPI_Bcast(&len, 1, MPI_INT, 0,
m_com);
605 stat = nc_get_att_double(
m_file_id, varid, att_name.c_str(), result.data());
607 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
612 MPI_Bcast(result.data(), len, MPI_DOUBLE, 0,
m_com);
616static int get_att_text(
int ncid,
int varid,
const std::string &att_name, std::string &result) {
620 stat = nc_inq_attlen(ncid, varid, att_name.c_str(), &attlen);
621 if (stat != NC_NOERR) {
626 std::vector<char> buffer(attlen + 1, 0);
627 stat = nc_get_att_text(ncid, varid, att_name.c_str(), buffer.data());
628 if (stat == NC_NOERR) {
629 result = buffer.data();
616static int get_att_text(
int ncid,
int varid,
const std::string &att_name, std::string &result) {
…}
639static int get_att_string(
int ncid,
int varid,
const std::string &att_name, std::string &result) {
643 stat = nc_inq_attlen(ncid, varid, att_name.c_str(), &attlen);
644 if (stat != NC_NOERR) {
649 std::vector<char *> buffer(attlen + 1, 0);
650 stat = nc_get_att_string(ncid, varid, att_name.c_str(), buffer.data());
651 if (stat == NC_NOERR) {
652 std::vector<std::string> strings(attlen);
653 for (
size_t k = 0;
k < attlen; ++
k) {
654 strings[
k] = buffer[
k];
656 result =
join(strings,
",");
660 stat = nc_free_string(attlen, buffer.data());
639static int get_att_string(
int ncid,
int varid,
const std::string &att_name, std::string &result) {
…}
671 std::string &result)
const {
679 if (varid >= NC_GLOBAL) {
680 nc_type nctype = NC_NAT;
681 stat = nc_inq_atttype(
m_file_id, varid, att_name.c_str(), &nctype);
683 if (stat == NC_NOERR) {
695 }
else if (stat == NC_ENOTATT) {
703 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
706 unsigned int len = result.size();
707 MPI_Bcast(&len, 1, MPI_UNSIGNED, 0,
m_com);
710 MPI_Bcast(&result[0], (
int)len, MPI_CHAR, 0,
m_com);
719 io::Type nctype,
const std::vector<double> &data)
const {
725 if (varid >= NC_GLOBAL) {
727 data.size(), data.data());
734 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
745 const std::string &value)
const {
751 if (varid >= NC_GLOBAL) {
752 stat = nc_put_att_text(
m_file_id, varid, att_name.c_str(), value.size(), value.c_str());
759 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
769 std::string &result)
const {
771 std::vector<char> name(NC_MAX_NAME + 1, 0);
776 if (varid >= NC_GLOBAL) {
777 stat = nc_inq_attname(
m_file_id, varid,
n, name.data());
784 MPI_Bcast(name.data(), NC_MAX_NAME, MPI_CHAR, 0,
m_com);
785 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
789 result = name.data();
803 if (varid >= NC_GLOBAL) {
805 nc_type nctype = NC_NAT;
806 stat = nc_inq_atttype(
m_file_id, varid, att_name.c_str(), &nctype);
807 if (stat == NC_ENOTATT) {
811 tmp =
static_cast<int>(nctype);
818 MPI_Bcast(&tmp, 1, MPI_INT, 0,
m_com);
820 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
832 stat = nc_set_fill(
m_file_id, fillmode, &old_modep);
836 MPI_Bcast(&old_modep, 1, MPI_INT, 0,
m_com);
837 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
846 int stat = nc_inq_format(
m_file_id, &format);
850 MPI_Bcast(&format, 1, MPI_INT, 0,
m_com);
853 case NC_FORMAT_CLASSIC:
854 case NC_FORMAT_64BIT_OFFSET:
856 case NC_FORMAT_64BIT_DATA:
858 case NC_FORMAT_NETCDF4:
859 case NC_FORMAT_NETCDF4_CLASSIC:
871 if (varid >= NC_GLOBAL) {
872 stat = nc_del_att(
m_file_id, varid, att_name.c_str());
877 MPI_Bcast(&stat, 1, MPI_INT, 0,
m_com);
889 if (variable_name ==
"PISM_GLOBAL") {
895 int stat = nc_inq_varid(
m_file_id, variable_name.c_str(), &varid);
897 if (stat == NC_NOERR) {
The PISM wrapper for a subset of the NetCDF C API.
void get_att_text_impl(const std::string &variable_name, const std::string &att_name, std::string &result) const
Gets a text attribute.
virtual void create_impl(const std::string &filename)
Create a NetCDF file.
void def_dim_impl(const std::string &name, size_t length) const
Define a dimension.
void inq_varname_impl(unsigned int j, std::string &result) const
void put_vara_double_impl(const std::string &variable_name, const std::vector< unsigned int > &start, const std::vector< unsigned int > &count, const double *op) const
void set_fill_impl(int fillmode, int &old_modep) const
Sets the fill mode.
std::string get_format() const
void enddef_impl() const
Exit define mode.
void get_vara_double_impl(const std::string &variable_name, const std::vector< unsigned int > &start, const std::vector< unsigned int > &count, double *ip) const
void inq_vardimid_impl(const std::string &variable_name, std::vector< std::string > &result) const
Get dimensions a variable depends on.
void put_att_text_impl(const std::string &variable_name, const std::string &att_name, const std::string &value) const
Writes a text attribute.
void get_att_double_impl(const std::string &variable_name, const std::string &att_name, std::vector< double > &result) const
Gets a double attribute.
void put_att_double_impl(const std::string &variable_name, const std::string &att_name, io::Type xtype, const std::vector< double > &data) const
Writes a double attribute.
void open_impl(const std::string &filename, io::Mode mode)
void inq_varnatts_impl(const std::string &variable_name, int &result) const
Get the number of attributes of a variable.
void inq_varid_impl(const std::string &variable_name, bool &exists) const
Finds a variable and sets the "exists" flag.
void get_var_double(const std::string &variable_name, const std::vector< unsigned int > &start, const std::vector< unsigned int > &count, double *ip) const
Get variable data.
void inq_attname_impl(const std::string &variable_name, unsigned int n, std::string &result) const
Gets the name of a numbered attribute.
virtual void def_var_impl(const std::string &name, io::Type nctype, const std::vector< std::string > &dims) const
Define a variable.
virtual void set_compression_level_impl(int level) const
void close_impl()
Close a NetCDF file.
void del_att_impl(const std::string &variable_name, const std::string &att_name) const
void inq_nvars_impl(int &result) const
Get the number of variables.
void inq_dimlen_impl(const std::string &dimension_name, unsigned int &result) const
Get a dimension length.
void redef_impl() const
Enter define mode.
void inq_dimid_impl(const std::string &dimension_name, bool &exists) const
int get_varid(const std::string &variable_name) const
void inq_unlimdim_impl(std::string &result) const
Get an unlimited dimension.
void inq_atttype_impl(const std::string &variable_name, const std::string &att_name, io::Type &result) const
Gets the type of an attribute.
#define PISM_ERROR_LOCATION
@ PISM_READONLY
open an existing file for reading only
static void check(const ErrorLocation &where, int return_code)
Prints an error message; for debugging.
static void check_and_abort(MPI_Comm com, const ErrorLocation &where, int return_code)
call MPI_Abort() if a NetCDF call failed
static void get_att_string(int ncid, int varid, const std::string &att_name, std::string &result)
static void get_att_text(int ncid, int varid, const std::string &att_name, std::string &result)
std::string join(const std::vector< std::string > &strings, const std::string &separator)
Concatenate strings, inserting separator between elements.
static pism::io::Type nc_type_to_pism_type(int input)
static nc_type pism_type_to_nc_type(pism::io::Type input)