ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
PlaneSampler Struct Reference

#include <ERF_SampleData.H>

Collaboration diagram for PlaneSampler:

Public Member Functions

 PlaneSampler ()
 
amrex::Box getIndexBox (const amrex::RealBox &real_box, const amrex::Geometry &geom)
 
void get_sample_data (amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
 
void write_sample_data (amrex::Vector< amrex::Real > &time, amrex::Vector< int > &level_steps, amrex::Vector< amrex::IntVect > &ref_ratio, amrex::Vector< amrex::Geometry > &geom)
 

Public Attributes

amrex::Vector< int > m_dir
 
amrex::Vector< int > m_lev
 
amrex::Vector< amrex::RealBox > m_bnd_rbx
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_ps_mf
 
amrex::Vector< std::string > m_name
 

Constructor & Destructor Documentation

◆ PlaneSampler()

PlaneSampler::PlaneSampler ( )
inline
568  {
569  amrex::ParmParse pp("erf");
570 
571  // Count number of lo and hi points define the plane
572  int n_plane_lo = pp.countval("sample_plane_lo") / AMREX_SPACEDIM;
573  int n_plane_hi = pp.countval("sample_plane_hi") / AMREX_SPACEDIM;
574  int n_plane_dir = pp.countval("sample_plane_dir");
575  AMREX_ALWAYS_ASSERT( (n_plane_lo==n_plane_hi ) &&
576  (n_plane_lo==n_plane_dir) );
577 
578  // Parse the data
579  if (n_plane_lo > 0) {
580  // Parse lo
581  amrex::Vector<amrex::Real> r_lo; r_lo.resize(n_plane_lo*AMREX_SPACEDIM);
582  amrex::Vector<amrex::Vector<amrex::Real>> rv_lo;
583  pp.queryarr("sample_plane_lo",r_lo,0,n_plane_lo*AMREX_SPACEDIM);
584  for (int i(0); i < n_plane_lo; i++) {
585  amrex::Vector<amrex::Real> rv = {r_lo[AMREX_SPACEDIM*i+0],
586  r_lo[AMREX_SPACEDIM*i+1],
587  r_lo[AMREX_SPACEDIM*i+2]};
588  rv_lo.push_back(rv);
589  }
590 
591  // Parse hi
592  amrex::Vector<amrex::Real> r_hi; r_hi.resize(n_plane_hi*AMREX_SPACEDIM);
593  amrex::Vector<amrex::Vector<amrex::Real>> rv_hi;
594  pp.queryarr("sample_plane_hi",r_hi,0,n_plane_hi*AMREX_SPACEDIM);
595  for (int i(0); i < n_plane_hi; i++) {
596  amrex::Vector<amrex::Real> rv = {r_hi[AMREX_SPACEDIM*i+0],
597  r_hi[AMREX_SPACEDIM*i+1],
598  r_hi[AMREX_SPACEDIM*i+2]};
599  rv_hi.push_back(rv);
600  }
601 
602  // Construct vector of bounding real boxes
603  m_bnd_rbx.resize(n_plane_lo);
604  for (int i(0); i < n_plane_hi; i++){
605  amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
606  m_bnd_rbx[i] = rbx;
607  }
608 
609  // Parse directionality
610  m_dir.resize(n_plane_dir);
611  pp.queryarr("sample_plane_dir",m_dir,0,n_plane_dir);
612 
613  // Parse names
614  std::string name_base = "plt_plane_";
615  m_name.resize(n_plane_lo);
616  int n_names = pp.countval("sample_plane_name");
617  if (n_names > 0) {
618  AMREX_ALWAYS_ASSERT( n_names==n_plane_lo );
619  pp.queryarr("sample_plane_name",m_name,0,n_names);
620  } else {
621  for (int iplane(0); iplane<n_plane_lo; ++iplane) {
622  m_name[iplane] = amrex::Concatenate(name_base, iplane , 5);
623  }
624  }
625 
626  // Allocate space for level indicator
627  m_lev.resize(n_plane_dir,0);
628 
629  // Allocate space for MF pointers
630  m_ps_mf.resize(n_plane_lo);
631  }
632  }
ParmParse pp("prob")
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:760
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:761
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:764
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_ps_mf
Definition: ERF_SampleData.H:763
amrex::Vector< amrex::RealBox > m_bnd_rbx
Definition: ERF_SampleData.H:762
Here is the call graph for this function:

Member Function Documentation

◆ get_sample_data()

void PlaneSampler::get_sample_data ( amrex::Vector< amrex::Geometry > &  geom,
amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars_new 
)
inline
654  {
655  int nlev = static_cast<int>(vars_new.size());
656  int nplane = static_cast<int>(m_bnd_rbx.size());
657  int ncomp = 2;
658  bool interpolate = true;
659 
660  // Loop over each plane
661  for (int iplane(0); iplane<nplane; ++iplane) {
662  int dir = m_dir[iplane];
663  amrex::RealBox bnd_rbx = m_bnd_rbx[iplane];
664  amrex::Real point = bnd_rbx.lo(dir);
665 
666  // Search each level to get the finest data possible
667  for (int ilev(nlev-1); ilev>=0; --ilev) {
668 
669  // Construct CC velocities
670  amrex::MultiFab mf_cc_vel;
671  auto ba = vars_new[ilev][Vars::cons].boxArray();
672  auto dm = vars_new[ilev][Vars::cons].DistributionMap();
673  mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
674  average_face_to_cellcenter(mf_cc_vel,0,
675  amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][Vars::xvel],
676  &vars_new[ilev][Vars::yvel],
677  &vars_new[ilev][Vars::zvel]});
678 
679  // Construct vector of MFs holding T and WSP
680  amrex::MultiFab mf_cc_data;
681  mf_cc_data.define(ba, dm, ncomp, 1);
682 #ifdef _OPENMP
683 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
684 #endif
685  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
686  const amrex::Box& tbx = mfi.tilebox();
687  auto const& dfab = mf_cc_data.array(mfi);
688  auto const& tfab = vars_new[ilev][Vars::cons].array(mfi);
689  auto const& wfab = mf_cc_vel.array(mfi);
690  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
691  {
692  dfab(i,j,k,0) = tfab(i,j,k,1)/tfab(i,j,k,0);
693  dfab(i,j,k,1) = std::sqrt(wfab(i,j,k,0)*wfab(i,j,k,0)
694  + wfab(i,j,k,1)*wfab(i,j,k,1)
695  + wfab(i,j,k,2)*wfab(i,j,k,2)) ;
696  });
697 
698  }
699 
700  m_lev[iplane] = ilev;
701  m_ps_mf[iplane] = get_slice_data(dir, point, mf_cc_data, geom[ilev],
702  0, ncomp, interpolate, bnd_rbx);
703 
704  // We can stop if we got the entire plane
705  auto min_bnd_bx = m_ps_mf[iplane]->boxArray().minimalBox();
706  amrex::Box bnd_bx = getIndexBox(bnd_rbx, geom[ilev]);
707  if (bnd_bx == min_bnd_bx) { break; }
708 
709  } // ilev
710  }// iplane
711  }
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+0.5) *dx[0];const Real y=(j+0.5) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - 1.0)/(2.0 *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(1.0+deltaT, inv_gm1);const Real T=(1.0+deltaT) *T_inf;const Real p=std::pow(rho_norm, Gamma)/Gamma *rho_0 *a_inf *a_inf;const Real rho_theta=rho_0 *rho_norm *(T *std::pow(p_0/p, rdOcp));state_pert(i, j, k, RhoTheta_comp)=rho_theta - getRhoThetagivenP(p_hse(i, j, k));const Real r2d_xy=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc));state_pert(i, j, k, RhoScalar_comp)=0.25 *(1.0+std::cos(PI *std::min(r2d_xy, R)/R));})
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ xvel
Definition: ERF_IndexDefines.H:159
@ cons
Definition: ERF_IndexDefines.H:158
@ zvel
Definition: ERF_IndexDefines.H:161
@ yvel
Definition: ERF_IndexDefines.H:160
amrex::Box getIndexBox(const amrex::RealBox &real_box, const amrex::Geometry &geom)
Definition: ERF_SampleData.H:636
Here is the call graph for this function:

◆ getIndexBox()

amrex::Box PlaneSampler::getIndexBox ( const amrex::RealBox &  real_box,
const amrex::Geometry &  geom 
)
inline
637  {
638  amrex::IntVect slice_lo, slice_hi;
639 
640  AMREX_D_TERM(slice_lo[0]=static_cast<int>(std::floor((real_box.lo(0) - geom.ProbLo(0))/geom.CellSize(0)));,
641  slice_lo[1]=static_cast<int>(std::floor((real_box.lo(1) - geom.ProbLo(1))/geom.CellSize(1)));,
642  slice_lo[2]=static_cast<int>(std::floor((real_box.lo(2) - geom.ProbLo(2))/geom.CellSize(2))););
643 
644  AMREX_D_TERM(slice_hi[0]=static_cast<int>(std::floor((real_box.hi(0) - geom.ProbLo(0))/geom.CellSize(0)));,
645  slice_hi[1]=static_cast<int>(std::floor((real_box.hi(1) - geom.ProbLo(1))/geom.CellSize(1)));,
646  slice_hi[2]=static_cast<int>(std::floor((real_box.hi(2) - geom.ProbLo(2))/geom.CellSize(2))););
647 
648  return amrex::Box(slice_lo, slice_hi) & geom.Domain();
649  }

Referenced by get_sample_data(), and write_sample_data().

Here is the caller graph for this function:

◆ write_sample_data()

void PlaneSampler::write_sample_data ( amrex::Vector< amrex::Real > &  time,
amrex::Vector< int > &  level_steps,
amrex::Vector< amrex::IntVect > &  ref_ratio,
amrex::Vector< amrex::Geometry > &  geom 
)
inline
718  {
719  amrex::Vector<std::string> varnames = {"T", "Wsp"};
720 
721  int nplane = m_ps_mf.size();
722  for (int iplane(0); iplane<nplane; ++iplane) {
723  // Data members that can be used as-is
724  int dir = m_dir[iplane];
725  int lev = m_lev[iplane];
726  amrex::Real m_time = time[lev];
727  amrex::Vector<int> m_level_steps = {level_steps[lev]};
728  amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
729 
730  // Create modified geometry object corresponding to the plane
731  amrex::RealBox m_rb = m_bnd_rbx[iplane];
732  amrex::Box m_dom = getIndexBox(m_rb, geom[lev]);
733  amrex::Real point = m_rb.hi(dir);
734  amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
735  for (int d(0); d<AMREX_SPACEDIM; ++d) {
736  if (d==dir) {
737  m_rb.setLo(d,point-myhalf*geom[lev].CellSize(d));
738  m_rb.setHi(d,point+myhalf*geom[lev].CellSize(d));
739  }
740  is_per[d] = geom[lev].isPeriodic(d);
741  }
742  amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
743  m_geom[0].define(m_dom, &m_rb, geom[lev].Coord(), is_per.data());
744 
745  // Create plotfile name
746  std::string name_plane = m_name[iplane];
747  name_plane += "_step_";
748  std::string plotfilename = amrex::Concatenate(name_plane, m_level_steps[0], 5);
749 
750  // Get the data
751  amrex::Vector<const amrex::MultiFab*> mf = {m_ps_mf[iplane].get()};
752 
753  // Write each plane
754  WriteMultiLevelPlotfile(plotfilename, 1, mf,
755  varnames, m_geom, m_time,
756  m_level_steps, m_ref_ratio);
757  } // iplane
758  }
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
Coord
Definition: ERF_DataStruct.H:92
Here is the call graph for this function:

Member Data Documentation

◆ m_bnd_rbx

amrex::Vector<amrex::RealBox> PlaneSampler::m_bnd_rbx

◆ m_dir

amrex::Vector<int> PlaneSampler::m_dir

◆ m_lev

amrex::Vector<int> PlaneSampler::m_lev

◆ m_name

amrex::Vector<std::string> PlaneSampler::m_name

Referenced by PlaneSampler(), and write_sample_data().

◆ m_ps_mf

amrex::Vector<std::unique_ptr<amrex::MultiFab> > PlaneSampler::m_ps_mf

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