ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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_plane_mfs (amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
 
void write_plane_mfs (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
398  {
399  amrex::ParmParse pp("erf");
400 
401  // Count number of lo and hi points define the plane
402  int n_plane_lo = pp.countval("sample_plane_lo") / AMREX_SPACEDIM;
403  int n_plane_hi = pp.countval("sample_plane_hi") / AMREX_SPACEDIM;
404  int n_plane_dir = pp.countval("sample_plane_dir");
405  AMREX_ALWAYS_ASSERT( (n_plane_lo==n_plane_hi ) &&
406  (n_plane_lo==n_plane_dir) );
407 
408  // Parse the data
409  if (n_plane_lo > 0) {
410  // Parse lo
411  amrex::Vector<amrex::Real> r_lo; r_lo.resize(n_plane_lo*AMREX_SPACEDIM);
412  amrex::Vector<amrex::Vector<amrex::Real>> rv_lo;
413  pp.queryarr("sample_plane_lo",r_lo,0,n_plane_lo*AMREX_SPACEDIM);
414  for (int i(0); i < n_plane_lo; i++) {
415  amrex::Vector<amrex::Real> rv = {r_lo[AMREX_SPACEDIM*i+0],
416  r_lo[AMREX_SPACEDIM*i+1],
417  r_lo[AMREX_SPACEDIM*i+2]};
418  rv_lo.push_back(rv);
419  }
420 
421  // Parse hi
422  amrex::Vector<amrex::Real> r_hi; r_hi.resize(n_plane_hi*AMREX_SPACEDIM);
423  amrex::Vector<amrex::Vector<amrex::Real>> rv_hi;
424  pp.queryarr("sample_plane_hi",r_hi,0,n_plane_hi*AMREX_SPACEDIM);
425  for (int i(0); i < n_plane_hi; i++) {
426  amrex::Vector<amrex::Real> rv = {r_hi[AMREX_SPACEDIM*i+0],
427  r_hi[AMREX_SPACEDIM*i+1],
428  r_hi[AMREX_SPACEDIM*i+2]};
429  rv_hi.push_back(rv);
430  }
431 
432  // Construct vector of bounding real boxes
433  m_bnd_rbx.resize(n_plane_lo);
434  for (int i(0); i < n_plane_hi; i++){
435  amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
436  m_bnd_rbx[i] = rbx;
437  }
438 
439  // Parse directionality
440  m_dir.resize(n_plane_dir);
441  pp.queryarr("sample_plane_dir",m_dir,0,n_plane_dir);
442 
443  // Parse names
444  std::string name_base = "plt_plane_";
445  m_name.resize(n_plane_lo);
446  int n_names = pp.countval("sample_plane_name");
447  if (n_names > 0) {
448  AMREX_ALWAYS_ASSERT( n_names==n_plane_lo );
449  pp.queryarr("sample_plane_name",m_name,0,n_names);
450  } else {
451  for (int iplane(0); iplane<n_plane_lo; ++iplane) {
452  m_name[iplane] = amrex::Concatenate(name_base, iplane , 5);
453  }
454  }
455 
456  // Allocate space for level indicator
457  m_lev.resize(n_plane_dir,0);
458 
459  // Allocate space for MF pointers
460  m_ps_mf.resize(n_plane_lo);
461  }
462  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:219
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:590
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:591
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:594
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_ps_mf
Definition: ERF_SampleData.H:593
amrex::Vector< amrex::RealBox > m_bnd_rbx
Definition: ERF_SampleData.H:592
Here is the call graph for this function:

Member Function Documentation

◆ get_plane_mfs()

void PlaneSampler::get_plane_mfs ( amrex::Vector< amrex::Geometry > &  geom,
amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars_new 
)
inline
484  {
485  int nlev = vars_new.size();
486  int nplane = m_bnd_rbx.size();
487  int ncomp = 2;
488  bool interpolate = true;
489 
490  // Loop over each plane
491  for (int iplane(0); iplane<nplane; ++iplane) {
492  int dir = m_dir[iplane];
493  amrex::RealBox bnd_rbx = m_bnd_rbx[iplane];
494  amrex::Real point = bnd_rbx.lo(dir);
495 
496  // Search each level to get the finest data possible
497  for (int ilev(nlev-1); ilev>=0; --ilev) {
498 
499  // Construct CC velocities
500  amrex::MultiFab mf_cc_vel;
501  auto ba = vars_new[ilev][Vars::cons].boxArray();
502  auto dm = vars_new[ilev][Vars::cons].DistributionMap();
503  mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
504  average_face_to_cellcenter(mf_cc_vel,0,
505  amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][Vars::xvel],
506  &vars_new[ilev][Vars::yvel],
507  &vars_new[ilev][Vars::zvel]});
508 
509  // Construct vector of MFs holding T and WSP
510  amrex::MultiFab mf_cc_data;
511  mf_cc_data.define(ba, dm, ncomp, 1);
512 #ifdef _OPENMP
513 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
514 #endif
515  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
516  const amrex::Box& tbx = mfi.tilebox();
517  auto const& dfab = mf_cc_data.array(mfi);
518  auto const& tfab = vars_new[ilev][Vars::cons].array(mfi);
519  auto const& wfab = mf_cc_vel.array(mfi);
520  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
521  {
522  dfab(i,j,k,0) = tfab(i,j,k,1)/tfab(i,j,k,0);
523  dfab(i,j,k,1) = std::sqrt(wfab(i,j,k,0)*wfab(i,j,k,0)
524  + wfab(i,j,k,1)*wfab(i,j,k,1)
525  + wfab(i,j,k,2)*wfab(i,j,k,2)) ;
526  });
527 
528  }
529 
530  m_lev[iplane] = ilev;
531  m_ps_mf[iplane] = get_slice_data(dir, point, mf_cc_data, geom[ilev],
532  0, ncomp, interpolate, bnd_rbx);
533 
534  // We can stop if we got the entire plane
535  auto min_bnd_bx = m_ps_mf[iplane]->boxArray().minimalBox();
536  amrex::Box bnd_bx = getIndexBox(bnd_rbx, geom[ilev]);
537  if (bnd_bx == min_bnd_bx) { continue; }
538 
539  } // ilev
540  }// iplane
541  }
@ 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:466
Here is the call graph for this function:

◆ getIndexBox()

amrex::Box PlaneSampler::getIndexBox ( const amrex::RealBox &  real_box,
const amrex::Geometry &  geom 
)
inline
467  {
468  amrex::IntVect slice_lo, slice_hi;
469 
470  AMREX_D_TERM(slice_lo[0]=static_cast<int>(std::floor((real_box.lo(0) - geom.ProbLo(0))/geom.CellSize(0)));,
471  slice_lo[1]=static_cast<int>(std::floor((real_box.lo(1) - geom.ProbLo(1))/geom.CellSize(1)));,
472  slice_lo[2]=static_cast<int>(std::floor((real_box.lo(2) - geom.ProbLo(2))/geom.CellSize(2))););
473 
474  AMREX_D_TERM(slice_hi[0]=static_cast<int>(std::floor((real_box.hi(0) - geom.ProbLo(0))/geom.CellSize(0)));,
475  slice_hi[1]=static_cast<int>(std::floor((real_box.hi(1) - geom.ProbLo(1))/geom.CellSize(1)));,
476  slice_hi[2]=static_cast<int>(std::floor((real_box.hi(2) - geom.ProbLo(2))/geom.CellSize(2))););
477 
478  return amrex::Box(slice_lo, slice_hi) & geom.Domain();
479  }

Referenced by get_plane_mfs(), and write_plane_mfs().

Here is the caller graph for this function:

◆ write_plane_mfs()

void PlaneSampler::write_plane_mfs ( amrex::Vector< amrex::Real > &  time,
amrex::Vector< int > &  level_steps,
amrex::Vector< amrex::IntVect > &  ref_ratio,
amrex::Vector< amrex::Geometry > &  geom 
)
inline
548  {
549  amrex::Vector<std::string> varnames = {"T", "Wsp"};
550 
551  int nplane = m_ps_mf.size();
552  for (int iplane(0); iplane<nplane; ++iplane) {
553  // Data members that can be used as-is
554  int dir = m_dir[iplane];
555  int lev = m_lev[iplane];
556  amrex::Real m_time = time[lev];
557  amrex::Vector<int> m_level_steps = {level_steps[lev]};
558  amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
559 
560  // Create modified geometry object corresponding to the plane
561  amrex::RealBox m_rb = m_bnd_rbx[iplane];
562  amrex::Box m_dom = getIndexBox(m_rb, geom[lev]);
563  amrex::Real point = m_rb.hi(dir);
564  amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
565  for (int d(0); d<AMREX_SPACEDIM; ++d) {
566  if (d==dir) {
567  m_rb.setLo(d,point-0.5*geom[lev].CellSize(d));
568  m_rb.setHi(d,point+0.5*geom[lev].CellSize(d));
569  }
570  is_per[d] = geom[lev].isPeriodic(d);
571  }
572  amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
573  m_geom[0].define(m_dom, &m_rb, geom[lev].Coord(), is_per.data());
574 
575  // Create plotfile name
576  std::string name_plane = m_name[iplane];
577  name_plane += "_step_";
578  std::string plotfilename = amrex::Concatenate(name_plane, m_level_steps[0], 5);
579 
580  // Get the data
581  amrex::Vector<const amrex::MultiFab*> mf = {m_ps_mf[iplane].get()};
582 
583  // Write each plane
584  WriteMultiLevelPlotfile(plotfilename, 1, mf,
585  varnames, m_geom, m_time,
586  m_level_steps, m_ref_ratio);
587  } // iplane
588  }
Coord
Definition: ERF_DataStruct.H:68
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_plane_mfs().

◆ 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: