120 amrex::Vector<
NDArray<DType> >& arrays,
bool one_time=
false,
int fill_time=0)
122 AMREX_ASSERT(arrays.size() == names.size());
125 ncmpi_begin_indep_data(ncf.ncid);
126 if (amrex::ParallelDescriptor::IOProcessor())
143 for (
auto n=0; n<arrays.size(); ++n) {
145 std::string vname_to_write = names[n];
146 std::string vname_to_read = names[n];
148 if (!ncf.has_var(vname_to_read)) {
149 amrex::Print() <<
"Variable name " << vname_to_read <<
" not found!" << std::endl;
155 std::vector<MPI_Offset> shape = ncf.var(vname_to_read).shape();
157 DType* dataPtr = arrays[n].get_data();
159 std::vector<MPI_Offset> start(shape.size(), 0);
162 start[0] = fill_time;
174 ncf.var(vname_to_read).get(dataPtr, start, shape);
201 const std::string& var_name,
208 nsy = nc_arrays[iv].get_vshape()[0];
209 nsx = nc_arrays[iv].get_vshape()[1];
212 nsz = nc_arrays[iv].get_vshape()[0];
213 nsy = nc_arrays[iv].get_vshape()[1];
214 nsx = nc_arrays[iv].get_vshape()[2];
218 nsy = nc_arrays[iv].get_vshape()[1];
219 nsx = nc_arrays[iv].get_vshape()[2];
222 nsz = nc_arrays[iv].get_vshape()[1];
223 nsy = nc_arrays[iv].get_vshape()[2];
224 nsx = nc_arrays[iv].get_vshape()[3];
227 nsz = nc_arrays[iv].get_vshape()[1];
229 nsx = nc_arrays[iv].get_vshape()[2];
233 nsx = nc_arrays[iv].get_vshape()[1];
235 nsz = nc_arrays[iv].get_vshape()[0];
237 nsx = nc_arrays[iv].get_vshape()[1];
239 amrex::Abort(
"Dont know this NC_Data_Dims_Type");
244 if (var_name ==
"u" || var_name ==
"ubar" || var_name ==
"mask_u" || var_name ==
"x_u" ||
245 var_name ==
"y_u" || var_name ==
"sustr")
247 my_box.setSmall(amrex::IntVect(0,-1,0));
248 my_box.setBig(amrex::IntVect(nsx-1,nsy-2,nsz-1));
249 my_box.setType(amrex::IndexType(amrex::IntVect(1,0,0)));
251 else if (var_name ==
"v" || var_name ==
"vbar" || var_name ==
"mask_v" || var_name ==
"x_v" ||
252 var_name ==
"y_v" || var_name ==
"svstr")
254 my_box.setSmall(amrex::IntVect(-1,0,0));
255 my_box.setBig(amrex::IntVect(nsx-2,nsy-1,nsz-1));
256 my_box.setType(amrex::IndexType(amrex::IntVect(0,1,0)));
258 else if (var_name ==
"mask_psi" || var_name ==
"x_psi" || var_name ==
"y_psi")
260 my_box.setSmall(amrex::IntVect(0,0,0));
261 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
262 my_box.setType(amrex::IndexType(amrex::IntVect(1,1,0)));
267 my_box.setSmall(amrex::IntVect(0,0,0));
268 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
272 my_box.setSmall(amrex::IntVect(-1,-1,0));
273 my_box.setBig(amrex::IntVect(nsx-2,nsy-2,nsz-1));
280 temp.resize(my_box,1,amrex::The_Pinned_Arena());
282 temp.resize(my_box,1);
284 amrex::Array4<DType> fab_arr = temp.array();
286 int ioff = temp.box().smallEnd()[0];
287 int joff = temp.box().smallEnd()[1];
289 auto num_pts_in_box = my_box.numPts();
290 auto num_pts_in_src = nsx*nsy*nsz;
293 AMREX_ALWAYS_ASSERT(num_pts_in_box >= num_pts_in_src);
295 for (
int n(0); n < num_pts_in_src; ++n) {
296 int nplane = nsx*nsy;
298 int j = (n - k*nplane) / nsx + joff;
299 int i = n - k*nplane - (j-joff) * nsx + ioff;
300 fab_arr(i,j,k,0) =
static_cast<DType
>(*(nc_arrays[iv].get_data()+n));
317 const std::string &fname,
318 amrex::Vector<std::string> nc_var_names,
319 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types,
320 amrex::Vector<FAB*> fab_vars,
321 bool one_time=
false,
int fill_time=0)
323 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
325 amrex::Vector<NDArray<amrex::Real>> nc_arrays(nc_var_names.size());
328 ReadNetCDFFile(fname, nc_var_names, nc_arrays, one_time, fill_time);
330 for (
int iv = 0; iv < nc_var_names.size(); iv++)
333 if (amrex::ParallelDescriptor::IOProcessor()) {
334 fill_fab_from_arrays<FAB,DType>(iv, nc_arrays, nc_var_names[iv], NC_dim_types[iv], tmp);
337 int ncomp = tmp.nComp();
338 amrex::Box box = tmp.box();
340 amrex::ParallelDescriptor::Bcast(&box, 1, ioproc);
341 amrex::ParallelDescriptor::Bcast(&ncomp, 1, ioproc);
343 if (!amrex::ParallelDescriptor::IOProcessor()) {
345 tmp.resize(box,ncomp,amrex::The_Pinned_Arena());
347 tmp.resize(box,ncomp);
351 amrex::ParallelDescriptor::Bcast(tmp.dataPtr(), tmp.size(), ioproc);
354 amrex::Box fab_bx = tmp.box();
355 amrex::Dim3 dom_lb = lbound(domain);
356 fab_bx += amrex::IntVect(dom_lb.x,dom_lb.y,dom_lb.z);
358 fab_vars[iv]->resize(fab_bx,1);
360 amrex::Gpu::copy(amrex::Gpu::hostToDevice,
361 tmp.dataPtr(), tmp.dataPtr() + tmp.size(),
362 fab_vars[iv]->dataPtr());
365 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)
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.