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
484  {
485  amrex::ParmParse pp("erf");
486 
487  // Count number of lo and hi points define the plane
488  int n_plane_lo = pp.countval("sample_plane_lo") / AMREX_SPACEDIM;
489  int n_plane_hi = pp.countval("sample_plane_hi") / AMREX_SPACEDIM;
490  int n_plane_dir = pp.countval("sample_plane_dir");
491  AMREX_ALWAYS_ASSERT( (n_plane_lo==n_plane_hi ) &&
492  (n_plane_lo==n_plane_dir) );
493 
494  // Parse the data
495  if (n_plane_lo > 0) {
496  // Parse lo
497  amrex::Vector<amrex::Real> r_lo; r_lo.resize(n_plane_lo*AMREX_SPACEDIM);
498  amrex::Vector<amrex::Vector<amrex::Real>> rv_lo;
499  pp.queryarr("sample_plane_lo",r_lo,0,n_plane_lo*AMREX_SPACEDIM);
500  for (int i(0); i < n_plane_lo; i++) {
501  amrex::Vector<amrex::Real> rv = {r_lo[AMREX_SPACEDIM*i+0],
502  r_lo[AMREX_SPACEDIM*i+1],
503  r_lo[AMREX_SPACEDIM*i+2]};
504  rv_lo.push_back(rv);
505  }
506 
507  // Parse hi
508  amrex::Vector<amrex::Real> r_hi; r_hi.resize(n_plane_hi*AMREX_SPACEDIM);
509  amrex::Vector<amrex::Vector<amrex::Real>> rv_hi;
510  pp.queryarr("sample_plane_hi",r_hi,0,n_plane_hi*AMREX_SPACEDIM);
511  for (int i(0); i < n_plane_hi; i++) {
512  amrex::Vector<amrex::Real> rv = {r_hi[AMREX_SPACEDIM*i+0],
513  r_hi[AMREX_SPACEDIM*i+1],
514  r_hi[AMREX_SPACEDIM*i+2]};
515  rv_hi.push_back(rv);
516  }
517 
518  // Construct vector of bounding real boxes
519  m_bnd_rbx.resize(n_plane_lo);
520  for (int i(0); i < n_plane_hi; i++){
521  amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
522  m_bnd_rbx[i] = rbx;
523  }
524 
525  // Parse directionality
526  m_dir.resize(n_plane_dir);
527  pp.queryarr("sample_plane_dir",m_dir,0,n_plane_dir);
528 
529  // Parse names
530  std::string name_base = "plt_plane_";
531  m_name.resize(n_plane_lo);
532  int n_names = pp.countval("sample_plane_name");
533  if (n_names > 0) {
534  AMREX_ALWAYS_ASSERT( n_names==n_plane_lo );
535  pp.queryarr("sample_plane_name",m_name,0,n_names);
536  } else {
537  for (int iplane(0); iplane<n_plane_lo; ++iplane) {
538  m_name[iplane] = amrex::Concatenate(name_base, iplane , 5);
539  }
540  }
541 
542  // Allocate space for level indicator
543  m_lev.resize(n_plane_dir,0);
544 
545  // Allocate space for MF pointers
546  m_ps_mf.resize(n_plane_lo);
547  }
548  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:233
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:676
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:677
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:680
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_ps_mf
Definition: ERF_SampleData.H:679
amrex::Vector< amrex::RealBox > m_bnd_rbx
Definition: ERF_SampleData.H:678
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
570  {
571  int nlev = static_cast<int>(vars_new.size());
572  int nplane = static_cast<int>(m_bnd_rbx.size());
573  int ncomp = 2;
574  bool interpolate = true;
575 
576  // Loop over each plane
577  for (int iplane(0); iplane<nplane; ++iplane) {
578  int dir = m_dir[iplane];
579  amrex::RealBox bnd_rbx = m_bnd_rbx[iplane];
580  amrex::Real point = bnd_rbx.lo(dir);
581 
582  // Search each level to get the finest data possible
583  for (int ilev(nlev-1); ilev>=0; --ilev) {
584 
585  // Construct CC velocities
586  amrex::MultiFab mf_cc_vel;
587  auto ba = vars_new[ilev][Vars::cons].boxArray();
588  auto dm = vars_new[ilev][Vars::cons].DistributionMap();
589  mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
590  average_face_to_cellcenter(mf_cc_vel,0,
591  amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][Vars::xvel],
592  &vars_new[ilev][Vars::yvel],
593  &vars_new[ilev][Vars::zvel]});
594 
595  // Construct vector of MFs holding T and WSP
596  amrex::MultiFab mf_cc_data;
597  mf_cc_data.define(ba, dm, ncomp, 1);
598 #ifdef _OPENMP
599 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
600 #endif
601  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
602  const amrex::Box& tbx = mfi.tilebox();
603  auto const& dfab = mf_cc_data.array(mfi);
604  auto const& tfab = vars_new[ilev][Vars::cons].array(mfi);
605  auto const& wfab = mf_cc_vel.array(mfi);
606  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
607  {
608  dfab(i,j,k,0) = tfab(i,j,k,1)/tfab(i,j,k,0);
609  dfab(i,j,k,1) = std::sqrt(wfab(i,j,k,0)*wfab(i,j,k,0)
610  + wfab(i,j,k,1)*wfab(i,j,k,1)
611  + wfab(i,j,k,2)*wfab(i,j,k,2)) ;
612  });
613 
614  }
615 
616  m_lev[iplane] = ilev;
617  m_ps_mf[iplane] = get_slice_data(dir, point, mf_cc_data, geom[ilev],
618  0, ncomp, interpolate, bnd_rbx);
619 
620  // We can stop if we got the entire plane
621  auto min_bnd_bx = m_ps_mf[iplane]->boxArray().minimalBox();
622  amrex::Box bnd_bx = getIndexBox(bnd_rbx, geom[ilev]);
623  if (bnd_bx == min_bnd_bx) { break; }
624 
625  } // ilev
626  }// iplane
627  }
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
amrex::Box getIndexBox(const amrex::RealBox &real_box, const amrex::Geometry &geom)
Definition: ERF_SampleData.H:552
Here is the call graph for this function:

◆ getIndexBox()

amrex::Box PlaneSampler::getIndexBox ( const amrex::RealBox &  real_box,
const amrex::Geometry &  geom 
)
inline
553  {
554  amrex::IntVect slice_lo, slice_hi;
555 
556  AMREX_D_TERM(slice_lo[0]=static_cast<int>(std::floor((real_box.lo(0) - geom.ProbLo(0))/geom.CellSize(0)));,
557  slice_lo[1]=static_cast<int>(std::floor((real_box.lo(1) - geom.ProbLo(1))/geom.CellSize(1)));,
558  slice_lo[2]=static_cast<int>(std::floor((real_box.lo(2) - geom.ProbLo(2))/geom.CellSize(2))););
559 
560  AMREX_D_TERM(slice_hi[0]=static_cast<int>(std::floor((real_box.hi(0) - geom.ProbLo(0))/geom.CellSize(0)));,
561  slice_hi[1]=static_cast<int>(std::floor((real_box.hi(1) - geom.ProbLo(1))/geom.CellSize(1)));,
562  slice_hi[2]=static_cast<int>(std::floor((real_box.hi(2) - geom.ProbLo(2))/geom.CellSize(2))););
563 
564  return amrex::Box(slice_lo, slice_hi) & geom.Domain();
565  }

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
634  {
635  amrex::Vector<std::string> varnames = {"T", "Wsp"};
636 
637  int nplane = m_ps_mf.size();
638  for (int iplane(0); iplane<nplane; ++iplane) {
639  // Data members that can be used as-is
640  int dir = m_dir[iplane];
641  int lev = m_lev[iplane];
642  amrex::Real m_time = time[lev];
643  amrex::Vector<int> m_level_steps = {level_steps[lev]};
644  amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
645 
646  // Create modified geometry object corresponding to the plane
647  amrex::RealBox m_rb = m_bnd_rbx[iplane];
648  amrex::Box m_dom = getIndexBox(m_rb, geom[lev]);
649  amrex::Real point = m_rb.hi(dir);
650  amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
651  for (int d(0); d<AMREX_SPACEDIM; ++d) {
652  if (d==dir) {
653  m_rb.setLo(d,point-0.5*geom[lev].CellSize(d));
654  m_rb.setHi(d,point+0.5*geom[lev].CellSize(d));
655  }
656  is_per[d] = geom[lev].isPeriodic(d);
657  }
658  amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
659  m_geom[0].define(m_dom, &m_rb, geom[lev].Coord(), is_per.data());
660 
661  // Create plotfile name
662  std::string name_plane = m_name[iplane];
663  name_plane += "_step_";
664  std::string plotfilename = amrex::Concatenate(name_plane, m_level_steps[0], 5);
665 
666  // Get the data
667  amrex::Vector<const amrex::MultiFab*> mf = {m_ps_mf[iplane].get()};
668 
669  // Write each plane
670  WriteMultiLevelPlotfile(plotfilename, 1, mf,
671  varnames, m_geom, m_time,
672  m_level_steps, m_ref_ratio);
673  } // iplane
674  }
Coord
Definition: ERF_DataStruct.H:85
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: