ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_ProbCommon.H
Go to the documentation of this file.
1 #ifndef ERF_PROBCOMMON_H_
2 #define ERF_PROBCOMMON_H_
3 
4 #include <AMReX_ParmParse.H>
5 #include <AMReX_Geometry.H>
6 #include <AMReX_FArrayBox.H>
7 #include <AMReX_MultiFab.H>
8 
9 #include "ERF_DataStruct.H"
10 #include "ERF_EOS.H"
11 #include "ERF_HSEUtils.H"
12 #include "ERF_TileNoZ.H"
13 
17 };
18 
19 template <typename T>
20 T erf_get_single_value (std::istream& is, int n)
21 {
22  T r = 0;
23  std::string line;
24  if (std::getline(is, line)) {
25  std::stringstream ss(line);
26  ss >> r;
27  if (ss.fail()) {
28  amrex::Abort("Failed to read");
29  }
30  if (ss >> r) {
31  amrex::Print() << "Trying to read line " << n << " in this file " << std::endl;
32  amrex::Abort("Wrong format: more than one number in this line");
33  }
34  } else {
35  amrex::Abort("Unable to read this line");
36  }
37  return r;
38 }
39 
40 void
41 init_my_custom_terrain (const amrex::Geometry& geom,
42  amrex::FArrayBox& terrain_fab,
43  const amrex::Real& time);
44 
45 /**
46  * Class to hold problem-specific routines
47 */
49 {
50 public:
51 
52  /**
53  * Virtual destructor to avoid data leakage with derived class
54  */
55  virtual ~ProblemBase () = default;
56 
57  /**
58  * Function to initialize the hydrostatic reference density
59  *
60  * @param[out] rho_hse hydrostatic reference density
61  * @param[in] z_phys_nd height coordinate at nodes
62  * @param[in] z_phys_cc height coordinate at cell centers
63  * @param[in] geom container for geometric information
64  */
65  virtual void
66  erf_init_dens_hse_dry (amrex::MultiFab& /*rho_hse*/,
67  std::unique_ptr<amrex::MultiFab>& /*z_phys_nd*/,
68  std::unique_ptr<amrex::MultiFab>& /*z_phys_cc*/,
69  amrex::Geometry const& /*geom*/,
70  const amrex::Vector<amrex::Real>& /*stretched_dz_h*/,
71  bool /*is_constant_dz*/, bool /*is_stretched_dz*/)
72  {
73  amrex::Print() << "Hydrostatically balanced density was NOT set"
74  << " -- an appropriate init_type should probably have been specified"
75  << " (e.g., input_sounding, WRFInput, or Metgrid)"
76  << std::endl;
77  amrex::Error("Should never call this version of erf_init_dens_hse_dry for "+name()+" problem");
78  }
79 
80  virtual void
81  erf_init_const_dens_hse(amrex::MultiFab& /*rho_hse*/)
82  {
83  amrex::Error("Should never call this version of erf_init_const_dens_hse for "+name()+" problem");
84  }
85 
86  virtual void
87  erf_init_const_dens_and_th_hse (amrex::MultiFab& /*rho_hse*/, amrex::MultiFab& /* p_hse*/,
88  amrex::MultiFab& /* pi_hse*/, amrex::MultiFab& /*th_hse*/,
89  amrex::MultiFab& /* qv_hse*/, amrex::Real /*l_rdOcp*/)
90  {
91  amrex::Error("Should never call this version of erf_init_const_dens_and_th_hse for "+name()+" problem");
92  }
93 
94  virtual void
95  erf_init_const_dens_and_linear_th_hse (amrex::MultiFab& /*rho_hse*/, amrex::MultiFab& /* p_hse*/,
96  amrex::MultiFab& /* pi_hse*/, amrex::MultiFab& /*th_hse*/,
97  amrex::MultiFab& /* qv_hse*/, amrex::Real /*l_rdOcp*/,
98  std::unique_ptr<amrex::MultiFab>& /*z_phys_cc*/)
99  {
100  amrex::Error("Should never call this version of erf_init_const_dens_and_linear_th_hse for "+name()+" problem");
101  }
102 
103  virtual void
104  erf_init_dens_hse_moist (amrex::MultiFab& /*rho_hse*/,
105  std::unique_ptr<amrex::MultiFab>& /*z_phys_nd*/,
106  amrex::Geometry const& /*geom*/)
107  {
108 
109  }
110 
111  /**
112  * Function to perform custom initialization of a test problem
113  *
114  * @param[in] bx cell-centered box on which to initialize scalars
115  * @param[out] state cell-centered variables to be filled in this routine
116  * @param[out] r_hse hydrostatic reference density
117  * @param[out] p_hse hydrostatic reference pressure
118  * @param[in] z_nd height coordinate at nodes
119  * @param[in] z_cc height coordinate at cell centers
120  * @param[in] mf_m map factor on cell centers
121  * @param[in] mf_u map factor on x-faces
122  * @param[in] mf_v map factor on y-faces
123  * @param[in] sc SolverChoice structure that carries parameters
124  * @param[in] lev AMR level
125  */
126  virtual void
127  init_custom_pert (const amrex::Box& /*bx*/,
128  amrex::Array4<amrex::Real const> const& /*state*/,
129  amrex::Array4<amrex::Real > const& /*state_pert*/,
130  amrex::Array4<amrex::Real > const& /*r_hse*/,
131  amrex::Array4<amrex::Real > const& /*p_hse*/,
132  amrex::Array4<amrex::Real const> const& /*z_nd*/,
133  amrex::Array4<amrex::Real const> const& /*z_cc*/,
134  amrex::GeometryData const& /*geomdata*/,
135  amrex::Array4<amrex::Real const> const& /*mf_m*/,
136  const SolverChoice& /*sc*/, const int /*lev*/)
137  {
138  amrex::Print() << "No perturbation to background fields supplied for "
139  << name() << " problem" << std::endl;
140  }
141 
142  /**
143  * Function to perform custom initialization of the velocities
144  *
145  * @param[in] xbx box on which to initialize x_vel_pert
146  * @param[in] ybx box on which to initialize y_vel_pert
147  * @param[in] zbx box on which to initialize z_vel_pert
148  * @param[out] x_vel_pert x-component of velocity perturbation to be filled in this routine
149  * @param[out] y_vel_pert y-component of velocity perturbation to be filled in this routine
150  * @param[out] z_vel_pert z-component of velocity perturbation to be filled in this routine
151  * @param[in] z_nd height coordinate at nodes
152  * @param[in] mf_m map factor on cell centers
153  * @param[in] mf_u map factor on x-faces
154  * @param[in] mf_v map factor on y-faces
155  * @param[in] sc SolverChoice structure that carries parameters
156  * @param[in] lev AMR level
157  */
158  virtual void
159  init_custom_pert_vels (const amrex::Box& /*xbx*/,
160  const amrex::Box& /*ybx*/,
161  const amrex::Box& /*zbx*/,
162  amrex::Array4<amrex::Real > const& /*x_vel_pert*/,
163  amrex::Array4<amrex::Real > const& /*y_vel_pert*/,
164  amrex::Array4<amrex::Real > const& /*z_vel_pert*/,
165  amrex::Array4<amrex::Real const> const& /*z_nd*/,
166  amrex::GeometryData const& /*geomdata*/,
167  amrex::Array4<amrex::Real const> const& /*mf_u*/,
168  amrex::Array4<amrex::Real const> const& /*mf_v*/,
169  const SolverChoice& /*sc*/, const int /*lev*/)
170  {
171  amrex::Print() << "No perturbation velocities supplied for " << name() << " problem" << std::endl;
172  }
173 
174  /**
175  * Function to update user-specified temperature source terms that can
176  * vary with time and height.
177  *
178  * @param[in] time current time
179  * @param[out] rhotheta_source forcing profile
180  * @param[in] geom container for geometric information
181  * @param[in] z_phys_cc height coordinate at cell centers
182  */
183  virtual void
185  amrex::MultiFab* src,
186  const amrex::Geometry& /*geom*/,
187  std::unique_ptr<amrex::MultiFab>& /*z_phys_cc*/)
188  {
189  if (src->empty()) return;
190 
191  amrex::Warning("Temperature forcing not defined for "+name()+" problem");
192  for ( amrex::MFIter mfi(*src, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
193  {
194  const auto &box = mfi.tilebox();
195  const amrex::Array4<amrex::Real>& src_arr = src->array(mfi);
196  // src is a spatial function if erf.spatial_rhotheta_forcing = true
197  // otherwise, qsrc_arr is defined only in over Z (box x and y dimensions are 1)
198  ParallelFor(box, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
199  // const amrex::Real z_cc = prob_lo[2] + (k+myhalf)* dx[2];
200  // set RHS term of RhoTheta equation based on time, z_cc here
201  src_arr(i, j, k) = zero;
202  });
203  }
204  }
205 
206  /**
207  * Function to update user-specified moisture source terms that can
208  * vary with time and height.
209  *
210  * @param[in] time current time
211  * @param[out] rhoqt_source moisture forcing profile
212  * @param[in] geom container for geometric information
213  * @param[in] z_phys_cc height coordinate at cell centers
214  */
215  virtual void
217  amrex::MultiFab* qsrc,
218  const amrex::Geometry& /*geom*/,
219  std::unique_ptr<amrex::MultiFab>& /*z_phys_cc*/)
220  {
221  if (qsrc->empty()) return;
222 
223  amrex::Warning("Moisture forcing not defined for "+name()+" problem");
224  for ( amrex::MFIter mfi(*qsrc, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
225  {
226  const auto &box = mfi.tilebox();
227  const amrex::Array4<amrex::Real>& qsrc_arr = qsrc->array(mfi);
228  // src is a spatial function if erf.spatial_moisture_forcing = true
229  // otherwise, qsrc_arr is defined only in over Z (box x and y dimensions are 1)
230  ParallelFor(box, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
231  // const amrex::Real z_cc = prob_lo[2] + (k+myhalf)* dx[2];
232  // set RHS term of RhoQ1 equation based on time, z_cc here
233  qsrc_arr(i, j, k) = zero;
234  });
235  }
236  }
237 
238  /**
239  * Function to update the vertical velocity profile, used to add subsidence
240  * source terms for x-mom, y-mom, rho*theta, rho*Q1, and rho*Q2.
241  *
242  * TODO: Currently, this is only called by InitData, so there is no time
243  * dependence.
244  *
245  * @param[in] time current time
246  * @param[out] wbar w vel forcing profile
247  * @param[in] geom container for geometric information
248  * @param[in] z_phys_cc height coordinate at cell centers
249  */
250  virtual void
252  amrex::Vector<amrex::Real>& wbar,
253  amrex::Gpu::DeviceVector<amrex::Real>& d_wbar,
254  const amrex::MultiFab& /*state*/,
255  const amrex::Geometry& geom,
256  std::unique_ptr<amrex::MultiFab>& /*z_phys_cc*/)
257  {
258  if (wbar.empty()) return;
259 
260  amrex::Warning("Moisture forcing not defined for "+name()+" problem");
261 
262  const int khi = geom.Domain().bigEnd()[2];
263  // const amrex::Real* prob_lo = geom.ProbLo();
264  // const auto dx = geom.CellSize();
265  for (int k = 0; k <= khi; k++)
266  {
267  // const amrex::Real z_cc = prob_lo[2] + (k+myhalf)* dx[2];
268  // set vertical velocity profile based on time, z_cc here
269  wbar[k] = zero;
270  }
271 
272  // Copy from host version to device version
273  amrex::Gpu::copy(amrex::Gpu::hostToDevice, wbar.begin(), wbar.end(), d_wbar.begin());
274  }
275 
276  /**
277  * Function to update user-specified geostrophic wind profile.
278  *
279  * @param[in] time current time
280  * @param[out] u_geos geostrophic wind profile
281  * @param[out] v_geos geostrophic wind profile
282  * @param[in] geom container for geometric information
283  * @param[in] z_phys_cc height coordinate at cell centers
284  */
285  virtual void
287  amrex::Vector<amrex::Real>& u_geos,
288  amrex::Gpu::DeviceVector<amrex::Real>& d_u_geos,
289  amrex::Vector<amrex::Real>& v_geos,
290  amrex::Gpu::DeviceVector<amrex::Real>& d_v_geos,
291  const amrex::Geometry& geom,
292  std::unique_ptr<amrex::MultiFab>& /*z_phys_cc*/)
293  {
294  if (u_geos.empty()) return;
295 
296  amrex::Warning("Geostrophic wind profile not defined for "+name()+" problem");
297 
298  const int khi = geom.Domain().bigEnd()[2];
299  // const amrex::Real* prob_lo = geom.ProbLo();
300  // const auto dx = geom.CellSize();
301  for (int k = 0; k <= khi; k++)
302  {
303  // const amrex::Real z_cc = prob_lo[2] + (k+myhalf)* dx[2];
304  // set RHS term of RhoTheta equation based on time, z_cc here
305  u_geos[k] = zero;
306  v_geos[k] = zero;
307  }
308 
309  // Copy from host version to device version
310  amrex::Gpu::copy(amrex::Gpu::hostToDevice, u_geos.begin(), u_geos.end(), d_u_geos.begin());
311  amrex::Gpu::copy(amrex::Gpu::hostToDevice, v_geos.begin(), v_geos.end(), d_v_geos.begin());
312  }
313 
314  /**
315  * Function to perform custom initialization of buildings
316  *
317  * @param[in] geom container for geometric information
318  * @param[out] z_phys_nd height coordinate at nodes
319  * @param[in] time current time
320  */
321  void
322  init_buildings_surface (const amrex::Geometry& geom,
323  amrex::FArrayBox& buildings_fab,
324  const amrex::Real& time)
325  {
326  // Check if a valid text file exists for the buildings
327  std::string fname;
328  amrex::ParmParse pp("erf");
329  auto valid_fname = pp.query("buildings_file_name",fname);
330 
331  if (valid_fname) {
332  read_custom_terrain(fname,false,geom,buildings_fab,time);
333  } else {
334  init_my_custom_terrain(geom, buildings_fab, time);
335  }
336  }
337 
338 
339  /**
340  * Function to perform custom initialization of terrain
341  *
342  * @param[in] geom container for geometric information
343  * @param[out] z_phys_nd height coordinate at nodes
344  * @param[in] time current time
345  */
346  void
347  init_terrain_surface (const amrex::Geometry& geom,
348  amrex::FArrayBox& terrain_fab,
349  const amrex::Real& time)
350  {
351  // Check if a valid text file exists for the terrain
352  std::string fname, fname_usgs;
353  amrex::ParmParse pp("erf");
354  auto valid_fname = pp.query("terrain_file_name",fname);
355  auto valid_fname_USGS = pp.query("terrain_file_name_USGS",fname_usgs);
356 
357  bool is_usgs;
358 
359  if (valid_fname) {
360  is_usgs = false;
361  read_custom_terrain(fname,is_usgs,geom,terrain_fab,time);
362 
363  } else if (valid_fname_USGS) {
364  is_usgs = true;
365  read_custom_terrain(fname_usgs,is_usgs,geom,terrain_fab,time);
366 
367  } else {
368  init_my_custom_terrain (geom, terrain_fab, time);
369  }
370  }
371 
372  void
373  read_custom_terrain(const std::string& fname,
374  const bool is_usgs,
375  const amrex::Geometry& geom,
376  amrex::FArrayBox& terrain_fab,
377  const amrex::Real& /*time*/)
378  {
379  amrex::Vector<amrex::Real> m_xterrain,m_yterrain,m_zterrain;
380 
381  int nx = 0; int ny = 0;
382 
383  if (amrex::ParallelDescriptor::IOProcessor()) {
384  // Read terrain file
385  amrex::Print()<<"Reading terrain file: "<< fname<< std::endl;
386  std::ifstream file(fname);
387 
388  if (!file.is_open()) {
389  amrex::Abort("Error: Could not open the file " + fname+ "\n");
390  }
391 
392  // Check if file is empty
393  if (file.peek() == std::ifstream::traits_type::eof()) {
394  amrex::Abort("Error: The file " + fname+ " is empty.\n");
395  }
396 
397  amrex::Real value1,value2,value3;
398 
399  if (is_usgs) {
400  amrex::Real lat_min, lon_min;
401 
402  file >> lon_min >> lat_min;
403  if(std::fabs(lon_min) > amrex::Real(180.0)) {
404  amrex::Error("The value of longitude for entry in the first line in " + fname
405  + " should not exceed amrex::Real(180.) It is " + std::to_string(lon_min));
406  }
407  if(std::fabs(lat_min) > amrex::Real(90.0)) {
408  amrex::Error("The value of latitude for entry in the first line in " + fname
409  + " should not exceed amrex::Real(90.) It is " + std::to_string(lat_min));
410  }
411 
412  file >> nx >> ny;
413 
414  int counter = 0;
415  while (file >> value1 >> value2 >> value3) {
416  m_xterrain.push_back(value1);
417  if(counter%nx==0) {
418  m_yterrain.push_back(value2);
419  }
420  m_zterrain.push_back(value3);
421  counter += 1;
422  }
423  AMREX_ASSERT(m_xterrain.size() == static_cast<long int>(nx*ny));
424  AMREX_ASSERT(m_yterrain.size() == static_cast<long int>(ny));
425  AMREX_ASSERT(m_zterrain.size() == static_cast<long int>(nx*ny));
426 
427  } else {
428  int cnt = 1;
429  nx = erf_get_single_value<int>(file,cnt); cnt++;
430  ny = erf_get_single_value<int>(file,cnt); cnt++;
431  amrex::Print()<<"Expecting " << nx << " values of x, " <<
432  ny << " values of y, and " <<
433  nx*ny << " values of z" << std::endl;
434  AMREX_ALWAYS_ASSERT(nx > 0);
435  AMREX_ALWAYS_ASSERT(ny > 0);
436  m_xterrain.resize(nx);
437  m_yterrain.resize(ny);
438  m_zterrain.resize(nx * ny);
439  for (int n = 0; n < nx; n++) {
440  m_xterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
441  cnt++;
442  }
443  for (int n = 0; n < ny; n++) {
444  m_yterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
445  cnt++;
446  }
447  for (int n = 0; n < nx * ny; n++) {
448  m_zterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
449  cnt++;
450  }
451  }
452 
453  // Close the file after reading
454  file.close();
455  }
456 
457  amrex::ParallelDescriptor::Bcast(&nx,1,amrex::ParallelDescriptor::IOProcessorNumber());
458  amrex::ParallelDescriptor::Bcast(&ny,1,amrex::ParallelDescriptor::IOProcessorNumber());
459 
460  int nz = nx * ny;
461 
462  m_xterrain.resize(nx);
463  m_yterrain.resize(ny);
464  m_zterrain.resize(nz);
465 
466  amrex::ParallelDescriptor::Bcast(m_xterrain.data(),nx,amrex::ParallelDescriptor::IOProcessorNumber());
467  amrex::ParallelDescriptor::Bcast(m_yterrain.data(),ny,amrex::ParallelDescriptor::IOProcessorNumber());
468  amrex::ParallelDescriptor::Bcast(m_zterrain.data(),nz,amrex::ParallelDescriptor::IOProcessorNumber());
469 
470  // Copy data to the GPU
471  amrex::Gpu::DeviceVector<amrex::Real> d_xterrain(nx),d_yterrain(ny),d_zterrain(nz);
472  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_xterrain.begin(), m_xterrain.end(), d_xterrain.begin());
473  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_yterrain.begin(), m_yterrain.end(), d_yterrain.begin());
474  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_zterrain.begin(), m_zterrain.end(), d_zterrain.begin());
475 
476  amrex::Real* d_xt = d_xterrain.data();
477  amrex::Real* d_yt = d_yterrain.data();
478  amrex::Real* d_zt = d_zterrain.data();
479 
480  auto dx = geom.CellSizeArray();
481  auto ProbLoArr = geom.ProbLoArray();
482 
483  int ilo = geom.Domain().smallEnd(0);
484  int jlo = geom.Domain().smallEnd(1);
485  int klo = geom.Domain().smallEnd(2);
486  int ihi = geom.Domain().bigEnd(0) + 1;
487  int jhi = geom.Domain().bigEnd(1) + 1;
488 
489  amrex::Box zbx = terrain_fab.box();
490  amrex::Array4<amrex::Real> const& z_arr = terrain_fab.array();
491 
492  amrex::ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/)
493  {
494  // Clip indices for ghost cells
495  int ii = amrex::min(amrex::max(i,ilo),ihi);
496  int jj = amrex::min(amrex::max(j,jlo),jhi);
497 
498  // Location of nodes
499  amrex::Real x = ProbLoArr[0] + ii * dx[0] + amrex::Real(1e-3);
500  amrex::Real y = ProbLoArr[1] + jj * dx[1] + amrex::Real(1e-3);
501 
502  int ind11, ind12, ind21, ind22;
503  amrex::Real x1, x2, y1, y2;
504 
505  int iindex_terrain=-1;
506  int jindex_terrain=-1;
507  // *******************************************************************
508  // NOTE: usgs-format is contiguous in x
509  // *******************************************************************
510  if (is_usgs) {
511 
512  for (int it = 0; it < ny && jindex_terrain == -1; it++) {
513  if (d_yt[it] > y) {
514  jindex_terrain = it-1;
515  }
516  }
517  if (jindex_terrain == -1) {
518  jindex_terrain = ny-1;
519  }
520 
521  int gstart = (jindex_terrain )*nx;
522  int gend = (jindex_terrain+1)*nx-1;
523  for (int it = gstart; it <= gend && iindex_terrain == -1; it++) {
524  if (d_xt[it] > x) {
525  iindex_terrain = it-gstart;
526  }
527  }
528 
529  // Define the four values to interpolate between
530  ind11 = jindex_terrain*nx + iindex_terrain; // (x1,y1)
531  ind12 = ind11+nx; // (x1,y2)
532  ind21 = ind11+1; // (x2,y1)
533  ind22 = ind12+1; // (x2,y2)
534 
535  x1 = d_xt[ind11];
536  x2 = d_xt[ind21];
537  y1 = d_yt[jindex_terrain];
538  y2 = d_yt[jindex_terrain+1];
539 
540  amrex::Real denom = (x2-x1)*(y2-y1);
541  amrex::Real w_11 = (x2-x)*(y2-y)/denom; // (x1,y1)
542  amrex::Real w_12 = (x2-x)*(y-y1)/denom; // (x1,y2)
543  amrex::Real w_21 = (x-x1)*(y2-y)/denom; // (x2,y1)
544  amrex::Real w_22 = (x-x1)*(y-y1)/denom; // (x2,y2)
545 
546  // Do bilinear interpolation
547  z_arr(i,j,klo) = w_11*d_zt[ind11] + w_12*d_zt[ind12] + w_21*d_zt[ind21] + w_22*d_zt[ind22];
548 
549  } else {
550 
551  for (int it = 0; it < ny && jindex_terrain == -1; it++) {
552  if (d_yt[it] > y) {
553  jindex_terrain = it-1;
554  }
555  }
556  if (jindex_terrain == -1) {
557  jindex_terrain = ny-1;
558  }
559 
560  for (int it = 0; it < nx && iindex_terrain == -1; it++) {
561  if (d_xt[it] > x) {
562  iindex_terrain = it-1;
563  }
564  }
565  if (iindex_terrain == -1) {
566  iindex_terrain = nx-1;
567  }
568 
569  // Define the four values to interpolate between
570  x1 = d_xt[iindex_terrain];
571  x2 = d_xt[iindex_terrain+1];
572  y1 = d_yt[jindex_terrain];
573  y2 = d_yt[jindex_terrain+1];
574 
575 #if 1
576  // *******************************************************************
577  // NOTE: this format is contiguous in y to match the AMR-Wind read
578  // *******************************************************************
579  ind11 = iindex_terrain * ny + jindex_terrain; // (x1,y1)
580  ind21 = std::min(iindex_terrain+1,nx-1) * ny + jindex_terrain; // (x2,y1)
581 
582  ind12 = iindex_terrain * ny + std::min(jindex_terrain+1,ny-1); // (x1,y2)
583  ind22 = std::min(iindex_terrain+1,nx-1) * ny + std::min(jindex_terrain+1,ny-1); // (x1,y2)
584 #else
585  // *******************************************************************
586  // NOTE: this format is contiguous in x as an alternative
587  // *******************************************************************
588 
589  ind11 = jindex_terrain * nx + iindex_terrain; // (x1,y1)
590  ind12 = std::min(jindex_terrain+1,ny-1) * nx + iindex_terrain; // (x1,y2)
591 
592  ind21 = jindex_terrain * nx + std::min(iindex_terrain+1,nx-1); // (x2,y1)
593  ind22 = std::min(jindex_terrain+1,ny-1) * nx + std::min(iindex_terrain+1,nx-1); // (x2,y2)
594 #endif
595 
596  if (iindex_terrain == nx-1 && jindex_terrain == ny-1)
597  {
598  z_arr(i,j,klo) = d_zt[ind11];
599  }
600  else if (iindex_terrain != nx-1 && jindex_terrain == ny-1)
601  {
602  amrex::Real w_11 = (x2-x); // (x1,y1)
603  amrex::Real w_21 = (x-x1); // (x2,y1)
604  amrex::Real denom = (x2-x1);
605  z_arr(i,j,klo) = (w_11*d_zt[ind11] + w_21*d_zt[ind21])/denom;
606  }
607  else if (iindex_terrain == nx-1 && jindex_terrain != ny-1)
608  {
609  amrex::Real w_11 = (y2-y); // (x1,y1)
610  amrex::Real w_12 = (y-y1); // (x1,y2)
611  amrex::Real denom = (y2-y1);
612  z_arr(i,j,klo) = (w_11*d_zt[ind11] + w_12*d_zt[ind12])/denom;
613  }
614  else
615  {
616  amrex::Real w_11 = (x2-x)*(y2-y); // (x1,y1)
617  amrex::Real w_21 = (x-x1)*(y2-y); // (x2,y1)
618  amrex::Real w_12 = (x2-x)*(y-y1); // (x1,y2)
619  amrex::Real w_22 = (x-x1)*(y-y1); // (x2,y2)
620  amrex::Real denom = (x2-x1)*(y2-y1);
621  z_arr(i,j,klo) = (w_11*d_zt[ind11] + w_12*d_zt[ind12] + w_21*d_zt[ind21] + w_22*d_zt[ind22]) / denom;
622  }
623  } // usgs?
624  });
625  }
626 
627  /**
628  * Function to perform custom initialization of terrain
629  *
630  * Note: Terrain functionality can also be used to provide grid stretching.
631  *
632  * @param[in] geom container for geometric information
633  * @param[out] terrain_mf height coordinate at nodes
634  * @param[in] time current time
635  */
636  virtual void
637  init_custom_terrain (const amrex::Geometry& geom,
638  amrex::FArrayBox& terrain_fab,
639  const amrex::Real& time)
640  {
641  std::string custom_terrain_type = "None";
642  amrex::ParmParse pp_prob("prob"); pp_prob.query("custom_terrain_type",custom_terrain_type);
643 
644  if (custom_terrain_type != "None")
645  {
646  amrex::Print() << "Calling custom terrain initialization" << std::endl;
647  init_my_custom_terrain(geom,terrain_fab,time);
648 
649  } else {
650  amrex::Print() << "Initializing flat terrain" << std::endl;
651  terrain_fab.template setVal<amrex::RunOn::Device>(0);
652 
653  if (SolverChoice::mesh_type == MeshType::VariableDz) {
654  SolverChoice::set_mesh_type(MeshType::StretchedDz);
655  amrex::Print() << "Resetting mesh type to StretchedDz" << std::endl;
656  }
657  }
658  }
659 
660 #ifdef ERF_USE_TERRAIN_VELOCITY
661  virtual amrex::Real compute_terrain_velocity(const amrex::Real /*time*/)
662  {
663  amrex::Error("Should never call compute_terrain_velocity for "+name()+" problem");
664  }
665 #endif
666 
667  /**
668  * Function to define the quantities needed to impose Rayleigh damping
669  *
670  * @param[out] rayleigh_ptrs = {strength of Rayleigh damping, reference values for xvel/yvel/zvel/theta used to define Rayleigh damping}
671  * @param[in] geom container for geometric information
672  * @param[in] z_phys_cc height coordinate at cell centers
673  */
674  virtual void
675  erf_init_rayleigh (amrex::Vector<amrex::Vector<amrex::Real> >& /*rayleigh_ptrs*/,
676  amrex::Geometry const& /*geom*/,
677  std::unique_ptr<amrex::MultiFab>& /*z_phys_nd*/,
678  amrex::Real /*zdamp*/)
679  {
680  // Default which does no harm
681  }
682 
683 protected:
684  // Struct to store problem parameters
686 
687  /**
688  * Function to update default base parameters, currently only used for
689  * init_type == InitType::Uniform
690  */
693  base_parms.T_0 = T_0;
694  }
695 
696  std::string name() {
697  std::string prob_name;
698  amrex::ParmParse pp("erf");
699  pp.get("prob_name",prob_name); return prob_name;
700  }
701 };
702 
703 
704 /**
705  * Function to init the physical bounds of the domain
706  * and instantiate a Problem derived from ProblemBase
707 */
708 extern std::unique_ptr<ProblemBase> amrex_probinit (const amrex_real* problo,
709  const amrex_real* probhi) AMREX_ATTRIBUTE_WEAK;
710 #endif
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
@ wbar
Definition: ERF_DataStruct.H:98
auto probhi
Definition: ERF_InitCustomPertVels_ABL.H:21
auto problo
Definition: ERF_InitCustomPertVels_ABL.H:20
Real rho_0
Definition: ERF_InitCustomPert_ABL.H:4
ParmParse pp("prob")
Real T_0
Definition: ERF_InitCustomPert_ABL.H:5
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
ParmParse pp_prob("prob")
Real T
Definition: ERF_InitCustomPert_Bubble.H:105
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
T erf_get_single_value(std::istream &is, int n)
Definition: ERF_ProbCommon.H:20
void init_my_custom_terrain(const amrex::Geometry &geom, amrex::FArrayBox &terrain_fab, const amrex::Real &time)
std::unique_ptr< ProblemBase > amrex_probinit(const amrex_real *problo, const amrex_real *probhi) AMREX_ATTRIBUTE_WEAK
const Box zbx
Definition: ERF_SetupDiff.H:9
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_ProbCommon.H:49
virtual void init_custom_terrain(const amrex::Geometry &geom, amrex::FArrayBox &terrain_fab, const amrex::Real &time)
Definition: ERF_ProbCommon.H:637
virtual void update_w_subsidence(const amrex::Real &, amrex::Vector< amrex::Real > &wbar, amrex::Gpu::DeviceVector< amrex::Real > &d_wbar, const amrex::MultiFab &, const amrex::Geometry &geom, std::unique_ptr< amrex::MultiFab > &)
Definition: ERF_ProbCommon.H:251
virtual void update_rhoqt_sources(const amrex::Real &, amrex::MultiFab *qsrc, const amrex::Geometry &, std::unique_ptr< amrex::MultiFab > &)
Definition: ERF_ProbCommon.H:216
virtual void erf_init_const_dens_and_linear_th_hse(amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::Real, std::unique_ptr< amrex::MultiFab > &)
Definition: ERF_ProbCommon.H:95
void init_terrain_surface(const amrex::Geometry &geom, amrex::FArrayBox &terrain_fab, const amrex::Real &time)
Definition: ERF_ProbCommon.H:347
virtual void init_custom_pert_vels(const amrex::Box &, const amrex::Box &, const amrex::Box &, amrex::Array4< amrex::Real > const &, amrex::Array4< amrex::Real > const &, amrex::Array4< amrex::Real > const &, amrex::Array4< amrex::Real const > const &, amrex::GeometryData const &, amrex::Array4< amrex::Real const > const &, amrex::Array4< amrex::Real const > const &, const SolverChoice &, const int)
Definition: ERF_ProbCommon.H:159
virtual void erf_init_const_dens_hse(amrex::MultiFab &)
Definition: ERF_ProbCommon.H:81
virtual void erf_init_dens_hse_dry(amrex::MultiFab &, std::unique_ptr< amrex::MultiFab > &, std::unique_ptr< amrex::MultiFab > &, amrex::Geometry const &, const amrex::Vector< amrex::Real > &, bool, bool)
Definition: ERF_ProbCommon.H:66
void read_custom_terrain(const std::string &fname, const bool is_usgs, const amrex::Geometry &geom, amrex::FArrayBox &terrain_fab, const amrex::Real &)
Definition: ERF_ProbCommon.H:373
void init_buildings_surface(const amrex::Geometry &geom, amrex::FArrayBox &buildings_fab, const amrex::Real &time)
Definition: ERF_ProbCommon.H:322
virtual void erf_init_rayleigh(amrex::Vector< amrex::Vector< amrex::Real > > &, amrex::Geometry const &, std::unique_ptr< amrex::MultiFab > &, amrex::Real)
Definition: ERF_ProbCommon.H:675
void init_base_parms(amrex::Real rho_0, amrex::Real T_0)
Definition: ERF_ProbCommon.H:691
virtual void update_rhotheta_sources(const amrex::Real &, amrex::MultiFab *src, const amrex::Geometry &, std::unique_ptr< amrex::MultiFab > &)
Definition: ERF_ProbCommon.H:184
virtual ~ProblemBase()=default
virtual void erf_init_dens_hse_moist(amrex::MultiFab &, std::unique_ptr< amrex::MultiFab > &, amrex::Geometry const &)
Definition: ERF_ProbCommon.H:104
virtual void init_custom_pert(const amrex::Box &, amrex::Array4< amrex::Real const > const &, amrex::Array4< amrex::Real > const &, amrex::Array4< amrex::Real > const &, amrex::Array4< amrex::Real > const &, amrex::Array4< amrex::Real const > const &, amrex::Array4< amrex::Real const > const &, amrex::GeometryData const &, amrex::Array4< amrex::Real const > const &, const SolverChoice &, const int)
Definition: ERF_ProbCommon.H:127
std::string name()
Definition: ERF_ProbCommon.H:696
virtual void erf_init_const_dens_and_th_hse(amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::Real)
Definition: ERF_ProbCommon.H:87
ProbParmDefaults base_parms
Definition: ERF_ProbCommon.H:685
virtual void update_geostrophic_profile(const amrex::Real &, amrex::Vector< amrex::Real > &u_geos, amrex::Gpu::DeviceVector< amrex::Real > &d_u_geos, amrex::Vector< amrex::Real > &v_geos, amrex::Gpu::DeviceVector< amrex::Real > &d_v_geos, const amrex::Geometry &geom, std::unique_ptr< amrex::MultiFab > &)
Definition: ERF_ProbCommon.H:286
Definition: ERF_ProbCommon.H:14
amrex::Real T_0
Definition: ERF_ProbCommon.H:16
amrex::Real rho_0
Definition: ERF_ProbCommon.H:15
Definition: ERF_DataStruct.H:141
static MeshType mesh_type
Definition: ERF_DataStruct.H:1134
static void set_mesh_type(MeshType new_mesh_type)
Definition: ERF_DataStruct.H:1137