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 (amrex::MultiFab &, std::unique_ptr< amrex::MultiFab > &, std::unique_ptr< amrex::MultiFab > &, amrex::Geometry const &)
 
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
97  {
98  amrex::Error("Should never call this version of erf_init_const_dens_and_linear_th_hse for "+name()+" problem");
99  }
std::string name()
Definition: ERF_ProbCommon.H:694
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
88  {
89  amrex::Error("Should never call this version of erf_init_const_dens_and_th_hse for "+name()+" problem");
90  }
Here is the call graph for this function:

◆ erf_init_const_dens_hse()

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

◆ erf_init_dens_hse()

virtual void ProblemBase::erf_init_dens_hse ( amrex::MultiFab &  ,
std::unique_ptr< amrex::MultiFab > &  ,
std::unique_ptr< amrex::MultiFab > &  ,
amrex::Geometry const &   
)
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
70  {
71  amrex::Print() << "Hydrostatically balanced density was NOT set"
72  << " -- an appropriate init_type should probably have been specified"
73  << " (e.g., input_sounding, WRFInput, or Metgrid)"
74  << std::endl;
75  amrex::Error("Should never call this version of erf_init_dens_hse for "+name()+" problem");
76  }
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
105  {
106 
107  }

◆ 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
677  {
678  // Default which does no harm
679  }

◆ 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

689  {
691  base_parms.T_0 = T_0;
692  }
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:683
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
323  {
324  // Check if a valid text file exists for the buildings
325  std::string fname;
326  amrex::ParmParse pp("erf");
327  auto valid_fname = pp.query("buildings_file_name",fname);
328 
329  if (valid_fname) {
330  read_custom_terrain(fname,false,geom,buildings_fab,time);
331  } else {
332  init_my_custom_terrain(geom, buildings_fab, time);
333  }
334  }
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:371
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
135  {
136  amrex::Print() << "No perturbation to background fields supplied for "
137  << name() << " problem" << std::endl;
138  }
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
168  {
169  amrex::Print() << "No perturbation velocities supplied for " << name() << " problem" << std::endl;
170  }
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
638  {
639  std::string custom_terrain_type = "None";
640  amrex::ParmParse pp_prob("prob"); pp_prob.query("custom_terrain_type",custom_terrain_type);
641 
642  if (custom_terrain_type != "None")
643  {
644  amrex::Print() << "Calling custom terrain initialization" << std::endl;
645  init_my_custom_terrain(geom,terrain_fab,time);
646 
647  } else {
648  amrex::Print() << "Initializing flat terrain" << std::endl;
649  terrain_fab.template setVal<amrex::RunOn::Device>(0.0);
650 
651  if (SolverChoice::mesh_type == MeshType::VariableDz) {
652  SolverChoice::set_mesh_type(MeshType::StretchedDz);
653  amrex::Print() << "Resetting mesh type to StretchedDz" << std::endl;
654  }
655  }
656  }
ParmParse pp_prob("prob")
static MeshType mesh_type
Definition: ERF_DataStruct.H:1068
static void set_mesh_type(MeshType new_mesh_type)
Definition: ERF_DataStruct.H:1071
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
348  {
349  // Check if a valid text file exists for the terrain
350  std::string fname, fname_usgs;
351  amrex::ParmParse pp("erf");
352  auto valid_fname = pp.query("terrain_file_name",fname);
353  auto valid_fname_USGS = pp.query("terrain_file_name_USGS",fname_usgs);
354 
355  bool is_usgs;
356 
357  if (valid_fname) {
358  is_usgs = false;
359  read_custom_terrain(fname,is_usgs,geom,terrain_fab,time);
360 
361  } else if (valid_fname_USGS) {
362  is_usgs = true;
363  read_custom_terrain(fname_usgs,is_usgs,geom,terrain_fab,time);
364 
365  } else {
366  init_my_custom_terrain (geom, terrain_fab, time);
367  }
368  }
Here is the call graph for this function:

◆ name()

std::string ProblemBase::name ( )
inlineprotected
694 { std::string prob_name; amrex::ParmParse pp("erf"); pp.query("prob_name",prob_name); return prob_name; }

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(), 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
376  {
377  amrex::Vector<amrex::Real> m_xterrain,m_yterrain,m_zterrain;
378 
379  int nx = 0; int ny = 0;
380 
381  if (amrex::ParallelDescriptor::IOProcessor()) {
382  // Read terrain file
383  amrex::Print()<<"Reading terrain file: "<< fname<< std::endl;
384  std::ifstream file(fname);
385 
386  if (!file.is_open()) {
387  amrex::Abort("Error: Could not open the file " + fname+ "\n");
388  }
389 
390  // Check if file is empty
391  if (file.peek() == std::ifstream::traits_type::eof()) {
392  amrex::Abort("Error: The file " + fname+ " is empty.\n");
393  }
394 
395  amrex::Real value1,value2,value3;
396 
397  if (is_usgs) {
398  amrex::Real lat_min, lon_min;
399 
400  file >> lon_min >> lat_min;
401  if(std::fabs(lon_min) > 180.0) {
402  amrex::Error("The value of longitude for entry in the first line in " + fname
403  + " should not exceed 180. It is " + std::to_string(lon_min));
404  }
405  if(std::fabs(lat_min) > 90.0) {
406  amrex::Error("The value of latitude for entry in the first line in " + fname
407  + " should not exceed 90. It is " + std::to_string(lat_min));
408  }
409 
410  file >> nx >> ny;
411 
412  int counter = 0;
413  while (file >> value1 >> value2 >> value3) {
414  m_xterrain.push_back(value1);
415  if(counter%nx==0) {
416  m_yterrain.push_back(value2);
417  }
418  m_zterrain.push_back(value3);
419  counter += 1;
420  }
421  AMREX_ASSERT(m_xterrain.size() == static_cast<long int>(nx*ny));
422  AMREX_ASSERT(m_yterrain.size() == static_cast<long int>(ny));
423  AMREX_ASSERT(m_zterrain.size() == static_cast<long int>(nx*ny));
424 
425  } else {
426  int cnt = 1;
427  nx = erf_get_single_value<int>(file,cnt); cnt++;
428  ny = erf_get_single_value<int>(file,cnt); cnt++;
429  amrex::Print()<<"Expecting " << nx << " values of x, " <<
430  ny << " values of y, and " <<
431  nx*ny << " values of z" << std::endl;
432  AMREX_ALWAYS_ASSERT(nx > 0);
433  AMREX_ALWAYS_ASSERT(ny > 0);
434  m_xterrain.resize(nx);
435  m_yterrain.resize(ny);
436  m_zterrain.resize(nx * ny);
437  for (int n = 0; n < nx; n++) {
438  m_xterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
439  cnt++;
440  }
441  for (int n = 0; n < ny; n++) {
442  m_yterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
443  cnt++;
444  }
445  for (int n = 0; n < nx * ny; n++) {
446  m_zterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
447  cnt++;
448  }
449  }
450 
451  // Close the file after reading
452  file.close();
453  }
454 
455  amrex::ParallelDescriptor::Bcast(&nx,1,amrex::ParallelDescriptor::IOProcessorNumber());
456  amrex::ParallelDescriptor::Bcast(&ny,1,amrex::ParallelDescriptor::IOProcessorNumber());
457 
458  int nz = nx * ny;
459 
460  m_xterrain.resize(nx);
461  m_yterrain.resize(ny);
462  m_zterrain.resize(nz);
463 
464  amrex::ParallelDescriptor::Bcast(m_xterrain.data(),nx,amrex::ParallelDescriptor::IOProcessorNumber());
465  amrex::ParallelDescriptor::Bcast(m_yterrain.data(),ny,amrex::ParallelDescriptor::IOProcessorNumber());
466  amrex::ParallelDescriptor::Bcast(m_zterrain.data(),nz,amrex::ParallelDescriptor::IOProcessorNumber());
467 
468  // Copy data to the GPU
469  amrex::Gpu::DeviceVector<amrex::Real> d_xterrain(nx),d_yterrain(ny),d_zterrain(nz);
470  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_xterrain.begin(), m_xterrain.end(), d_xterrain.begin());
471  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_yterrain.begin(), m_yterrain.end(), d_yterrain.begin());
472  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_zterrain.begin(), m_zterrain.end(), d_zterrain.begin());
473 
474  amrex::Real* d_xt = d_xterrain.data();
475  amrex::Real* d_yt = d_yterrain.data();
476  amrex::Real* d_zt = d_zterrain.data();
477 
478  auto dx = geom.CellSizeArray();
479  auto ProbLoArr = geom.ProbLoArray();
480 
481  int ilo = geom.Domain().smallEnd(0);
482  int jlo = geom.Domain().smallEnd(1);
483  int klo = geom.Domain().smallEnd(2);
484  int ihi = geom.Domain().bigEnd(0) + 1;
485  int jhi = geom.Domain().bigEnd(1) + 1;
486 
487  amrex::Box zbx = terrain_fab.box();
488  amrex::Array4<amrex::Real> const& z_arr = terrain_fab.array();
489 
490  amrex::ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/)
491  {
492  // Clip indices for ghost cells
493  int ii = amrex::min(amrex::max(i,ilo),ihi);
494  int jj = amrex::min(amrex::max(j,jlo),jhi);
495 
496  // Location of nodes
497  amrex::Real x = ProbLoArr[0] + ii * dx[0] + 1e-3;
498  amrex::Real y = ProbLoArr[1] + jj * dx[1] + 1e-3;
499 
500  int ind11, ind12, ind21, ind22;
501  amrex::Real x1, x2, y1, y2;
502 
503  int iindex_terrain=-1;
504  int jindex_terrain=-1;
505  // *******************************************************************
506  // NOTE: usgs-format is contiguous in x
507  // *******************************************************************
508  if (is_usgs) {
509 
510  for (int it = 0; it < ny && jindex_terrain == -1; it++) {
511  if (d_yt[it] > y) {
512  jindex_terrain = it-1;
513  }
514  }
515  if (jindex_terrain == -1) {
516  jindex_terrain = ny-1;
517  }
518 
519  int gstart = (jindex_terrain )*nx;
520  int gend = (jindex_terrain+1)*nx-1;
521  for (int it = gstart; it <= gend && iindex_terrain == -1; it++) {
522  if (d_xt[it] > x) {
523  iindex_terrain = it-gstart;
524  }
525  }
526 
527  // Define the four values to interpolate between
528  ind11 = jindex_terrain*nx + iindex_terrain; // (x1,y1)
529  ind12 = ind11+nx; // (x1,y2)
530  ind21 = ind11+1; // (x2,y1)
531  ind22 = ind12+1; // (x2,y2)
532 
533  x1 = d_xt[ind11];
534  x2 = d_xt[ind21];
535  y1 = d_yt[jindex_terrain];
536  y2 = d_yt[jindex_terrain+1];
537 
538  amrex::Real denom = (x2-x1)*(y2-y1);
539  amrex::Real w_11 = (x2-x)*(y2-y)/denom; // (x1,y1)
540  amrex::Real w_12 = (x2-x)*(y-y1)/denom; // (x1,y2)
541  amrex::Real w_21 = (x-x1)*(y2-y)/denom; // (x2,y1)
542  amrex::Real w_22 = (x-x1)*(y-y1)/denom; // (x2,y2)
543 
544  // Do bilinear interpolation
545  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];
546 
547  } else {
548 
549  for (int it = 0; it < ny && jindex_terrain == -1; it++) {
550  if (d_yt[it] > y) {
551  jindex_terrain = it-1;
552  }
553  }
554  if (jindex_terrain == -1) {
555  jindex_terrain = ny-1;
556  }
557 
558  for (int it = 0; it < nx && iindex_terrain == -1; it++) {
559  if (d_xt[it] > x) {
560  iindex_terrain = it-1;
561  }
562  }
563  if (iindex_terrain == -1) {
564  iindex_terrain = nx-1;
565  }
566 
567  // Define the four values to interpolate between
568  x1 = d_xt[iindex_terrain];
569  x2 = d_xt[iindex_terrain+1];
570  y1 = d_yt[jindex_terrain];
571  y2 = d_yt[jindex_terrain+1];
572 
573 #if 1
574  // *******************************************************************
575  // NOTE: this format is contiguous in y to match the AMR-Wind read
576  // *******************************************************************
577  ind11 = iindex_terrain * ny + jindex_terrain; // (x1,y1)
578  ind21 = std::min(iindex_terrain+1,nx-1) * ny + jindex_terrain; // (x2,y1)
579 
580  ind12 = iindex_terrain * ny + std::min(jindex_terrain+1,ny-1); // (x1,y2)
581  ind22 = std::min(iindex_terrain+1,nx-1) * ny + std::min(jindex_terrain+1,ny-1); // (x1,y2)
582 #else
583  // *******************************************************************
584  // NOTE: this format is contiguous in x as an alternative
585  // *******************************************************************
586 
587  ind11 = jindex_terrain * nx + iindex_terrain; // (x1,y1)
588  ind12 = std::min(jindex_terrain+1,ny-1) * nx + iindex_terrain; // (x1,y2)
589 
590  ind21 = jindex_terrain * nx + std::min(iindex_terrain+1,nx-1); // (x2,y1)
591  ind22 = std::min(jindex_terrain+1,ny-1) * nx + std::min(iindex_terrain+1,nx-1); // (x2,y2)
592 #endif
593 
594  if (iindex_terrain == nx-1 && jindex_terrain == ny-1)
595  {
596  z_arr(i,j,klo) = d_zt[ind11];
597  }
598  else if (iindex_terrain != nx-1 && jindex_terrain == ny-1)
599  {
600  amrex::Real w_11 = (x2-x); // (x1,y1)
601  amrex::Real w_21 = (x-x1); // (x2,y1)
602  amrex::Real denom = (x2-x1);
603  z_arr(i,j,klo) = (w_11*d_zt[ind11] + w_21*d_zt[ind21])/denom;
604  }
605  else if (iindex_terrain == nx-1 && jindex_terrain != ny-1)
606  {
607  amrex::Real w_11 = (y2-y); // (x1,y1)
608  amrex::Real w_12 = (y-y1); // (x1,y2)
609  amrex::Real denom = (y2-y1);
610  z_arr(i,j,klo) = (w_11*d_zt[ind11] + w_12*d_zt[ind12])/denom;
611  }
612  else
613  {
614  amrex::Real w_11 = (x2-x)*(y2-y); // (x1,y1)
615  amrex::Real w_21 = (x-x1)*(y2-y); // (x2,y1)
616  amrex::Real w_12 = (x2-x)*(y-y1); // (x1,y2)
617  amrex::Real w_22 = (x-x1)*(y-y1); // (x2,y2)
618  amrex::Real denom = (x2-x1)*(y2-y1);
619  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;
620  }
621  } // usgs?
622  });
623  }
const int klo
Definition: ERF_InitCustomPertVels_ParticleTests.H:4
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+0.5) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=1.0) { Real dT=T_pert *(std::cos(PI *L)+1.0)/2.0;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
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
291  {
292  if (u_geos.empty()) return;
293 
294  amrex::Warning("Geostrophic wind profile not defined for "+name()+" problem");
295 
296  const int khi = geom.Domain().bigEnd()[2];
297  // const amrex::Real* prob_lo = geom.ProbLo();
298  // const auto dx = geom.CellSize();
299  for (int k = 0; k <= khi; k++)
300  {
301  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
302  // set RHS term of RhoTheta equation based on time, z_cc here
303  u_geos[k] = 0.0;
304  v_geos[k] = 0.0;
305  }
306 
307  // Copy from host version to device version
308  amrex::Gpu::copy(amrex::Gpu::hostToDevice, u_geos.begin(), u_geos.end(), d_u_geos.begin());
309  amrex::Gpu::copy(amrex::Gpu::hostToDevice, v_geos.begin(), v_geos.end(), d_v_geos.begin());
310  }
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
218  {
219  if (qsrc->empty()) return;
220 
221  amrex::Warning("Moisture forcing not defined for "+name()+" problem");
222  for ( amrex::MFIter mfi(*qsrc, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
223  {
224  const auto &box = mfi.tilebox();
225  const amrex::Array4<amrex::Real>& qsrc_arr = qsrc->array(mfi);
226  // src is a spatial function if erf.spatial_moisture_forcing = true
227  // otherwise, qsrc_arr is defined only in over Z (box x and y dimensions are 1)
228  ParallelFor(box, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
229  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
230  // set RHS term of RhoQ1 equation based on time, z_cc here
231  qsrc_arr(i, j, k) = 0.0;
232  });
233  }
234  }
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
186  {
187  if (src->empty()) return;
188 
189  amrex::Warning("Temperature forcing not defined for "+name()+" problem");
190  for ( amrex::MFIter mfi(*src, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
191  {
192  const auto &box = mfi.tilebox();
193  const amrex::Array4<amrex::Real>& src_arr = src->array(mfi);
194  // src is a spatial function if erf.spatial_rhotheta_forcing = true
195  // otherwise, qsrc_arr is defined only in over Z (box x and y dimensions are 1)
196  ParallelFor(box, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
197  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
198  // set RHS term of RhoTheta equation based on time, z_cc here
199  src_arr(i, j, k) = 0.0;
200  });
201  }
202  }
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
255  {
256  if (wbar.empty()) return;
257 
258  amrex::Warning("Moisture forcing not defined for "+name()+" problem");
259 
260  const int khi = geom.Domain().bigEnd()[2];
261  // const amrex::Real* prob_lo = geom.ProbLo();
262  // const auto dx = geom.CellSize();
263  for (int k = 0; k <= khi; k++)
264  {
265  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
266  // set vertical velocity profile based on time, z_cc here
267  wbar[k] = 0.0;
268  }
269 
270  // Copy from host version to device version
271  amrex::Gpu::copy(amrex::Gpu::hostToDevice, wbar.begin(), wbar.end(), d_wbar.begin());
272  }
@ 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: