9 #include "AMReX_FArrayBox.H"
10 #include "AMReX_IArrayBox.H"
29 template <
typename DataType>
32 using DType =
typename std::remove_const<DataType>::type;
35 explicit NDArray (
const std::string vname,
const std::vector<size_t>& vshape)
49 ref_counted.fetch_add(1, std::memory_order_relaxed);
58 ref_counted.fetch_add(1, std::memory_order_relaxed);
64 ref_counted.fetch_sub(1, std::memory_order_acq_rel);
70 ref_counted.fetch_add(1, std::memory_order_relaxed);
87 int isize =
static_cast<int>(
shape.size());
88 for (
auto i=0; i<isize; ++i) num *=
shape[i];
105 template<
typename DType>
109 AMREX_ASSERT(arrays.size() == names.size());
111 if (amrex::ParallelDescriptor::IOProcessor())
128 amrex::Print() <<
"Reading the dimensions from the netcdf file " <<
"\n";
130 for (
auto n=0; n<arrays.size(); ++n) {
132 std::string vname_to_write = names[n];
133 std::string vname_to_read = names[n];
138 std::vector<size_t> shape = ncf.var(vname_to_read).shape();
140 DType* dataPtr = arrays[n].get_data();
142 std::vector<size_t> start(shape.size(), 0);
152 ncf.var(vname_to_read).get(dataPtr, start, shape);
167 const std::string& var_name,
168 const std::string& attr_name);
180 template<
class FAB,
typename DType>
184 const std::string& var_name,
191 nsy = nc_arrays[iv].get_vshape()[0];
192 nsx = nc_arrays[iv].get_vshape()[1];
196 nsy = nc_arrays[iv].get_vshape()[1];
197 nsx = nc_arrays[iv].get_vshape()[2];
200 nsz = nc_arrays[iv].get_vshape()[1];
201 nsy = nc_arrays[iv].get_vshape()[2];
202 nsx = nc_arrays[iv].get_vshape()[3];
205 amrex::Abort(
"Dont know this NC_Data_Dims_Type");
210 if (var_name ==
"u" || var_name ==
"ubar")
212 my_box.setSmall(amrex::IntVect(0,-1,0));
213 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
214 my_box.setType(amrex::IndexType(amrex::IntVect(1,0,0)));
216 else if (var_name ==
"v" || var_name ==
"vbar")
218 my_box.setSmall(amrex::IntVect(-1,0,0));
219 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
220 my_box.setType(amrex::IndexType(amrex::IntVect(0,1,0)));
224 my_box.setSmall(amrex::IntVect(-1,-1,0));
225 my_box.setBig(amrex::IntVect(nsx-2,nsy-2,nsz-1));
232 temp.resize(my_box,1,amrex::The_Pinned_Arena());
234 temp.resize(my_box,1);
236 amrex::Array4<DType> fab_arr = temp.array();
238 int ioff = temp.box().smallEnd()[0];
239 int joff = temp.box().smallEnd()[1];
241 auto num_pts_in_box = my_box.numPts();
242 auto num_pts_in_src = nsx*nsy*nsz;
245 AMREX_ALWAYS_ASSERT(num_pts_in_box >= num_pts_in_src);
247 for (
int n(0); n < num_pts_in_src; ++n) {
248 int nplane = nsx*nsy;
250 int j = (n - k*nplane) / nsx + joff;
251 int i = n - k*nplane - (j-joff) * nsx + ioff;
252 fab_arr(i,j,k,0) =
static_cast<DType
>(*(nc_arrays[iv].get_data()+n));
264 template<
class FAB,
typename DType>
267 const std::string &fname,
268 amrex::Vector<std::string> nc_var_names,
269 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types,
270 amrex::Vector<FAB*> fab_vars)
272 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
274 amrex::Vector<NDArray<float>> nc_arrays(nc_var_names.size());
276 if (amrex::ParallelDescriptor::IOProcessor())
281 for (
int iv = 0; iv < nc_var_names.size(); iv++)
284 if (amrex::ParallelDescriptor::IOProcessor()) {
285 fill_fab_from_arrays<FAB,DType>(iv, nc_arrays, nc_var_names[iv], NC_dim_types[iv], tmp);
288 int ncomp = tmp.nComp();
289 amrex::Box box = tmp.box();
291 amrex::ParallelDescriptor::Bcast(&box, 1, ioproc);
292 amrex::ParallelDescriptor::Bcast(&ncomp, 1, ioproc);
294 if (!amrex::ParallelDescriptor::IOProcessor()) {
296 tmp.resize(box,ncomp,amrex::The_Pinned_Arena());
298 tmp.resize(box,ncomp);
302 amrex::ParallelDescriptor::Bcast(tmp.dataPtr(), tmp.size(), ioproc);
305 amrex::Box fab_bx = tmp.box();
306 amrex::Dim3 dom_lb = lbound(domain);
307 fab_bx += amrex::IntVect(dom_lb.x,dom_lb.y,dom_lb.z);
309 fab_vars[iv]->resize(fab_bx,1);
311 amrex::Gpu::copy(amrex::Gpu::hostToDevice,
312 tmp.dataPtr(), tmp.dataPtr() + tmp.size(),
313 fab_vars[iv]->dataPtr());
316 fab_vars[iv]->copy(tmp,tmp.box(),0,fab_bx,0,1);
NC_Data_Dims_Type
Definition: NCFile.H:15
void fill_fab_from_arrays(int iv, amrex::Vector< NDArray< float >> &nc_arrays, const std::string &var_name, NC_Data_Dims_Type &NC_dim_type, FAB &temp)
Definition: NCFile.H:182
void ReadNetCDFFile(const std::string &fname, amrex::Vector< std::string > names, amrex::Vector< NDArray< DType > > &arrays)
Definition: NCFile.H:106
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)
Definition: NCFile.H:266
std::string ReadNetCDFVarAttrStr(const std::string &fname, const std::string &var_name, const std::string &attr_name)
Definition: NCFile.cpp:7
amrex::Vector< amrex::FArrayBox > PlaneVector
Definition: REMORA_PhysBCFunct.H:19
static NCFile open(const std::string &name, const int cmode=NC_NOWRITE)
Definition: NCInterface.cpp:577
bool owned
Definition: NCFile.H:101
std::string name
Definition: NCFile.H:98
size_t ndim()
Definition: NCFile.H:85
NDArray()
Definition: NCFile.H:41
std::string get_vname()
Definition: NCFile.H:75
NDArray(const NDArray &array)
Definition: NCFile.H:44
NDArray & operator=(const NDArray &array)
Definition: NCFile.H:53
void set_vshape(std::vector< size_t > vshape)
Definition: NCFile.H:93
NDArray(const std::string vname, const std::vector< size_t > &vshape)
Definition: NCFile.H:35
std::atomic< size_t > ref_counted
Definition: NCFile.H:100
typename std::remove_const< DataType >::type DType
Definition: NCFile.H:32
std::vector< size_t > shape
Definition: NCFile.H:99
DType * data
Definition: NCFile.H:102
decltype(auto) get_data()
Definition: NCFile.H:69
std::vector< size_t > get_vshape()
Definition: NCFile.H:80
~NDArray()
Definition: NCFile.H:63