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_dens_hse_moist (amrex::MultiFab &, std::unique_ptr< amrex::MultiFab > &, amrex::Geometry const &)
 
virtual void init_custom_pert (const amrex::Box &, const amrex::Box &, const amrex::Box &, 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 &, 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 &, 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)
 
void init_uniform (const amrex::Box &bx, amrex::Array4< amrex::Real > const &state)
 

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_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 erf_init_dens_hse for "+name()+" problem");
76  }
std::string name()
Definition: ERF_ProbCommon.H:675
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
82  {
83 
84  }

◆ 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
639  {
640  amrex::Error("Should never call erf_init_rayleigh for "+name()+" problem");
641  }
Here is the call graph for this function:

◆ 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

670  {
671  base_parms.rho_0 = rho_0;
672  base_parms.T_0 = T_0;
673  }
ProbParmDefaults base_parms
Definition: ERF_ProbCommon.H:664
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
285  {
286  // Check if a valid text file exists for the buildings
287  std::string fname;
288  amrex::ParmParse pp("erf");
289  auto valid_fname = pp.query("buildings_file_name",fname);
290 
291  if (valid_fname) {
292  read_custom_terrain(fname,false,geom,buildings_fab,time);
293  } else {
294  init_my_custom_terrain(geom, buildings_fab, time);
295  }
296  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:233
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:333
Here is the call graph for this function:

◆ init_custom_pert()

virtual void ProblemBase::init_custom_pert ( const amrex::Box &  ,
const amrex::Box &  ,
const amrex::Box &  ,
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 &  ,
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 &  ,
amrex::Array4< amrex::Real const > 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
[in]xbxbox on which to initialize x-component of velocity
[in]ybxbox on which to initialize y-component of velocity
[in]zbxbox on which to initialize z-component of velocity
[out]statecell-centered variables to be filled in this routine
[out]x_velocityx-component of velocity to be filled in this routine
[out]y_velocityy-component of velocity to be filled in this routine
[out]z_velocityz-component of velocity 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]qvwater vapor
[in]qccloud water
[in]qicloud ice
[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
129  {
130  amrex::Print() << "No perturbation to background fields supplied for "
131  << name() << " problem" << std::endl;
132  }
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
600  {
601  std::string custom_terrain_type = "None";
602  amrex::ParmParse pp_prob("prob"); pp_prob.query("custom_terrain_type",custom_terrain_type);
603 
604  if (custom_terrain_type != "None")
605  {
606  amrex::Print() << "Calling custom terrain initialization" << std::endl;
607  init_my_custom_terrain(geom,terrain_fab,time);
608 
609  } else {
610  amrex::Print() << "Initializing flat terrain" << std::endl;
611  terrain_fab.template setVal<amrex::RunOn::Device>(0.0);
612 
613  if (SolverChoice::mesh_type == MeshType::VariableDz) {
614  SolverChoice::set_mesh_type(MeshType::StretchedDz);
615  amrex::Print() << "Resetting mesh type to StretchedDz" << std::endl;
616  }
617  }
618  }
static MeshType mesh_type
Definition: ERF_DataStruct.H:1048
static void set_mesh_type(MeshType new_mesh_type)
Definition: ERF_DataStruct.H:1051
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
310  {
311  // Check if a valid text file exists for the terrain
312  std::string fname, fname_usgs;
313  amrex::ParmParse pp("erf");
314  auto valid_fname = pp.query("terrain_file_name",fname);
315  auto valid_fname_USGS = pp.query("terrain_file_name_USGS",fname_usgs);
316 
317  bool is_usgs;
318 
319  if (valid_fname) {
320  is_usgs = false;
321  read_custom_terrain(fname,is_usgs,geom,terrain_fab,time);
322 
323  } else if (valid_fname_USGS) {
324  is_usgs = true;
325  read_custom_terrain(fname_usgs,is_usgs,geom,terrain_fab,time);
326 
327  } else {
328  init_my_custom_terrain (geom, terrain_fab, time);
329  }
330  }
Here is the call graph for this function:

◆ init_uniform()

void ProblemBase::init_uniform ( const amrex::Box &  bx,
amrex::Array4< amrex::Real > const &  state 
)
inline

Function to set uniform background density and potential temperature fields

648  {
649  amrex::Real rho_0 = base_parms.rho_0;
650  amrex::Real T_0 = base_parms.T_0;
651  amrex::Print() << "Initializing uniform fields"
652  << " rho=" << rho_0 << " theta=" << T_0
653  << std::endl;
654 
655  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
656  {
657  state(i, j, k, Rho_comp) = rho_0;
658  state(i, j, k, RhoTheta_comp) = rho_0 * T_0;
659  });
660  }
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
amrex::Real Real
Definition: ERF_ShocInterface.H:19

◆ name()

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

Referenced by erf_init_dens_hse(), erf_init_rayleigh(), init_custom_pert(), 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
338  {
339  amrex::Vector<amrex::Real> m_xterrain,m_yterrain,m_zterrain;
340 
341  int nx = 0; int ny = 0;
342 
343  if (amrex::ParallelDescriptor::IOProcessor()) {
344  // Read terrain file
345  amrex::Print()<<"Reading terrain file: "<< fname<< std::endl;
346  std::ifstream file(fname);
347 
348  if (!file.is_open()) {
349  amrex::Abort("Error: Could not open the file " + fname+ "\n");
350  }
351 
352  // Check if file is empty
353  if (file.peek() == std::ifstream::traits_type::eof()) {
354  amrex::Abort("Error: The file " + fname+ " is empty.\n");
355  }
356 
357  amrex::Real value1,value2,value3;
358 
359  if (is_usgs) {
360  amrex::Real lat_min, lon_min;
361 
362  file >> lon_min >> lat_min;
363  if(std::fabs(lon_min) > 180.0) {
364  amrex::Error("The value of longitude for entry in the first line in " + fname
365  + " should not exceed 180. It is " + std::to_string(lon_min));
366  }
367  if(std::fabs(lat_min) > 90.0) {
368  amrex::Error("The value of latitude for entry in the first line in " + fname
369  + " should not exceed 90. It is " + std::to_string(lat_min));
370  }
371 
372  file >> nx >> ny;
373 
374  int counter = 0;
375  while (file >> value1 >> value2 >> value3) {
376  m_xterrain.push_back(value1);
377  if(counter%nx==0) {
378  m_yterrain.push_back(value2);
379  }
380  m_zterrain.push_back(value3);
381  counter += 1;
382  }
383  AMREX_ASSERT(m_xterrain.size() == static_cast<long int>(nx*ny));
384  AMREX_ASSERT(m_yterrain.size() == static_cast<long int>(ny));
385  AMREX_ASSERT(m_zterrain.size() == static_cast<long int>(nx*ny));
386 
387  } else {
388  int cnt = 1;
389  nx = erf_get_single_value<int>(file,cnt); cnt++;
390  ny = erf_get_single_value<int>(file,cnt); cnt++;
391  amrex::Print()<<"Expecting " << nx << " values of x, " <<
392  ny << " values of y, and " <<
393  nx*ny << " values of z" << std::endl;
394  AMREX_ALWAYS_ASSERT(nx > 0);
395  AMREX_ALWAYS_ASSERT(ny > 0);
396  m_xterrain.resize(nx);
397  m_yterrain.resize(ny);
398  m_zterrain.resize(nx * ny);
399  for (int n = 0; n < nx; n++) {
400  m_xterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
401  cnt++;
402  }
403  for (int n = 0; n < ny; n++) {
404  m_yterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
405  cnt++;
406  }
407  for (int n = 0; n < nx * ny; n++) {
408  m_zterrain[n] = erf_get_single_value<amrex::Real>(file,cnt);
409  cnt++;
410  }
411  }
412 
413  // Close the file after reading
414  file.close();
415  }
416 
417  amrex::ParallelDescriptor::Bcast(&nx,1,amrex::ParallelDescriptor::IOProcessorNumber());
418  amrex::ParallelDescriptor::Bcast(&ny,1,amrex::ParallelDescriptor::IOProcessorNumber());
419 
420  int nz = nx * ny;
421 
422  m_xterrain.resize(nx);
423  m_yterrain.resize(ny);
424  m_zterrain.resize(nz);
425 
426  amrex::ParallelDescriptor::Bcast(m_xterrain.data(),nx,amrex::ParallelDescriptor::IOProcessorNumber());
427  amrex::ParallelDescriptor::Bcast(m_yterrain.data(),ny,amrex::ParallelDescriptor::IOProcessorNumber());
428  amrex::ParallelDescriptor::Bcast(m_zterrain.data(),nz,amrex::ParallelDescriptor::IOProcessorNumber());
429 
430  // Copy data to the GPU
431  amrex::Gpu::DeviceVector<amrex::Real> d_xterrain(nx),d_yterrain(ny),d_zterrain(nz);
432  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_xterrain.begin(), m_xterrain.end(), d_xterrain.begin());
433  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_yterrain.begin(), m_yterrain.end(), d_yterrain.begin());
434  amrex::Gpu::copy(amrex::Gpu::hostToDevice, m_zterrain.begin(), m_zterrain.end(), d_zterrain.begin());
435 
436  amrex::Real* d_xt = d_xterrain.data();
437  amrex::Real* d_yt = d_yterrain.data();
438  amrex::Real* d_zt = d_zterrain.data();
439 
440  auto dx = geom.CellSizeArray();
441  auto ProbLoArr = geom.ProbLoArray();
442 
443  int ilo = geom.Domain().smallEnd(0);
444  int jlo = geom.Domain().smallEnd(1);
445  int klo = geom.Domain().smallEnd(2);
446  int ihi = geom.Domain().bigEnd(0) + 1;
447  int jhi = geom.Domain().bigEnd(1) + 1;
448 
449  amrex::Box zbx = terrain_fab.box();
450  amrex::Array4<amrex::Real> const& z_arr = terrain_fab.array();
451 
452  amrex::ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/)
453  {
454  // Clip indices for ghost cells
455  int ii = amrex::min(amrex::max(i,ilo),ihi);
456  int jj = amrex::min(amrex::max(j,jlo),jhi);
457 
458  // Location of nodes
459  amrex::Real x = ProbLoArr[0] + ii * dx[0] + 1e-3;
460  amrex::Real y = ProbLoArr[1] + jj * dx[1] + 1e-3;
461 
462  int ind11, ind12, ind21, ind22;
463  amrex::Real x1, x2, y1, y2;
464 
465  int iindex_terrain=-1;
466  int jindex_terrain=-1;
467  // *******************************************************************
468  // NOTE: usgs-format is contiguous in x
469  // *******************************************************************
470  if (is_usgs) {
471 
472  for (int it = 0; it < ny && jindex_terrain == -1; it++) {
473  if (d_yt[it] > y) {
474  jindex_terrain = it-1;
475  }
476  }
477  if (jindex_terrain == -1) {
478  jindex_terrain = ny-1;
479  }
480 
481  int gstart = (jindex_terrain )*nx;
482  int gend = (jindex_terrain+1)*nx-1;
483  for (int it = gstart; it <= gend && iindex_terrain == -1; it++) {
484  if (d_xt[it] > x) {
485  iindex_terrain = it-gstart;
486  }
487  }
488 
489  // Define the four values to interpolate between
490  ind11 = jindex_terrain*nx + iindex_terrain; // (x1,y1)
491  ind12 = ind11+nx; // (x1,y2)
492  ind21 = ind11+1; // (x2,y1)
493  ind22 = ind12+1; // (x2,y2)
494 
495  x1 = d_xt[ind11];
496  x2 = d_xt[ind21];
497  y1 = d_yt[jindex_terrain];
498  y2 = d_yt[jindex_terrain+1];
499 
500  amrex::Real denom = (x2-x1)*(y2-y1);
501  amrex::Real w_11 = (x2-x)*(y2-y)/denom; // (x1,y1)
502  amrex::Real w_12 = (x2-x)*(y-y1)/denom; // (x1,y2)
503  amrex::Real w_21 = (x-x1)*(y2-y)/denom; // (x2,y1)
504  amrex::Real w_22 = (x-x1)*(y-y1)/denom; // (x2,y2)
505 
506  // Do bilinear interpolation
507  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];
508 
509  } else {
510 
511  for (int it = 0; it < ny && jindex_terrain == -1; it++) {
512  if (d_yt[it] > y) {
513  jindex_terrain = it-1;
514  }
515  }
516  if (jindex_terrain == -1) {
517  jindex_terrain = ny-1;
518  }
519 
520  for (int it = 0; it < nx && iindex_terrain == -1; it++) {
521  if (d_xt[it] > x) {
522  iindex_terrain = it-1;
523  }
524  }
525  if (iindex_terrain == -1) {
526  iindex_terrain = nx-1;
527  }
528 
529  // Define the four values to interpolate between
530  x1 = d_xt[iindex_terrain];
531  x2 = d_xt[iindex_terrain+1];
532  y1 = d_yt[jindex_terrain];
533  y2 = d_yt[jindex_terrain+1];
534 
535 #if 1
536  // *******************************************************************
537  // NOTE: this format is contiguous in y to match the AMR-Wind read
538  // *******************************************************************
539  ind11 = iindex_terrain * ny + jindex_terrain; // (x1,y1)
540  ind21 = std::min(iindex_terrain+1,nx-1) * ny + jindex_terrain; // (x2,y1)
541 
542  ind12 = iindex_terrain * ny + std::min(jindex_terrain+1,ny-1); // (x1,y2)
543  ind22 = std::min(iindex_terrain+1,nx-1) * ny + std::min(jindex_terrain+1,ny-1); // (x1,y2)
544 #else
545  // *******************************************************************
546  // NOTE: this format is contiguous in x as an alternative
547  // *******************************************************************
548 
549  ind11 = jindex_terrain * nx + iindex_terrain; // (x1,y1)
550  ind12 = std::min(jindex_terrain+1,ny-1) * nx + iindex_terrain; // (x1,y2)
551 
552  ind21 = jindex_terrain * nx + std::min(iindex_terrain+1,nx-1); // (x2,y1)
553  ind22 = std::min(jindex_terrain+1,ny-1) * nx + std::min(iindex_terrain+1,nx-1); // (x2,y2)
554 #endif
555 
556  if (iindex_terrain == nx-1 && jindex_terrain == ny-1)
557  {
558  z_arr(i,j,klo) = d_zt[ind11];
559  }
560  else if (iindex_terrain != nx-1 && jindex_terrain == ny-1)
561  {
562  amrex::Real w_11 = (x2-x); // (x1,y1)
563  amrex::Real w_21 = (x-x1); // (x2,y1)
564  amrex::Real denom = (x2-x1);
565  z_arr(i,j,klo) = (w_11*d_zt[ind11] + w_21*d_zt[ind21])/denom;
566  }
567  else if (iindex_terrain == nx-1 && jindex_terrain != ny-1)
568  {
569  amrex::Real w_11 = (y2-y); // (x1,y1)
570  amrex::Real w_12 = (y-y1); // (x1,y2)
571  amrex::Real denom = (y2-y1);
572  z_arr(i,j,klo) = (w_11*d_zt[ind11] + w_12*d_zt[ind12])/denom;
573  }
574  else
575  {
576  amrex::Real w_11 = (x2-x)*(y2-y); // (x1,y1)
577  amrex::Real w_21 = (x-x1)*(y2-y); // (x2,y1)
578  amrex::Real w_12 = (x2-x)*(y-y1); // (x1,y2)
579  amrex::Real w_22 = (x-x1)*(y-y1); // (x2,y2)
580  amrex::Real denom = (x2-x1)*(y2-y1);
581  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;
582  }
583  } // usgs?
584  });
585  }
const Box zbx
Definition: ERF_SetupDiff.H:9

Referenced by init_buildings_surface(), and init_terrain_surface().

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
253  {
254  if (u_geos.empty()) return;
255 
256  amrex::Warning("Geostrophic wind profile not defined for "+name()+" problem");
257 
258  const int khi = geom.Domain().bigEnd()[2];
259  // const amrex::Real* prob_lo = geom.ProbLo();
260  // const auto dx = geom.CellSize();
261  for (int k = 0; k <= khi; k++)
262  {
263  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
264  // set RHS term of RhoTheta equation based on time, z_cc here
265  u_geos[k] = 0.0;
266  v_geos[k] = 0.0;
267  }
268 
269  // Copy from host version to device version
270  amrex::Gpu::copy(amrex::Gpu::hostToDevice, u_geos.begin(), u_geos.end(), d_u_geos.begin());
271  amrex::Gpu::copy(amrex::Gpu::hostToDevice, v_geos.begin(), v_geos.end(), d_v_geos.begin());
272  }
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
180  {
181  if (qsrc->empty()) return;
182 
183  amrex::Warning("Moisture forcing not defined for "+name()+" problem");
184  for ( amrex::MFIter mfi(*qsrc, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
185  {
186  const auto &box = mfi.tilebox();
187  const amrex::Array4<amrex::Real>& qsrc_arr = qsrc->array(mfi);
188  // src is a spatial function if erf.spatial_moisture_forcing = true
189  // otherwise, qsrc_arr is defined only in over Z (box x and y dimensions are 1)
190  ParallelFor(box, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
191  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
192  // set RHS term of RhoQ1 equation based on time, z_cc here
193  qsrc_arr(i, j, k) = 0.0;
194  });
195  }
196  }
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
148  {
149  if (src->empty()) return;
150 
151  amrex::Warning("Temperature forcing not defined for "+name()+" problem");
152  for ( amrex::MFIter mfi(*src, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
153  {
154  const auto &box = mfi.tilebox();
155  const amrex::Array4<amrex::Real>& src_arr = src->array(mfi);
156  // src is a spatial function if erf.spatial_rhotheta_forcing = true
157  // otherwise, qsrc_arr is defined only in over Z (box x and y dimensions are 1)
158  ParallelFor(box, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
159  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
160  // set RHS term of RhoTheta equation based on time, z_cc here
161  src_arr(i, j, k) = 0.0;
162  });
163  }
164  }
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
217  {
218  if (wbar.empty()) return;
219 
220  amrex::Warning("Moisture forcing not defined for "+name()+" problem");
221 
222  const int khi = geom.Domain().bigEnd()[2];
223  // const amrex::Real* prob_lo = geom.ProbLo();
224  // const auto dx = geom.CellSize();
225  for (int k = 0; k <= khi; k++)
226  {
227  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
228  // set vertical velocity profile based on time, z_cc here
229  wbar[k] = 0.0;
230  }
231 
232  // Copy from host version to device version
233  amrex::Gpu::copy(amrex::Gpu::hostToDevice, wbar.begin(), wbar.end(), d_wbar.begin());
234  }
@ 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(), and init_uniform().


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