REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_NCFile.H
Go to the documentation of this file.
1#ifndef _REMORA_NCFILE_H_
2#define _REMORA_NCFILE_H_
3
4#include <sstream>
5#include <string>
6#include <ctime>
7#include <memory>
8
9#include "AMReX_FArrayBox.H"
10#include "AMReX_IArrayBox.H"
11// this includes pnetcdf.h, and mpi.h
12#include "REMORA_NCInterface.H"
13
14using PlaneVector = amrex::Vector<amrex::FArrayBox>;
15
20 SN_WE,
23 BT_Riv,
24};
25
26/**
27 * \brief NDArray is the datatype designed to hold any data, including scalars, multidimensional
28 * arrays, that read from the NetCDF file.
29 *
30 * The data read from NetCDF file are stored in a continuous memory, and the data layout is described
31 * by using a vector (shape). AMRex Box can be constructed using the data shape information, and MultiFab
32 * data array can be setup using the data that stored in the NDArray.
33 */
34template <typename DataType>
35struct NDArray
36{
37 using DType = typename std::remove_const<DataType>::type;
38
39 /** \brief constructor */
40 explicit NDArray (const std::string vname, const std::vector<MPI_Offset>& vshape)
41 : name{vname}, shape{vshape},
42 data{std::shared_ptr<DType[]>(new DType[this->ndim()],
43 std::default_delete<DType[]>())} {}
44
45 /** \brief default constructor */
46 NDArray () : name{"null"}, data{nullptr} {}
47
48 /** \brief get the data pointer */
49 decltype(auto) get_data () {
50 return data.get();
51 }
52
53 /** \brief get the variable name */
54 std::string get_vname () {
55 return name;
56 }
57
58 /** \brief get the variable data shape */
59 std::vector<MPI_Offset> get_vshape () {
60 return shape;
61 }
62
63 /** \brief return the total number of data */
64 size_t ndim () {
65 size_t num = 1;
66 int isize = static_cast<int>(shape.size());
67 for (auto i=0; i<isize; ++i) num *= shape[i];
68 return num;
69 }
70
71 /** \brief set the data shape information */
72 void set_vshape (std::vector<MPI_Offset> vshape) {
73 shape = vshape;
74 }
75
76 private:
77 std::string name;
78 std::vector<MPI_Offset> shape;
79 std::shared_ptr<DType[]> data;
80};
81
82/** \brief Read in data from netcdf file and save to data arrays
83 *
84 * @param[in ] fname file name to read from
85 * @param[in ] names field names to read
86 * @param[inout] arrays vector of data arrays to fill
87 * @param[in ] one_time whether to read in a single time step of data
88 * @param[in ] fill_time if reading a single time, what time step to read from the file?
89 */
90template<typename DType>
91void ReadNetCDFFile (const std::string& fname, amrex::Vector<std::string> names,
92 amrex::Vector<NDArray<DType> >& arrays, bool one_time=false, int fill_time=0)
93{
94 AMREX_ASSERT(arrays.size() == names.size());
95
96 auto ncf = ncutils::NCFile::open(fname, NC_NOCLOBBER );
97 ncmpi_begin_indep_data(ncf.ncid);
98 if (amrex::ParallelDescriptor::IOProcessor())
99 {
100
101 /*
102 // get the dimension information
103 int Time = static_cast<int>(ncf.dim("Time").len());
104 int DateStrLen = static_cast<int>(ncf.dim("DateStrLen").len());
105 int west_east = static_cast<int>(ncf.dim("west_east").len());
106 int south_north = static_cast<int>(ncf.dim("south_north").len());
107 int bottom_top = static_cast<int>(ncf.dim("bottom_top").len());
108 int bottom_top_stag = static_cast<int>(ncf.dim("bottom_top_stag").len());
109 int west_east_stag = static_cast<int>(ncf.dim("west_east_stag").len());
110 int south_north_stag = static_cast<int>(ncf.dim("south_north_stag").len());
111 */
112
113// amrex::Print() << "Reading the dimensions from the netcdf file " << "\n";
114
115 for (auto n=0; n<arrays.size(); ++n) {
116 // read the data from NetCDF file
117 std::string vname_to_write = names[n];
118 std::string vname_to_read = names[n];
119
120 if (!ncf.has_var(vname_to_read)) {
121 amrex::Print() << "Variable name " << vname_to_read << " not found!" << std::endl;
122 }
123
124 // amrex::AllPrint() << "About to read " << vname_to_read
125 // << " while filling the array for " << vname_to_write << std::endl;
126
127 std::vector<MPI_Offset> shape = ncf.var(vname_to_read).shape();
128 arrays[n] = NDArray<DType>(vname_to_read,shape);
129 DType* dataPtr = arrays[n].get_data();
130
131 std::vector<MPI_Offset> start(shape.size(), 0);
132 // If only reading one time step, adjust start and shape accordingly
133 if (one_time) {
134 start[0] = fill_time;
135 shape[0] = 1;
136 }
137
138 // auto numPts = arrays[n].ndim();
139 // amrex::Print() << "NetCDF Variable name = " << vname_to_read << std::endl;
140 // amrex::Print() << "numPts read from NetCDF file/var = " << numPts << std::endl;
141 // amrex::Print() << "Dims of the variable = (";
142 // for (auto &dim:shape)
143 // amrex::Print() << dim << ", " ;
144 // amrex::Print() << ")" << std::endl;
145
146 ncf.var(vname_to_read).get(dataPtr, start, shape);
147 }
148 }
149 ncf.close();
150}
151
152/**
153 * \brief Helper function for reading a single variable attribute
154 */
155std::string ReadNetCDFVarAttrStr (const std::string& fname,
156 const std::string& var_name,
157 const std::string& attr_name);
158
159/**
160 * \brief Helper function for reading data from NetCDF file into a
161 * provided FAB.
162 *
163 * @param iv Index for which variable we are going to fill
164 * @param nc_arrays Arrays of data from NetCDF file
165 * @param var_name Variable name
166 * @param NC_dim_type Dimension type for the variable as stored in the NetCDF file
167 * @param temp FAB where we store the variable data from the NetCDF Arrays
168 */
169template<class FAB,typename DType>
170void
172 amrex::Vector<NDArray<amrex::Real>>& nc_arrays,
173 const std::string& var_name,
174 NC_Data_Dims_Type& NC_dim_type,
175 FAB& temp)
176{
177 int nsx, nsy, nsz;
178 if (NC_dim_type == NC_Data_Dims_Type::SN_WE) {
179 nsz = 1;
180 nsy = nc_arrays[iv].get_vshape()[0];
181 nsx = nc_arrays[iv].get_vshape()[1];
182 // amrex::Print() << "TYPE SN WE " << nsy << " " << nsx << std::endl;
183 } else if (NC_dim_type == NC_Data_Dims_Type::BT_SN_WE) {
184 nsz = nc_arrays[iv].get_vshape()[0];
185 nsy = nc_arrays[iv].get_vshape()[1];
186 nsx = nc_arrays[iv].get_vshape()[2];
187 // amrex::Print() << "TYPE BT SN WE " << nz << " " << nsy << " " << nsx << std::endl;
188 } else if (NC_dim_type == NC_Data_Dims_Type::Time_SN_WE) {
189 nsz = 1;
190 nsy = nc_arrays[iv].get_vshape()[1];
191 nsx = nc_arrays[iv].get_vshape()[2];
192 // amrex::Print() << "TYPE TIME SN WE " << nsy << " " << nsx << std::endl;
193 } else if (NC_dim_type == NC_Data_Dims_Type::Time_BT_SN_WE) {
194 nsz = nc_arrays[iv].get_vshape()[1];
195 nsy = nc_arrays[iv].get_vshape()[2];
196 nsx = nc_arrays[iv].get_vshape()[3];
197 // amrex::Print() << "TYPE TIME BT SN WE " << nsx << " " << nsy << " " << nsz << std::endl;
198 } else if (NC_dim_type == NC_Data_Dims_Type::Time_BT_Riv) {
199 nsz = nc_arrays[iv].get_vshape()[1]; // number of z-coords
200 nsy = 1; // dummy
201 nsx = nc_arrays[iv].get_vshape()[2]; // number of rivers
202 } else if (NC_dim_type == NC_Data_Dims_Type::Time_Riv) {
203 nsz = 1; // number of z-coords
204 nsy = 1; // dummy
205 nsx = nc_arrays[iv].get_vshape()[1]; // number of rivers
206 } else if (NC_dim_type == NC_Data_Dims_Type::BT_Riv) {
207 nsz = nc_arrays[iv].get_vshape()[0]; // number of z-coords
208 nsy = 1; // dummy
209 nsx = nc_arrays[iv].get_vshape()[1]; // number of rivers
210 } else {
211 amrex::Abort("Dont know this NC_Data_Dims_Type");
212 }
213
214 amrex::Box my_box;
215
216 if (var_name == "u" || var_name == "ubar" || var_name == "mask_u" || var_name == "x_u" ||
217 var_name == "y_u" || var_name == "sustr")
218 {
219 my_box.setSmall(amrex::IntVect(0,-1,0));
220 my_box.setBig(amrex::IntVect(nsx-1,nsy-2,nsz-1));
221 my_box.setType(amrex::IndexType(amrex::IntVect(1,0,0)));
222 }
223 else if (var_name == "v" || var_name == "vbar" || var_name == "mask_v" || var_name == "x_v" ||
224 var_name == "y_v" || var_name == "svstr")
225 {
226 my_box.setSmall(amrex::IntVect(-1,0,0));
227 my_box.setBig(amrex::IntVect(nsx-2,nsy-1,nsz-1));
228 my_box.setType(amrex::IndexType(amrex::IntVect(0,1,0)));
229 }
230 else if (var_name == "mask_psi" || var_name == "x_psi" || var_name == "y_psi")
231 {
232 my_box.setSmall(amrex::IntVect(0,0,0));
233 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
234 my_box.setType(amrex::IndexType(amrex::IntVect(1,1,0)));
235 // River data has no ghost cells
236 } else if (NC_dim_type == NC_Data_Dims_Type::Time_BT_Riv ||
237 NC_dim_type == NC_Data_Dims_Type::Time_Riv ||
238 NC_dim_type == NC_Data_Dims_Type::BT_Riv) {
239 my_box.setSmall(amrex::IntVect(0,0,0));
240 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
241 }
242 else // everything cell-centered -- these have one ghost cell
243 {
244 my_box.setSmall(amrex::IntVect(-1,-1,0));
245 my_box.setBig(amrex::IntVect(nsx-2,nsy-2,nsz-1));
246 }
247
248 // amrex::Print() <<" NAME and BOX for this VAR " << var_name << " " << my_box << std::endl;
249
250#ifdef AMREX_USE_GPU
251 // Make sure temp lives on CPU since nc_arrays lives on CPU only
252 temp.resize(my_box,1,amrex::The_Pinned_Arena());
253#else
254 temp.resize(my_box,1);
255#endif
256 amrex::Array4<DType> fab_arr = temp.array();
257
258 int ioff = temp.box().smallEnd()[0];
259 int joff = temp.box().smallEnd()[1];
260
261 auto num_pts_in_box = my_box.numPts();
262 auto num_pts_in_src = nsx*nsy*nsz;
263 // amrex::Print() <<" NUM PTS IN BOX " << num_pts_in_box << std::endl;
264 // amrex::Print() <<" NUM PTS IN SRC " << num_pts_in_src << std::endl;
265 AMREX_ALWAYS_ASSERT(num_pts_in_box >= num_pts_in_src);
266
267 for (int n(0); n < num_pts_in_src; ++n) {
268 int nplane = nsx*nsy;
269 int k = n / nplane;
270 int j = (n - k*nplane) / nsx + joff;
271 int i = n - k*nplane - (j-joff) * nsx + ioff;
272 fab_arr(i,j,k,0) = static_cast<DType>(*(nc_arrays[iv].get_data()+n));
273 }
274}
275
276/**
277 * \brief Function to read NetCDF variables and fill the corresponding Array4's
278 *
279 * @param fname Name of the NetCDF file to be read
280 * @param nc_var_names Variable names in the NetCDF file
281 * @param NC_dim_types NetCDF data dimension types
282 * @param fab_vars Fab data we are to fill
283 * @param one_time Whether to read in a single time step from file
284 * @param fill_time What time step to read from file, if only reading one step
285 */
286template<class FAB,typename DType>
287void
288BuildFABsFromNetCDFFile (const amrex::Box& domain,
289 const std::string &fname,
290 amrex::Vector<std::string> nc_var_names,
291 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types,
292 amrex::Vector<FAB*> fab_vars,
293 bool one_time=false, int fill_time=0)
294{
295 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber(); // I/O rank
296
297 amrex::Vector<NDArray<amrex::Real>> nc_arrays(nc_var_names.size());
298
299
300 ReadNetCDFFile(fname, nc_var_names, nc_arrays, one_time, fill_time); // this is filled only on proc 0
301
302 for (int iv = 0; iv < nc_var_names.size(); iv++)
303 {
304 FAB tmp;
305 if (amrex::ParallelDescriptor::IOProcessor()) {
306 fill_fab_from_arrays<FAB,DType>(iv, nc_arrays, nc_var_names[iv], NC_dim_types[iv], tmp);
307 }
308
309 int ncomp = tmp.nComp();
310 amrex::Box box = tmp.box();
311
312 amrex::ParallelDescriptor::Bcast(&box, 1, ioproc);
313 amrex::ParallelDescriptor::Bcast(&ncomp, 1, ioproc);
314
315 if (!amrex::ParallelDescriptor::IOProcessor()) {
316#ifdef AMREX_USE_GPU
317 tmp.resize(box,ncomp,amrex::The_Pinned_Arena());
318#else
319 tmp.resize(box,ncomp);
320#endif
321 }
322
323 amrex::ParallelDescriptor::Bcast(tmp.dataPtr(), tmp.size(), ioproc);
324
325 // Shift box by the domain lower corner
326 amrex::Box fab_bx = tmp.box();
327 amrex::Dim3 dom_lb = lbound(domain);
328 fab_bx += amrex::IntVect(dom_lb.x,dom_lb.y,dom_lb.z);
329 // fab_vars points to data on device
330 fab_vars[iv]->resize(fab_bx,1);
331#ifdef AMREX_USE_GPU
332 amrex::Gpu::copy(amrex::Gpu::hostToDevice,
333 tmp.dataPtr(), tmp.dataPtr() + tmp.size(),
334 fab_vars[iv]->dataPtr());
335#else
336 // Provided by BaseFab inheritance through FArrayBox
337 fab_vars[iv]->copy(tmp,tmp.box(),0,fab_bx,0,1);
338#endif
339 }
340}
341
342#endif
NC_Data_Dims_Type
void ReadNetCDFFile(const std::string &fname, amrex::Vector< std::string > names, amrex::Vector< NDArray< DType > > &arrays, bool one_time=false, int fill_time=0)
Read in data from netcdf file and save to data arrays.
void BuildFABsFromNetCDFFile(const amrex::Box &domain, const std::string &fname, amrex::Vector< std::string > nc_var_names, amrex::Vector< enum NC_Data_Dims_Type > NC_dim_types, amrex::Vector< FAB * > fab_vars, bool one_time=false, int fill_time=0)
Function to read NetCDF variables and fill the corresponding Array4's.
void fill_fab_from_arrays(int iv, amrex::Vector< NDArray< amrex::Real > > &nc_arrays, const std::string &var_name, NC_Data_Dims_Type &NC_dim_type, FAB &temp)
Helper function for reading data from NetCDF file into a provided FAB.
std::string ReadNetCDFVarAttrStr(const std::string &fname, const std::string &var_name, const std::string &attr_name)
Helper function for reading a single variable attribute.
amrex::Vector< amrex::FArrayBox > PlaneVector
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.
NDArray is the datatype designed to hold any data, including scalars, multidimensional arrays,...
std::string name
std::shared_ptr< DType[]> data
size_t ndim()
return the total number of data
NDArray()
default constructor
std::string get_vname()
get the variable name
std::vector< MPI_Offset > shape
NDArray(const std::string vname, const std::vector< MPI_Offset > &vshape)
constructor
typename std::remove_const< DataType >::type DType
decltype(auto) get_data()
get the data pointer
void set_vshape(std::vector< MPI_Offset > vshape)
set the data shape information
std::vector< MPI_Offset > get_vshape()
get the variable data shape