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

#include <ERF_SampleData.H>

Collaboration diagram for LineSampler:

Public Member Functions

 LineSampler ()
 
void get_line_mfs (amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
 
void write_line_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::Box > m_bnd_bx
 
amrex::Vector< amrex::MultiFab > m_ls_mf
 

Constructor & Destructor Documentation

◆ LineSampler()

LineSampler::LineSampler ( )
inline
16  {
17  amrex::ParmParse pp("erf");
18 
19  // Count number of lo and hi points define the line
20  int n_line_lo = pp.countval("sample_line_lo") / AMREX_SPACEDIM;
21  int n_line_hi = pp.countval("sample_line_hi") / AMREX_SPACEDIM;
22  int n_line_dir = pp.countval("sample_line_dir");
23  AMREX_ALWAYS_ASSERT( (n_line_lo==n_line_hi ) &&
24  (n_line_lo==n_line_dir) );
25 
26  // Parse the data
27  if (n_line_lo > 0) {
28  // Parse lo
29  amrex::Vector<int> idx_lo; idx_lo.resize(n_line_lo*AMREX_SPACEDIM);
30  amrex::Vector<amrex::IntVect> iv_lo; iv_lo.resize(n_line_lo);
31  pp.queryarr("sample_line_lo",idx_lo,0,n_line_lo*AMREX_SPACEDIM);
32  for (int i(0); i < n_line_lo; i++) {
33  amrex::IntVect iv(idx_lo[AMREX_SPACEDIM*i+0],
34  idx_lo[AMREX_SPACEDIM*i+1],
35  idx_lo[AMREX_SPACEDIM*i+2]);
36  iv_lo[i] = iv;
37  }
38 
39  // Parse hi
40  amrex::Vector<int> idx_hi; idx_hi.resize(n_line_hi*AMREX_SPACEDIM);
41  amrex::Vector<amrex::IntVect> iv_hi; iv_hi.resize(n_line_hi);
42  pp.queryarr("sample_line_hi",idx_hi,0,n_line_hi*AMREX_SPACEDIM);
43  for (int i(0); i < n_line_hi; i++) {
44  amrex::IntVect iv(idx_hi[AMREX_SPACEDIM*i+0],
45  idx_hi[AMREX_SPACEDIM*i+1],
46  idx_hi[AMREX_SPACEDIM*i+2]);
47  iv_hi[i] = iv;
48  }
49 
50  // Construct vector of bounding boxes
51  m_bnd_bx.resize(n_line_lo);
52  for (int i = 0; i < n_line_hi; i++){
53  amrex::Box lbx(iv_lo[i],iv_hi[i]);
54  m_bnd_bx[i] = lbx;
55  }
56 
57  // Parse directionality
58  m_dir.resize(n_line_dir);
59  pp.queryarr("sample_line_dir",m_dir,0,n_line_dir);
60 
61  // Allocate space for level indicator
62  m_lev.resize(n_line_dir,0);
63 
64  // Allocate space for MF pointers
65  m_ls_mf.resize(n_line_lo);
66  }
67  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:219
amrex::Vector< amrex::MultiFab > m_ls_mf
Definition: ERF_SampleData.H:182
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:179
amrex::Vector< amrex::Box > m_bnd_bx
Definition: ERF_SampleData.H:181
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:180
Here is the call graph for this function:

Member Function Documentation

◆ get_line_mfs()

void LineSampler::get_line_mfs ( amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars_new)
inline
71  {
72  int nlev = vars_new.size();
73  int nline = m_bnd_bx.size();
74  int ncomp = 2;
75 
76  // Loop over each line
77  for (int iline(0); iline<nline; ++iline) {
78  int dir = m_dir[iline];
79  amrex::Box bnd_bx = m_bnd_bx[iline];
80  amrex::IntVect cell = bnd_bx.smallEnd();
81 
82  // Search each level to get the finest data possible
83  for (int ilev(nlev-1); ilev>=0; --ilev) {
84 
85  // Construct CC velocities
86  amrex::MultiFab mf_cc_vel;
87  auto ba = vars_new[ilev][Vars::cons].boxArray();
88  auto dm = vars_new[ilev][Vars::cons].DistributionMap();
89  mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
90  average_face_to_cellcenter(mf_cc_vel,0,
91  amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][Vars::xvel],
92  &vars_new[ilev][Vars::yvel],
93  &vars_new[ilev][Vars::zvel]});
94 
95  // Construct vector of MFs holding T and WSP
96  amrex::MultiFab mf_cc_data;
97  mf_cc_data.define(ba, dm, ncomp, 1);
98 #ifdef _OPENMP
99 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
100 #endif
101  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
102  const amrex::Box& tbx = mfi.tilebox();
103  auto const& dfab = mf_cc_data.array(mfi);
104  auto const& tfab = vars_new[ilev][Vars::cons].array(mfi);
105  auto const& wfab = mf_cc_vel.array(mfi);
106  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
107  {
108  dfab(i,j,k,0) = tfab(i,j,k,1)/tfab(i,j,k,0);
109  dfab(i,j,k,1) = std::sqrt(wfab(i,j,k,0)*wfab(i,j,k,0)
110  + wfab(i,j,k,1)*wfab(i,j,k,1)
111  + wfab(i,j,k,2)*wfab(i,j,k,2)) ;
112  });
113 
114  }
115 
116  m_lev[iline] = ilev;
117  m_ls_mf[iline] = get_line_data(mf_cc_data, dir, cell, bnd_bx);
118 
119  // We can stop if we got the entire line
120  auto min_bnd_bx = m_ls_mf[iline].boxArray().minimalBox();
121  if (bnd_bx == min_bnd_bx) { continue; }
122 
123  } // ilev
124  }// iline
125  }
@ xvel
Definition: ERF_IndexDefines.H:130
@ cons
Definition: ERF_IndexDefines.H:129
@ zvel
Definition: ERF_IndexDefines.H:132
@ yvel
Definition: ERF_IndexDefines.H:131

◆ write_line_mfs()

void LineSampler::write_line_mfs ( amrex::Vector< amrex::Real > &  time,
amrex::Vector< int > &  level_steps,
amrex::Vector< amrex::IntVect > &  ref_ratio,
amrex::Vector< amrex::Geometry > &  geom 
)
inline
132  {
133  std::string name_base = "plt_line_";
134  amrex::Vector<std::string> varnames = {"T", "Wsp"};
135 
136  int nline = m_ls_mf.size();
137  for (int iline(0); iline<nline; ++iline) {
138  // Data members that can be used as-is
139  int dir = m_dir[iline];
140  int lev = m_lev[iline];
141  amrex::Real m_time = time[lev];
142  amrex::Vector<int> m_level_steps = {level_steps[lev]};
143  amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
144 
145  // Create modified geometry object corresponding to the line
146  auto plo = geom[lev].ProbLo();
147  auto dx = geom[lev].CellSize();
148  amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
149  amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
150  amrex::Box m_dom = m_bnd_bx[iline];
151  amrex::RealBox m_rb;
152  for (int d(0); d<AMREX_SPACEDIM; ++d) {
153  amrex::Real offset = (d==dir) ? 0 : 0.5;
154  amrex::Real lo = plo[d] + ( m_dom.smallEnd(d) - offset ) * dx[d];
155  amrex::Real hi = plo[d] + ( m_dom.bigEnd(d) + offset ) * dx[d];
156 
157  m_rb.setLo(d,lo);
158  m_rb.setHi(d,hi);
159 
160  is_per[d] = geom[lev].isPeriodic(d);
161  }
162  m_geom[0].define(m_dom, &m_rb, geom[lev].Coord(), is_per.data());
163 
164  // Create plotfile name
165  std::string name_line = amrex::Concatenate(name_base, iline , 5);
166  name_line += "_step_";
167  std::string plotfilename = amrex::Concatenate(name_line, m_level_steps[0], 5);
168 
169  // Get the data
170  amrex::Vector<const amrex::MultiFab*> mf = {&(m_ls_mf[iline])};
171 
172  // Write each line
173  WriteMultiLevelPlotfile(plotfilename, 1, mf,
174  varnames, m_geom, m_time,
175  m_level_steps, m_ref_ratio);
176  }
177  }
Coord
Definition: ERF_DataStruct.H:64
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
Here is the call graph for this function:

Member Data Documentation

◆ m_bnd_bx

amrex::Vector<amrex::Box> LineSampler::m_bnd_bx

◆ m_dir

amrex::Vector<int> LineSampler::m_dir

◆ m_lev

amrex::Vector<int> LineSampler::m_lev

◆ m_ls_mf

amrex::Vector<amrex::MultiFab> LineSampler::m_ls_mf

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