ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ProblemBase Class Reference

#include <ERF_ProbCommon.H>

Collaboration diagram for ProblemBase:

Public Member Functions

virtual ~ProblemBase ()=default
 
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)
 
virtual void erf_init_const_dens_hse (amrex::MultiFab &)
 
virtual void erf_init_const_dens_and_th_hse (amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::Real)
 
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 > &)
 
virtual void erf_init_dens_hse_moist (amrex::MultiFab &, std::unique_ptr< amrex::MultiFab > &, amrex::Geometry const &)
 
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)
 
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)
 
virtual void update_rhotheta_sources (const amrex::Real &, amrex::MultiFab *src, const amrex::Geometry &, std::unique_ptr< amrex::MultiFab > &)
 
virtual void update_rhoqt_sources (const amrex::Real &, amrex::MultiFab *qsrc, const amrex::Geometry &, std::unique_ptr< amrex::MultiFab > &)
 
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 > &)
 
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 > &)
 
void init_buildings_surface (const amrex::Geometry &geom, amrex::FArrayBox &buildings_fab, const amrex::Real &time)
 
void init_terrain_surface (const amrex::Geometry &geom, amrex::FArrayBox &terrain_fab, const amrex::Real &time)
 
void read_custom_terrain (const std::string &fname, const bool is_usgs, const amrex::Geometry &geom, amrex::FArrayBox &terrain_fab, const amrex::Real &)
 
virtual void init_custom_terrain (const amrex::Geometry &geom, amrex::FArrayBox &terrain_fab, const amrex::Real &time)
 
virtual void erf_init_rayleigh (amrex::Vector< amrex::Vector< amrex::Real > > &, amrex::Geometry const &, std::unique_ptr< amrex::MultiFab > &, amrex::Real)
 

Protected Member Functions

void init_base_parms (amrex::Real rho_0, amrex::Real T_0)
 
std::string name ()
 

Protected Attributes

ProbParmDefaults base_parms
 

Detailed Description

Class to hold problem-specific routines

Constructor & Destructor Documentation

◆ ~ProblemBase()

virtual ProblemBase::~ProblemBase ( )
virtualdefault

Virtual destructor to avoid data leakage with derived class

Member Function Documentation

◆ erf_init_const_dens_and_linear_th_hse()

virtual void ProblemBase::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 > &   
)
inlinevirtual
99  {
100  amrex::Error("Should never call this version of erf_init_const_dens_and_linear_th_hse for "+name()+" problem");
101  }
std::string name()
Definition: ERF_ProbCommon.H:696
Here is the call graph for this function:

◆ erf_init_const_dens_and_th_hse()

virtual void ProblemBase::erf_init_const_dens_and_th_hse ( amrex::MultiFab &  ,
amrex::MultiFab &  ,
amrex::MultiFab &  ,
amrex::MultiFab &  ,
amrex::MultiFab &  ,
amrex::Real   
)
inlinevirtual
90  {
91  amrex::Error("Should never call this version of erf_init_const_dens_and_th_hse for "+name()+" problem");
92  }
Here is the call graph for this function:

◆ erf_init_const_dens_hse()

virtual void ProblemBase::erf_init_const_dens_hse ( amrex::MultiFab &  )
inlinevirtual
82  {
83  amrex::Error("Should never call this version of erf_init_const_dens_hse for "+name()+" problem");
84  }
Here is the call graph for this function:

◆ erf_init_dens_hse_dry()

virtual void ProblemBase::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   
)
inlinevirtual

Function to initialize the hydrostatic reference density

Parameters
[out]rho_hsehydrostatic reference density
[in]z_phys_ndheight coordinate at nodes
[in]z_phys_ccheight coordinate at cell centers
[in]geomcontainer for geometric information
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  }
Here is the call graph for this function:

◆ erf_init_dens_hse_moist()

virtual void ProblemBase::erf_init_dens_hse_moist ( amrex::MultiFab &  ,
std::unique_ptr< amrex::MultiFab > &  ,
amrex::Geometry const &   
)
inlinevirtual
107  {
108 
109  }

◆ erf_init_rayleigh()

virtual void ProblemBase::erf_init_rayleigh ( amrex::Vector< amrex::Vector< amrex::Real > > &  ,
amrex::Geometry const &  ,
std::unique_ptr< amrex::MultiFab > &  ,
amrex::Real   
)
inlinevirtual

Function to define the quantities needed to impose Rayleigh damping

Parameters
[out]rayleigh_ptrs= {strength of Rayleigh damping, reference values for xvel/yvel/zvel/theta used to define Rayleigh damping}
[in]geomcontainer for geometric information
[in]z_phys_ccheight coordinate at cell centers
679  {
680  // Default which does no harm
681  }

◆ init_base_parms()

void ProblemBase::init_base_parms ( amrex::Real  rho_0,
amrex::Real  T_0 
)
inlineprotected

Function to update default base parameters, currently only used for init_type == InitType::Uniform

691  {
693  base_parms.T_0 = T_0;
694  }
Real rho_0
Definition: ERF_InitCustomPert_ABL.H:4
Real T_0
Definition: ERF_InitCustomPert_ABL.H:5
ProbParmDefaults base_parms
Definition: ERF_ProbCommon.H:685
amrex::Real T_0
Definition: ERF_ProbCommon.H:16
amrex::Real rho_0
Definition: ERF_ProbCommon.H:15

◆ init_buildings_surface()

void ProblemBase::init_buildings_surface ( const amrex::Geometry &  geom,
amrex::FArrayBox &  buildings_fab,
const amrex::Real time 
)
inline

Function to perform custom initialization of buildings

Parameters
[in]geomcontainer for geometric information
[out]z_phys_ndheight coordinate at nodes
[in]timecurrent 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  }
ParmParse pp("prob")
void init_my_custom_terrain(const amrex::Geometry &geom, amrex::FArrayBox &terrain_fab, const amrex::Real &time)
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
Here is the call graph for this function:

◆ init_custom_pert()

virtual void ProblemBase::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   
)
inlinevirtual

Function to perform custom initialization of a test problem

Parameters
[in]bxcell-centered box on which to initialize scalars
[out]statecell-centered variables to be filled in this routine
[out]r_hsehydrostatic reference density
[out]p_hsehydrostatic reference pressure
[in]z_ndheight coordinate at nodes
[in]z_ccheight coordinate at cell centers
[in]mf_mmap factor on cell centers
[in]mf_umap factor on x-faces
[in]mf_vmap factor on y-faces
[in]scSolverChoice structure that carries parameters
[in]levAMR level
137  {
138  amrex::Print() << "No perturbation to background fields supplied for "
139  << name() << " problem" << std::endl;
140  }
Here is the call graph for this function:

◆ init_custom_pert_vels()

virtual void ProblemBase::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   
)
inlinevirtual

Function to perform custom initialization of the velocities

Parameters
[in]xbxbox on which to initialize x_vel_pert
[in]ybxbox on which to initialize y_vel_pert
[in]zbxbox on which to initialize z_vel_pert
[out]x_vel_pertx-component of velocity perturbation to be filled in this routine
[out]y_vel_perty-component of velocity perturbation to be filled in this routine
[out]z_vel_pertz-component of velocity perturbation to be filled in this routine
[in]z_ndheight coordinate at nodes
[in]mf_mmap factor on cell centers
[in]mf_umap factor on x-faces
[in]mf_vmap factor on y-faces
[in]scSolverChoice structure that carries parameters
[in]levAMR level
170  {
171  amrex::Print() << "No perturbation velocities supplied for " << name() << " problem" << std::endl;
172  }
Here is the call graph for this function:

◆ init_custom_terrain()

virtual void ProblemBase::init_custom_terrain ( const amrex::Geometry &  geom,
amrex::FArrayBox &  terrain_fab,
const amrex::Real time 
)
inlinevirtual

Function to perform custom initialization of terrain

Note: Terrain functionality can also be used to provide grid stretching.

Parameters
[in]geomcontainer for geometric information
[out]terrain_mfheight coordinate at nodes
[in]timecurrent 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  }
ParmParse pp_prob("prob")
static MeshType mesh_type
Definition: ERF_DataStruct.H:1134
static void set_mesh_type(MeshType new_mesh_type)
Definition: ERF_DataStruct.H:1137
Here is the call graph for this function:

◆ init_terrain_surface()

void ProblemBase::init_terrain_surface ( const amrex::Geometry &  geom,
amrex::FArrayBox &  terrain_fab,
const amrex::Real time 
)
inline

Function to perform custom initialization of terrain

Parameters
[in]geomcontainer for geometric information
[out]z_phys_ndheight coordinate at nodes
[in]timecurrent 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  }
Here is the call graph for this function:

◆ name()

std::string ProblemBase::name ( )
inlineprotected
696  {
697  std::string prob_name;
698  amrex::ParmParse pp("erf");
699  pp.get("prob_name",prob_name); return prob_name;
700  }

Referenced by erf_init_const_dens_and_linear_th_hse(), erf_init_const_dens_and_th_hse(), erf_init_const_dens_hse(), erf_init_dens_hse_dry(), init_custom_pert(), init_custom_pert_vels(), update_geostrophic_profile(), update_rhoqt_sources(), update_rhotheta_sources(), and update_w_subsidence().

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

◆ read_custom_terrain()

void ProblemBase::read_custom_terrain ( const std::string &  fname,
const bool  is_usgs,
const amrex::Geometry &  geom,
amrex::FArrayBox &  terrain_fab,
const amrex::Real  
)
inline
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  }
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
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);})
const Box zbx
Definition: ERF_SetupDiff.H:9
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by init_buildings_surface(), and init_terrain_surface().

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

◆ update_geostrophic_profile()

virtual void ProblemBase::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 > &   
)
inlinevirtual

Function to update user-specified geostrophic wind profile.

Parameters
[in]timecurrent time
[out]u_geosgeostrophic wind profile
[out]v_geosgeostrophic wind profile
[in]geomcontainer for geometric information
[in]z_phys_ccheight coordinate at cell centers
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  }
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
Here is the call graph for this function:

◆ update_rhoqt_sources()

virtual void ProblemBase::update_rhoqt_sources ( const amrex::Real ,
amrex::MultiFab *  qsrc,
const amrex::Geometry &  ,
std::unique_ptr< amrex::MultiFab > &   
)
inlinevirtual

Function to update user-specified moisture source terms that can vary with time and height.

Parameters
[in]timecurrent time
[out]rhoqt_sourcemoisture forcing profile
[in]geomcontainer for geometric information
[in]z_phys_ccheight coordinate at cell centers
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  }
Here is the call graph for this function:

◆ update_rhotheta_sources()

virtual void ProblemBase::update_rhotheta_sources ( const amrex::Real ,
amrex::MultiFab *  src,
const amrex::Geometry &  ,
std::unique_ptr< amrex::MultiFab > &   
)
inlinevirtual

Function to update user-specified temperature source terms that can vary with time and height.

Parameters
[in]timecurrent time
[out]rhotheta_sourceforcing profile
[in]geomcontainer for geometric information
[in]z_phys_ccheight coordinate at cell centers
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  }
Here is the call graph for this function:

◆ update_w_subsidence()

virtual void ProblemBase::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 > &   
)
inlinevirtual

Function to update the vertical velocity profile, used to add subsidence source terms for x-mom, y-mom, rho*theta, rho*Q1, and rho*Q2.

TODO: Currently, this is only called by InitData, so there is no time dependence.

Parameters
[in]timecurrent time
[out]wbarw vel forcing profile
[in]geomcontainer for geometric information
[in]z_phys_ccheight coordinate at cell centers
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  }
@ wbar
Definition: ERF_DataStruct.H:98
Here is the call graph for this function:

Member Data Documentation

◆ base_parms

ProbParmDefaults ProblemBase::base_parms
protected

Referenced by init_base_parms().


The documentation for this class was generated from the following file: