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];
151 std::vector<MPI_Offset> shape = ncf.var(vname_to_read).shape();
153 DType* dataPtr = arrays[n].get_data();
155 std::vector<MPI_Offset> start(shape.size(), 0);
158 start[0] = fill_time;
170 ncf.var(vname_to_read).get(dataPtr, start, shape);
197 const std::string& var_name,
204 nsy = nc_arrays[iv].get_vshape()[0];
205 nsx = nc_arrays[iv].get_vshape()[1];
208 nsz = nc_arrays[iv].get_vshape()[0];
209 nsy = nc_arrays[iv].get_vshape()[1];
210 nsx = nc_arrays[iv].get_vshape()[2];
214 nsy = nc_arrays[iv].get_vshape()[1];
215 nsx = nc_arrays[iv].get_vshape()[2];
218 nsz = nc_arrays[iv].get_vshape()[1];
219 nsy = nc_arrays[iv].get_vshape()[2];
220 nsx = nc_arrays[iv].get_vshape()[3];
223 nsz = nc_arrays[iv].get_vshape()[1];
225 nsx = nc_arrays[iv].get_vshape()[2];
229 nsx = nc_arrays[iv].get_vshape()[1];
231 nsz = nc_arrays[iv].get_vshape()[0];
233 nsx = nc_arrays[iv].get_vshape()[1];
235 amrex::Abort(
"Dont know this NC_Data_Dims_Type");
240 if (var_name ==
"u" || var_name ==
"ubar" || var_name ==
"mask_u" || var_name ==
"x_u" ||
241 var_name ==
"y_u" || var_name ==
"sustr")
243 my_box.setSmall(amrex::IntVect(0,-1,0));
244 my_box.setBig(amrex::IntVect(nsx-1,nsy-2,nsz-1));
245 my_box.setType(amrex::IndexType(amrex::IntVect(1,0,0)));
247 else if (var_name ==
"v" || var_name ==
"vbar" || var_name ==
"mask_v" || var_name ==
"x_v" ||
248 var_name ==
"y_v" || var_name ==
"svstr")
250 my_box.setSmall(amrex::IntVect(-1,0,0));
251 my_box.setBig(amrex::IntVect(nsx-2,nsy-1,nsz-1));
252 my_box.setType(amrex::IndexType(amrex::IntVect(0,1,0)));
254 else if (var_name ==
"mask_psi" || var_name ==
"x_psi" || var_name ==
"y_psi")
256 my_box.setSmall(amrex::IntVect(0,0,0));
257 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
258 my_box.setType(amrex::IndexType(amrex::IntVect(1,1,0)));
263 my_box.setSmall(amrex::IntVect(0,0,0));
264 my_box.setBig(amrex::IntVect(nsx-1,nsy-1,nsz-1));
268 my_box.setSmall(amrex::IntVect(-1,-1,0));
269 my_box.setBig(amrex::IntVect(nsx-2,nsy-2,nsz-1));
276 temp.resize(my_box,1,amrex::The_Pinned_Arena());
278 temp.resize(my_box,1);
280 amrex::Array4<DType> fab_arr = temp.array();
282 int ioff = temp.box().smallEnd()[0];
283 int joff = temp.box().smallEnd()[1];
285 auto num_pts_in_box = my_box.numPts();
286 auto num_pts_in_src = nsx*nsy*nsz;
289 AMREX_ALWAYS_ASSERT(num_pts_in_box >= num_pts_in_src);
291 for (
int n(0); n < num_pts_in_src; ++n) {
292 int nplane = nsx*nsy;
294 int j = (n - k*nplane) / nsx + joff;
295 int i = n - k*nplane - (j-joff) * nsx + ioff;
296 fab_arr(i,j,k,0) =
static_cast<DType
>(*(nc_arrays[iv].get_data()+n));
313 const std::string &fname,
314 amrex::Vector<std::string> nc_var_names,
315 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types,
316 amrex::Vector<FAB*> fab_vars,
317 bool one_time=
false,
int fill_time=0)
319 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
321 amrex::Vector<NDArray<amrex::Real>> nc_arrays(nc_var_names.size());
324 ReadNetCDFFile(fname, nc_var_names, nc_arrays, one_time, fill_time);
326 for (
int iv = 0; iv < nc_var_names.size(); iv++)
329 if (amrex::ParallelDescriptor::IOProcessor()) {
330 fill_fab_from_arrays<FAB,DType>(iv, nc_arrays, nc_var_names[iv], NC_dim_types[iv], tmp);
333 int ncomp = tmp.nComp();
334 amrex::Box box = tmp.box();
336 amrex::ParallelDescriptor::Bcast(&box, 1, ioproc);
337 amrex::ParallelDescriptor::Bcast(&ncomp, 1, ioproc);
339 if (!amrex::ParallelDescriptor::IOProcessor()) {
341 tmp.resize(box,ncomp,amrex::The_Pinned_Arena());
343 tmp.resize(box,ncomp);
347 amrex::ParallelDescriptor::Bcast(tmp.dataPtr(), tmp.size(), ioproc);
350 amrex::Box fab_bx = tmp.box();
351 amrex::Dim3 dom_lb = lbound(domain);
352 fab_bx += amrex::IntVect(dom_lb.x,dom_lb.y,dom_lb.z);
354 fab_vars[iv]->resize(fab_bx,1);
356 amrex::Gpu::copy(amrex::Gpu::hostToDevice,
357 tmp.dataPtr(), tmp.dataPtr() + tmp.size(),
358 fab_vars[iv]->dataPtr());
361 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.