REMORA
Energy Research and Forecasting: An Atmospheric Modeling Code
NCFile.H File Reference
#include <sstream>
#include <string>
#include <ctime>
#include <atomic>
#include "AMReX_FArrayBox.H"
#include "AMReX_IArrayBox.H"
#include "NCInterface.H"
Include dependency graph for NCFile.H:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

struct  NDArray< DataType >
 

Typedefs

using PlaneVector = amrex::Vector< amrex::FArrayBox >
 

Enumerations

enum class  NC_Data_Dims_Type { Time_BT_SN_WE , Time_SN_WE , SN_WE }
 

Functions

template<typename DType >
void ReadNetCDFFile (const std::string &fname, amrex::Vector< std::string > names, amrex::Vector< NDArray< DType > > &arrays)
 
std::string ReadNetCDFVarAttrStr (const std::string &fname, const std::string &var_name, const std::string &attr_name)
 
template<class FAB , typename DType >
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)
 
template<class FAB , typename DType >
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)
 

Typedef Documentation

◆ PlaneVector

using PlaneVector = amrex::Vector<amrex::FArrayBox>

Enumeration Type Documentation

◆ NC_Data_Dims_Type

enum NC_Data_Dims_Type
strong
Enumerator
Time_BT_SN_WE 
Time_SN_WE 
SN_WE 

Function Documentation

◆ BuildFABsFromNetCDFFile()

template<class FAB , typename DType >
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 
)

Function to read NetCDF variables and fill the corresponding Array4's

Parameters
fnameName of the NetCDF file to be read
nc_var_namesVariable names in the NetCDF file
NC_dim_typesNetCDF data dimension types
fab_varsFab data we are to fill
271 {
272  int ioproc = amrex::ParallelDescriptor::IOProcessorNumber(); // I/O rank
273 
274  amrex::Vector<NDArray<float>> nc_arrays(nc_var_names.size());
275 
276  if (amrex::ParallelDescriptor::IOProcessor())
277  {
278  ReadNetCDFFile(fname, nc_var_names, nc_arrays);
279  }
280 
281  for (int iv = 0; iv < nc_var_names.size(); iv++)
282  {
283  FAB tmp;
284  if (amrex::ParallelDescriptor::IOProcessor()) {
285  fill_fab_from_arrays<FAB,DType>(iv, nc_arrays, nc_var_names[iv], NC_dim_types[iv], tmp);
286  }
287 
288  int ncomp = tmp.nComp();
289  amrex::Box box = tmp.box();
290 
291  amrex::ParallelDescriptor::Bcast(&box, 1, ioproc);
292  amrex::ParallelDescriptor::Bcast(&ncomp, 1, ioproc);
293 
294  if (!amrex::ParallelDescriptor::IOProcessor()) {
295 #ifdef AMREX_USE_GPU
296  tmp.resize(box,ncomp,amrex::The_Pinned_Arena());
297 #else
298  tmp.resize(box,ncomp);
299 #endif
300  }
301 
302  amrex::ParallelDescriptor::Bcast(tmp.dataPtr(), tmp.size(), ioproc);
303 
304  // Shift box by the domain lower corner
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);
308  // fab_vars points to data on device
309  fab_vars[iv]->resize(fab_bx,1);
310 #ifdef AMREX_USE_GPU
311  amrex::Gpu::copy(amrex::Gpu::hostToDevice,
312  tmp.dataPtr(), tmp.dataPtr() + tmp.size(),
313  fab_vars[iv]->dataPtr());
314 #else
315  // Provided by BaseFab inheritance through FArrayBox
316  fab_vars[iv]->copy(tmp,tmp.box(),0,fab_bx,0,1);
317 #endif
318  }
319 }
void ReadNetCDFFile(const std::string &fname, amrex::Vector< std::string > names, amrex::Vector< NDArray< DType > > &arrays)
Definition: NCFile.H:106
Here is the call graph for this function:

◆ fill_fab_from_arrays()

template<class FAB , typename DType >
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 
)

Helper function for reading data from NetCDF file into a provided FAB.

Parameters
ivIndex for which variable we are going to fill
nc_arraysArrays of data from NetCDF file
var_nameVariable name
NC_dim_typeDimension type for the variable as stored in the NetCDF file
tempFAB where we store the variable data from the NetCDF Arrays
187 {
188  int nsx, nsy, nsz;
189  if (NC_dim_type == NC_Data_Dims_Type::SN_WE) {
190  nsz = 1;
191  nsy = nc_arrays[iv].get_vshape()[0];
192  nsx = nc_arrays[iv].get_vshape()[1];
193  // amrex::Print() << "TYPE SN WE " << nsy << " " << nsx << std::endl;
194  } else if (NC_dim_type == NC_Data_Dims_Type::Time_SN_WE) {
195  nsz = 1;
196  nsy = nc_arrays[iv].get_vshape()[1];
197  nsx = nc_arrays[iv].get_vshape()[2];
198  // amrex::Print() << "TYPE TIME SN WE " << nsy << " " << nsx << std::endl;
199  } else if (NC_dim_type == NC_Data_Dims_Type::Time_BT_SN_WE) {
200  nsz = nc_arrays[iv].get_vshape()[1];
201  nsy = nc_arrays[iv].get_vshape()[2];
202  nsx = nc_arrays[iv].get_vshape()[3];
203  // amrex::Print() << "TYPE TIME BT SN WE " << nsx << " " << nsy << " " << nsz << std::endl;
204  } else {
205  amrex::Abort("Dont know this NC_Data_Dims_Type");
206  }
207 
208  amrex::Box my_box;
209 
210  if (var_name == "u" || var_name == "ubar")
211  {
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)));
215  }
216  else if (var_name == "v" || var_name == "vbar")
217  {
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)));
221  }
222  else // everything cell-centered -- these have one ghost cell
223  {
224  my_box.setSmall(amrex::IntVect(-1,-1,0));
225  my_box.setBig(amrex::IntVect(nsx-2,nsy-2,nsz-1));
226  }
227 
228  // amrex::Print() <<" NAME and BOX for this VAR " << var_name << " " << my_box << std::endl;
229 
230 #ifdef AMREX_USE_GPU
231  // Make sure temp lives on CPU since nc_arrays lives on CPU only
232  temp.resize(my_box,1,amrex::The_Pinned_Arena());
233 #else
234  temp.resize(my_box,1);
235 #endif
236  amrex::Array4<DType> fab_arr = temp.array();
237 
238  int ioff = temp.box().smallEnd()[0];
239  int joff = temp.box().smallEnd()[1];
240 
241  auto num_pts_in_box = my_box.numPts();
242  auto num_pts_in_src = nsx*nsy*nsz;
243  // amrex::Print() <<" NUM PTS IN BOX " << num_pts_in_box << std::endl;
244  // amrex::Print() <<" NUM PTS IN SRC " << num_pts_in_src << std::endl;
245  AMREX_ALWAYS_ASSERT(num_pts_in_box >= num_pts_in_src);
246 
247  for (int n(0); n < num_pts_in_src; ++n) {
248  int nplane = nsx*nsy;
249  int k = n / nplane;
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));
253  }
254 }

◆ ReadNetCDFFile()

template<typename DType >
void ReadNetCDFFile ( const std::string &  fname,
amrex::Vector< std::string >  names,
amrex::Vector< NDArray< DType > > &  arrays 
)
108 {
109  AMREX_ASSERT(arrays.size() == names.size());
110 
111  if (amrex::ParallelDescriptor::IOProcessor())
112  {
113  auto ncf = ncutils::NCFile::open(fname, NC_CLOBBER | NC_NETCDF4);
114 
115  /*
116  // get the dimension information
117  int Time = static_cast<int>(ncf.dim("Time").len());
118  int DateStrLen = static_cast<int>(ncf.dim("DateStrLen").len());
119  int west_east = static_cast<int>(ncf.dim("west_east").len());
120  int south_north = static_cast<int>(ncf.dim("south_north").len());
121  int bottom_top = static_cast<int>(ncf.dim("bottom_top").len());
122  int bottom_top_stag = static_cast<int>(ncf.dim("bottom_top_stag").len());
123  int west_east_stag = static_cast<int>(ncf.dim("west_east_stag").len());
124  int south_north_stag = static_cast<int>(ncf.dim("south_north_stag").len());
125  int bdy_width = static_cast<int>(ncf.dim("bdy_width").len());
126  */
127 
128  amrex::Print() << "Reading the dimensions from the netcdf file " << "\n";
129 
130  for (auto n=0; n<arrays.size(); ++n) {
131  // read the data from NetCDF file
132  std::string vname_to_write = names[n];
133  std::string vname_to_read = names[n];
134 
135  // amrex::AllPrint() << "About to read " << vname_to_read
136  // << " while filling the array for " << vname_to_write << std::endl;
137 
138  std::vector<size_t> shape = ncf.var(vname_to_read).shape();
139  arrays[n] = NDArray<DType>(vname_to_read,shape);
140  DType* dataPtr = arrays[n].get_data();
141 
142  std::vector<size_t> start(shape.size(), 0);
143 
144  // auto numPts = arrays[n].ndim();
145  // amrex::Print() << "NetCDF Variable name = " << vname_to_read << std::endl;
146  // amrex::Print() << "numPts read from NetCDF file/var = " << numPts << std::endl;
147  // amrex::Print() << "Dims of the variable = (";
148  // for (auto &dim:shape)
149  // amrex::Print() << dim << ", " ;
150  // amrex::Print() << ")" << std::endl;
151 
152  ncf.var(vname_to_read).get(dataPtr, start, shape);
153  }
154  ncf.close();
155  }
156 }
static NCFile open(const std::string &name, const int cmode=NC_NOWRITE)
Definition: NCInterface.cpp:577
Definition: NCFile.H:31

Referenced by BuildFABsFromNetCDFFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ReadNetCDFVarAttrStr()

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

Parameters
fnameName of NetCDF file
var_nameName of variable
attr_nameName of attribute to read
Returns
attribute value
10 {
11  std::string attr_val;
12  if (amrex::ParallelDescriptor::IOProcessor())
13  {
14  auto ncf = ncutils::NCFile::open(fname, NC_CLOBBER | NC_NETCDF4);
15  attr_val = ncf.var(var_name).get_attr(attr_name);
16  ncf.close();
17  }
18  return attr_val;
19 }
Here is the call graph for this function: