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_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_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  }
std::string name()
Definition: ERF_ProbCommon.H:685
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
96  {
97 
98  }

◆ 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
668  {
669  // Default which does no harm
670  }

◆ 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

680  {
682  base_parms.T_0 = T_0;
683  }
amrex::Real rho_0
Definition: ERF_InitCustomPert_IsentropicVortex.H:33
Real T_0
Definition: ERF_InitCustomPert_TaylorGreenVortex.H:4
ProbParmDefaults base_parms
Definition: ERF_ProbCommon.H:674
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
314  {
315  // Check if a valid text file exists for the buildings
316  std::string fname;
317  amrex::ParmParse pp("erf");
318  auto valid_fname = pp.query("buildings_file_name",fname);
319 
320  if (valid_fname) {
321  read_custom_terrain(fname,false,geom,buildings_fab,time);
322  } else {
323  init_my_custom_terrain(geom, buildings_fab, time);
324  }
325  }
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:362
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
126  {
127  amrex::Print() << "No perturbation to background fields supplied for "
128  << name() << " problem" << std::endl;
129  }
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
159  {
160  amrex::Print() << "No perturbation velocities supplied for " << name() << " problem" << std::endl;
161  }
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
629  {
630  std::string custom_terrain_type = "None";
631  amrex::ParmParse pp_prob("prob"); pp_prob.query("custom_terrain_type",custom_terrain_type);
632 
633  if (custom_terrain_type != "None")
634  {
635  amrex::Print() << "Calling custom terrain initialization" << std::endl;
636  init_my_custom_terrain(geom,terrain_fab,time);
637 
638  } else {
639  amrex::Print() << "Initializing flat terrain" << std::endl;
640  terrain_fab.template setVal<amrex::RunOn::Device>(0.0);
641 
642  if (SolverChoice::mesh_type == MeshType::VariableDz) {
643  SolverChoice::set_mesh_type(MeshType::StretchedDz);
644  amrex::Print() << "Resetting mesh type to StretchedDz" << std::endl;
645  }
646  }
647  }
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
339  {
340  // Check if a valid text file exists for the terrain
341  std::string fname, fname_usgs;
342  amrex::ParmParse pp("erf");
343  auto valid_fname = pp.query("terrain_file_name",fname);
344  auto valid_fname_USGS = pp.query("terrain_file_name_USGS",fname_usgs);
345 
346  bool is_usgs;
347 
348  if (valid_fname) {
349  is_usgs = false;
350  read_custom_terrain(fname,is_usgs,geom,terrain_fab,time);
351 
352  } else if (valid_fname_USGS) {
353  is_usgs = true;
354  read_custom_terrain(fname_usgs,is_usgs,geom,terrain_fab,time);
355 
356  } else {
357  init_my_custom_terrain (geom, terrain_fab, time);
358  }
359  }
Here is the call graph for this function:

◆ name()

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

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