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

#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::Vector< amrex::Real > &src, amrex::Gpu::DeviceVector< amrex::Real > &d_src, const amrex::Geometry &geom, std::unique_ptr< amrex::MultiFab > &)
 
virtual void update_rhoqt_sources (const amrex::Real &, amrex::Vector< amrex::Real > &qsrc, amrex::Gpu::DeviceVector< amrex::Real > &d_qsrc, const amrex::Geometry &geom, 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::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 &, amrex::FArrayBox &terrain_fab, const amrex::Real &)
 
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)
 
virtual std::string name ()=0
 

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
65  {
66  amrex::Print() << "Hydrostatically balanced density was NOT set"
67  << " -- an appropriate init_type should probably have been specified"
68  << " (e.g., input_sounding, WRFInput, or Metgrid)"
69  << std::endl;
70  amrex::Error("Should never call erf_init_dens_hse for "+name()+" problem");
71  }
virtual std::string name()=0
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
77  {
78 
79  }

◆ 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
628  {
629  amrex::Error("Should never call erf_init_rayleigh for "+name()+" problem");
630  }
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

659  {
660  base_parms.rho_0 = rho_0;
661  base_parms.T_0 = T_0;
662  }
ProbParmDefaults base_parms
Definition: ERF_ProbCommon.H:653
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
283  {
284  // Check if a valid text file exists for the buildings
285  std::string fname;
286  amrex::ParmParse pp("erf");
287  auto valid_fname = pp.query("buildings_file_name",fname);
288 
289  if (valid_fname) {
290  read_custom_terrain(fname,false,geom,buildings_fab,time);
291  } else {
292  init_custom_terrain(geom, buildings_fab, time);
293  }
294  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:233
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:331
virtual void init_custom_terrain(const amrex::Geometry &, amrex::FArrayBox &terrain_fab, const amrex::Real &)
Definition: ERF_ProbCommon.H:595
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
124  {
125  amrex::Print() << "No perturbation to background fields supplied for "
126  << name() << " problem" << std::endl;
127  }
Here is the call graph for this function:

◆ init_custom_terrain()

virtual void ProblemBase::init_custom_terrain ( const amrex::Geometry &  ,
amrex::FArrayBox &  terrain_fab,
const amrex::Real  
)
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
598  {
599  amrex::Print() << "Initializing flat terrain" << std::endl;
600 
601  if (SolverChoice::mesh_type == MeshType::VariableDz) {
602  SolverChoice::set_mesh_type(MeshType::StretchedDz);
603  amrex::Print() << "Resetting mesh type to StretchedDz" << std::endl;
604  }
605 
606  terrain_fab.template setVal<amrex::RunOn::Device>(0.0);
607  }
static MeshType mesh_type
Definition: ERF_DataStruct.H:983
static void set_mesh_type(MeshType new_mesh_type)
Definition: ERF_DataStruct.H:986

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:

◆ 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
308  {
309  // Check if a valid text file exists for the terrain
310  std::string fname, fname_usgs;
311  amrex::ParmParse pp("erf");
312  auto valid_fname = pp.query("terrain_file_name",fname);
313  auto valid_fname_USGS = pp.query("terrain_file_name_USGS",fname_usgs);
314 
315  bool is_usgs;
316 
317  if (valid_fname) {
318  is_usgs = false;
319  read_custom_terrain(fname,is_usgs,geom,terrain_fab,time);
320 
321  } else if (valid_fname_USGS) {
322  is_usgs = true;
323  read_custom_terrain(fname_usgs,is_usgs,geom,terrain_fab,time);
324 
325  } else {
326  init_custom_terrain (geom, terrain_fab, time);
327  }
328  }
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

637  {
638  amrex::Real rho_0 = base_parms.rho_0;
639  amrex::Real T_0 = base_parms.T_0;
640  amrex::Print() << "Initializing uniform fields"
641  << " rho=" << rho_0 << " theta=" << T_0
642  << std::endl;
643 
644  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
645  {
646  state(i, j, k, Rho_comp) = rho_0;
647  state(i, j, k, RhoTheta_comp) = rho_0 * T_0;
648  });
649  }
#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()

virtual std::string ProblemBase::name ( )
protectedpure virtual

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

◆ update_rhoqt_sources()

virtual void ProblemBase::update_rhoqt_sources ( const amrex::Real ,
amrex::Vector< amrex::Real > &  qsrc,
amrex::Gpu::DeviceVector< amrex::Real > &  d_qsrc,
const amrex::Geometry &  geom,
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
178  {
179  if (qsrc.empty()) return;
180 
181  amrex::Warning("Moisture forcing not defined for "+name()+" problem");
182 
183  const int khi = geom.Domain().bigEnd()[2];
184  // const amrex::Real* prob_lo = geom.ProbLo();
185  // const auto dx = geom.CellSize();
186  for (int k = 0; k <= khi; k++)
187  {
188  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
189  // set RHS term of RhoQ1 equation based on time, z_cc here
190  qsrc[k] = 0.0;
191  }
192 
193  // Copy from host version to device version
194  amrex::Gpu::copy(amrex::Gpu::hostToDevice, qsrc.begin(), qsrc.end(), d_qsrc.begin());
195  }
Here is the call graph for this function:

◆ update_rhotheta_sources()

virtual void ProblemBase::update_rhotheta_sources ( const amrex::Real ,
amrex::Vector< amrex::Real > &  src,
amrex::Gpu::DeviceVector< amrex::Real > &  d_src,
const amrex::Geometry &  geom,
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
144  {
145  if (src.empty()) return;
146 
147  amrex::Warning("Temperature forcing not defined for "+name()+" problem");
148 
149  const int khi = geom.Domain().bigEnd()[2];
150  // const amrex::Real* prob_lo = geom.ProbLo();
151  // const auto dx = geom.CellSize();
152  for (int k = 0; k <= khi; k++)
153  {
154  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
155  // set RHS term of RhoTheta equation based on time, z_cc here
156  src[k] = 0.0;
157  }
158 
159  // Copy from host version to device version
160  amrex::Gpu::copy(amrex::Gpu::hostToDevice, src.begin(), src.end(), d_src.begin());
161  }
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::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
215  {
216  if (wbar.empty()) return;
217 
218  amrex::Warning("Moisture forcing not defined for "+name()+" problem");
219 
220  const int khi = geom.Domain().bigEnd()[2];
221  // const amrex::Real* prob_lo = geom.ProbLo();
222  // const auto dx = geom.CellSize();
223  for (int k = 0; k <= khi; k++)
224  {
225  // const amrex::Real z_cc = prob_lo[2] + (k+0.5)* dx[2];
226  // set vertical velocity profile based on time, z_cc here
227  wbar[k] = 0.0;
228  }
229 
230  // Copy from host version to device version
231  amrex::Gpu::copy(amrex::Gpu::hostToDevice, wbar.begin(), wbar.end(), d_wbar.begin());
232  }
@ wbar
Definition: ERF_DataStruct.H:96
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: