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)
94 AMREX_ASSERT(arrays.size() == names.size());
97 ncmpi_begin_indep_data(ncf.ncid);
98 if (amrex::ParallelDescriptor::IOProcessor())
115 for (
auto n=0; n<arrays.size(); ++n) {
117 std::string vname_to_write = names[n];
118 std::string vname_to_read = names[n];
120 if (!ncf.has_var(vname_to_read)) {
121 amrex::Print() <<
"Variable name " << vname_to_read <<
" not found!" << std::endl;
127 std::vector<MPI_Offset> shape = ncf.var(vname_to_read).shape();
129 DType* dataPtr = arrays[n].get_data();
131 std::vector<MPI_Offset> start(shape.size(), 0);
134 start[0] = fill_time;
146 ncf.var(vname_to_read).get(dataPtr, start, shape);
180 const std::string& var_name,
182 FAB& temp, amrex::IntVect custom_ngrow)
187 nsy = nc_arrays[iv].get_vshape()[0];
188 nsx = nc_arrays[iv].get_vshape()[1];
191 nsz = nc_arrays[iv].get_vshape()[0];
192 nsy = nc_arrays[iv].get_vshape()[1];
193 nsx = nc_arrays[iv].get_vshape()[2];
197 nsy = nc_arrays[iv].get_vshape()[1];
198 nsx = nc_arrays[iv].get_vshape()[2];
201 nsz = nc_arrays[iv].get_vshape()[1];
202 nsy = nc_arrays[iv].get_vshape()[2];
203 nsx = nc_arrays[iv].get_vshape()[3];
206 nsz = nc_arrays[iv].get_vshape()[1];
208 nsx = nc_arrays[iv].get_vshape()[2];
212 nsx = nc_arrays[iv].get_vshape()[1];
214 nsz = nc_arrays[iv].get_vshape()[0];
216 nsx = nc_arrays[iv].get_vshape()[1];
218 amrex::Abort(
"Dont know this NC_Data_Dims_Type");
223 if (var_name ==
"u" || var_name ==
"ubar" || var_name ==
"mask_u" || var_name ==
"x_u" ||
224 var_name ==
"y_u" || var_name ==
"sustr")
226 my_box.setSmall(amrex::IntVect(0,-1,0));
227 my_box.setBig(amrex::IntVect(nsx-1,nsy-2,nsz-1));
228 my_box.setType(amrex::IndexType(amrex::IntVect(1,0,0)));
230 else if (var_name ==
"v" || var_name ==
"vbar" || var_name ==
"mask_v" || var_name ==
"x_v" ||
231 var_name ==
"y_v" || var_name ==
"svstr")
233 my_box.setSmall(amrex::IntVect(-1,0,0));
234 my_box.setBig(amrex::IntVect(nsx-2,nsy-1,nsz-1));
235 my_box.setType(amrex::IndexType(amrex::IntVect(0,1,0)));
237 else if (var_name ==
"mask_psi" || var_name ==
"x_psi" || var_name ==
"y_psi")
239 my_box.setSmall(amrex::IntVect(0,0,0));
240 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
241 my_box.setType(amrex::IndexType(amrex::IntVect(1,1,0)));
246 my_box.setSmall(amrex::IntVect(0,0,0));
247 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
251 my_box.setSmall(amrex::IntVect(-1,-1,0));
252 my_box.setBig(amrex::IntVect(nsx-2,nsy-2,nsz-1));
255 my_box.grow(custom_ngrow - amrex::IntVect(1,1,0));
261 temp.resize(my_box,1,amrex::The_Pinned_Arena());
263 temp.resize(my_box,1);
265 amrex::Array4<DType> fab_arr = temp.array();
267 int ioff = temp.box().smallEnd()[0];
268 int joff = temp.box().smallEnd()[1];
270 auto num_pts_in_box = my_box.numPts();
271 auto num_pts_in_src = nsx*nsy*nsz;
274 AMREX_ALWAYS_ASSERT(num_pts_in_box >= num_pts_in_src);
276 for (
int n(0); n < num_pts_in_src; ++n) {
277 int nplane = nsx*nsy;
279 int j = (n - k*nplane) / nsx + joff;
280 int i = n - k*nplane - (j-joff) * nsx + ioff;
281 fab_arr(i,j,k,0) =
static_cast<DType
>(*(nc_arrays[iv].get_data()+n));
299 const std::string &fname,
300 amrex::Vector<std::string> nc_var_names,
301 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types,
302 amrex::Vector<FAB*> fab_vars,
303 bool one_time=
false,
int fill_time=0,
304 amrex::IntVect custom_ngrow=amrex::IntVect(1,1,0))
306 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
308 amrex::Vector<NDArray<amrex::Real>> nc_arrays(nc_var_names.size());
311 ReadNetCDFFile(fname, nc_var_names, nc_arrays, one_time, fill_time);
313 for (
int iv = 0; iv < nc_var_names.size(); iv++)
316 if (amrex::ParallelDescriptor::IOProcessor()) {
317 fill_fab_from_arrays<FAB,DType>(iv, nc_arrays, nc_var_names[iv], NC_dim_types[iv], tmp, custom_ngrow);
320 int ncomp = tmp.nComp();
321 amrex::Box box = tmp.box();
323 amrex::ParallelDescriptor::Bcast(&box, 1, ioproc);
324 amrex::ParallelDescriptor::Bcast(&ncomp, 1, ioproc);
326 if (!amrex::ParallelDescriptor::IOProcessor()) {
328 tmp.resize(box,ncomp,amrex::The_Pinned_Arena());
330 tmp.resize(box,ncomp);
334 amrex::ParallelDescriptor::Bcast(tmp.dataPtr(), tmp.size(), ioproc);
337 amrex::Box fab_bx = tmp.box();
338 amrex::Dim3 dom_lb = lbound(domain);
339 fab_bx += amrex::IntVect(dom_lb.x,dom_lb.y,dom_lb.z);
341 fab_vars[iv]->resize(fab_bx,1);
343 amrex::Gpu::copy(amrex::Gpu::hostToDevice,
344 tmp.dataPtr(), tmp.dataPtr() + tmp.size(),
345 fab_vars[iv]->dataPtr());
348 fab_vars[iv]->copy(tmp,tmp.box(),0,fab_bx,0,1);
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, amrex::IntVect custom_ngrow=amrex::IntVect(1, 1, 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, amrex::IntVect custom_ngrow)
Helper function for reading data from NetCDF file into a provided FAB.