REMORA
Energy Research and Forecasting: An Atmospheric Modeling Code
NCFile.H
Go to the documentation of this file.
1 #ifndef _NCFILE_H_
2 #define _NCFILE_H_
3 
4 #include <sstream>
5 #include <string>
6 #include <ctime>
7 #include <atomic>
8 
9 #include "AMReX_FArrayBox.H"
10 #include "AMReX_IArrayBox.H"
11 #include "NCInterface.H"
12 
13 using PlaneVector = amrex::Vector<amrex::FArrayBox>;
14 
15 enum class NC_Data_Dims_Type {
17  Time_SN_WE,
18  SN_WE,
19 };
20 
21 //
22 // NDArray is the datatype designed to hold any data, including scalars, multidimensional
23 // arrays, that read from the NetCDF file.
24 //
25 // The data read from NetCDF file are stored in a continuous memory, and the data layout is descriped
26 // by using a vector (shape). AMRex Box can be constructed using the data shape information, and MultiFab
27 // data array can be setup using the data that stored in the NDArray.
28 //
29 template <typename DataType>
30 struct NDArray
31 {
32  using DType = typename std::remove_const<DataType>::type;
33 
34  // constructor
35  explicit NDArray (const std::string vname, const std::vector<size_t>& vshape)
36  : name{vname}, shape{vshape}, ref_counted{1}, owned{true} {
37  data = new DType [this->ndim()];
38  }
39 
40  // default constructor
41  NDArray () : name{"null"}, ref_counted{1}, owned{false}, data{nullptr} {}
42 
43  // copy constructor
44  NDArray (const NDArray& array) {
45  name = array.name;
46  shape = array.shape;
47  data = array.data;
48  owned = false;
49  ref_counted.fetch_add(1, std::memory_order_relaxed);
50  }
51 
52  // copy assignment
53  NDArray& operator=(const NDArray& array) {
54  name = array.name;
55  shape = array.shape;
56  data = array.data;
57  owned = false;
58  ref_counted.fetch_add(1, std::memory_order_relaxed);
59  return *this;
60  }
61 
62  // destructor
63  ~NDArray () {
64  ref_counted.fetch_sub(1, std::memory_order_acq_rel);
65  if (ref_counted == 1 && owned) delete[] data;
66  }
67 
68  // get the data pointer
69  decltype(auto) get_data () {
70  ref_counted.fetch_add(1, std::memory_order_relaxed);
71  return data;
72  }
73 
74  // get the variable name
75  std::string get_vname () {
76  return name;
77  }
78 
79  // get the variable data shape
80  std::vector<size_t> get_vshape () {
81  return shape;
82  }
83 
84  // return the total number of data
85  size_t ndim () {
86  size_t num = 1;
87  int isize = static_cast<int>(shape.size());
88  for (auto i=0; i<isize; ++i) num *= shape[i];
89  return num;
90  }
91 
92  // set the data shape information
93  void set_vshape (std::vector<size_t> vshape) {
94  shape = vshape;
95  }
96 
97  private:
98  std::string name;
99  std::vector<size_t> shape;
100  mutable std::atomic<size_t> ref_counted;
101  bool owned;
103 };
104 
105 template<typename DType>
106 void ReadNetCDFFile (const std::string& fname, amrex::Vector<std::string> names,
107  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 }
157 
158 /**
159  * Helper function for reading a single variable attribute
160  *
161  * @param fname Name of NetCDF file
162  * @param var_name Name of variable
163  * @param attr_name Name of attribute to read
164  * @returns attribute value
165  */
166 std::string ReadNetCDFVarAttrStr (const std::string& fname,
167  const std::string& var_name,
168  const std::string& attr_name);
169 
170 /**
171  * Helper function for reading data from NetCDF file into a
172  * provided FAB.
173  *
174  * @param iv Index for which variable we are going to fill
175  * @param nc_arrays Arrays of data from NetCDF file
176  * @param var_name Variable name
177  * @param NC_dim_type Dimension type for the variable as stored in the NetCDF file
178  * @param temp FAB where we store the variable data from the NetCDF Arrays
179  */
180 template<class FAB,typename DType>
181 void
183  amrex::Vector<NDArray<float>>& nc_arrays,
184  const std::string& var_name,
185  NC_Data_Dims_Type& NC_dim_type,
186  FAB& temp)
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 }
255 
256 /**
257  * Function to read NetCDF variables and fill the corresponding Array4's
258  *
259  * @param fname Name of the NetCDF file to be read
260  * @param nc_var_names Variable names in the NetCDF file
261  * @param NC_dim_types NetCDF data dimension types
262  * @param fab_vars Fab data we are to fill
263  */
264 template<class FAB,typename DType>
265 void
266 BuildFABsFromNetCDFFile (const amrex::Box& domain,
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)
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 }
320 
321 #endif
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
Definition: NCFile.H:31
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