1 #ifndef ERF_SAMPLEDATA_H
2 #define ERF_SAMPLEDATA_H
6 #include <AMReX_ParmParse.H>
7 #include <AMReX_MultiFab.H>
8 #include <AMReX_MultiFabUtil.H>
9 #include <AMReX_PlotFileUtil.H>
21 amrex::ParmParse
pp(
"erf");
24 int n_line_lo =
pp.countval(
"sample_line_lo") / AMREX_SPACEDIM;
25 int n_line_hi =
pp.countval(
"sample_line_hi") / AMREX_SPACEDIM;
26 int n_line_dir =
pp.countval(
"sample_line_dir");
27 AMREX_ALWAYS_ASSERT( (n_line_lo==n_line_hi ) &&
28 (n_line_lo==n_line_dir) );
33 amrex::Vector<int> idx_lo; idx_lo.resize(n_line_lo*AMREX_SPACEDIM);
34 amrex::Vector<amrex::IntVect> iv_lo; iv_lo.resize(n_line_lo);
35 pp.queryarr(
"sample_line_lo",idx_lo,0,n_line_lo*AMREX_SPACEDIM);
36 for (
int i(0); i < n_line_lo; i++) {
37 amrex::IntVect iv(idx_lo[AMREX_SPACEDIM*i+0],
38 idx_lo[AMREX_SPACEDIM*i+1],
39 idx_lo[AMREX_SPACEDIM*i+2]);
44 amrex::Vector<int> idx_hi; idx_hi.resize(n_line_hi*AMREX_SPACEDIM);
45 amrex::Vector<amrex::IntVect> iv_hi; iv_hi.resize(n_line_hi);
46 pp.queryarr(
"sample_line_hi",idx_hi,0,n_line_hi*AMREX_SPACEDIM);
47 for (
int i(0); i < n_line_hi; i++) {
48 amrex::IntVect iv(idx_hi[AMREX_SPACEDIM*i+0],
49 idx_hi[AMREX_SPACEDIM*i+1],
50 idx_hi[AMREX_SPACEDIM*i+2]);
56 for (
int i = 0; i < n_line_hi; i++){
57 amrex::Box lbx(iv_lo[i],iv_hi[i]);
62 m_dir.resize(n_line_dir);
63 pp.queryarr(
"sample_line_dir",
m_dir,0,n_line_dir);
66 std::string name_base =
"plt_line_";
68 int n_names =
pp.countval(
"sample_line_name");
70 AMREX_ALWAYS_ASSERT( n_names==n_line_lo );
71 pp.queryarr(
"sample_line_name",
m_name,0,n_names);
73 for (
int iline(0); iline<n_line_lo; ++iline) {
74 m_name[iline] = amrex::Concatenate(name_base, iline , 5);
79 m_lev.resize(n_line_dir,0);
85 if (
pp.countval(
"line_sampling_vars") > 0) {
87 amrex::Vector<std::string> requested_vars;
88 pp.queryarr(
"line_sampling_vars",requested_vars);
89 amrex::Print() <<
"Selected line sampling vars :";
92 amrex::Print() <<
" " <<
"density";
96 amrex::Print() <<
" " <<
"x_velocity";
100 amrex::Print() <<
" " <<
"y_velocity";
104 amrex::Print() <<
" " <<
"z_velocity";
108 amrex::Print() <<
" " <<
"magvel";
112 amrex::Print() <<
" " <<
"theta";
116 amrex::Print() <<
" " <<
"qv";
120 amrex::Print() <<
" " <<
"qc";
124 amrex::Print() <<
" " <<
"pressure";
126 amrex::Print() << std::endl;
132 if (
m_write_ascii && amrex::ParallelDescriptor::IOProcessor()) {
133 int nvar =
static_cast<int>(
m_varnames.size());
136 for (
int iline(0); iline<n_line_lo; ++iline) {
137 for (
int ivar(0); ivar<nvar; ++ivar) {
142 amrex::FileOpenFailed(filename);
152 write_coords (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& z_phys_cc)
156 amrex::Print() <<
"Writing out line coordinates to text" << std::endl;
158 for (
int lev(0); lev < z_phys_cc.size(); ++lev) {
160 std::ofstream outfile;
161 if (amrex::ParallelDescriptor::IOProcessor()) {
162 std::string fname = amrex::Concatenate(
"plt_line_lev", lev, 1);
166 if (!outfile.is_open()) {
167 amrex::AllPrint() <<
"Could not open " << fname << std::endl;
172 int nline =
static_cast<int>(
m_ls_mf.size());
173 for (
int iline(0); iline<nline; ++iline) {
174 int dir =
m_dir[iline];
175 amrex::Box bnd_bx =
m_bnd_bx[iline];
176 amrex::IntVect first_cell = bnd_bx.smallEnd();
179 amrex::MultiFab line_coords_mf = get_line_data(
180 *z_phys_cc[lev], dir, first_cell, bnd_bx
184 amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(
185 line_coords_mf, 0, 1, bnd_bx, dir
189 if (amrex::ParallelDescriptor::IOProcessor()) {
190 for (
const auto& zval : vec) {
191 outfile <<
" " << zval;
193 outfile << std::endl;
203 amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_new)
205 int nlev =
static_cast<int>(vars_new.size());
206 int nline =
static_cast<int>(
m_bnd_bx.size());
207 int ncomp =
static_cast<int>(
m_varnames.size());
212 for (
int iline(0); iline<nline; ++iline) {
213 int dir =
m_dir[iline];
214 amrex::Box bnd_bx =
m_bnd_bx[iline];
215 amrex::IntVect cell = bnd_bx.smallEnd();
218 for (
int ilev(nlev-1); ilev>=0; --ilev) {
221 amrex::MultiFab mf_cc_vel;
222 auto ba = vars_new[ilev][
Vars::cons].boxArray();
223 auto dm = vars_new[ilev][
Vars::cons].DistributionMap();
224 mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
225 average_face_to_cellcenter(mf_cc_vel,0,
226 amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][
Vars::xvel],
231 amrex::MultiFab mf_cc_data;
232 mf_cc_data.define(ba, dm, ncomp, 1);
237 amrex::MultiFab::Copy(mf_cc_data, vars_new[ilev][
Vars::cons],
Rho_comp, mf_comp, 1, 0);
242 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 0, mf_comp, 1, 0);
246 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 1, mf_comp, 1, 0);
250 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 2, mf_comp, 1, 0);
256 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
258 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
259 const amrex::Box& tbx = mfi.tilebox();
260 auto const& dfab = mf_cc_data.array(mfi);
261 auto const& vfab = mf_cc_vel.array(mfi);
263 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
265 dfab(i,j,k,mf_comp) = std::sqrt(vfab(i,j,k,0)*vfab(i,j,k,0)
266 + vfab(i,j,k,1)*vfab(i,j,k,1)
267 + vfab(i,j,k,2)*vfab(i,j,k,2)) ;
275 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
277 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
278 const amrex::Box& tbx = mfi.tilebox();
279 auto const& dfab = mf_cc_data.array(mfi);
280 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
282 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
292 if (qv_comp >= 0) AMREX_ALWAYS_ASSERT(qv_comp == mf_comp);
295 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
297 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
298 const amrex::Box& tbx = mfi.tilebox();
299 auto const& dfab = mf_cc_data.array(mfi);
300 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
302 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
311 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
313 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
314 const amrex::Box& tbx = mfi.tilebox();
315 auto const& dfab = mf_cc_data.array(mfi);
316 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
318 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
329 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
331 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
332 const amrex::Box& tbx = mfi.tilebox();
333 auto const& dfab = mf_cc_data.array(mfi);
334 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
336 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
339 dfab(i,j,k,qv_comp));
347 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
349 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
350 const amrex::Box& tbx = mfi.tilebox();
351 auto const& dfab = mf_cc_data.array(mfi);
352 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
354 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
363 m_ls_mf[iline] = get_line_data(mf_cc_data, dir, cell, bnd_bx);
366 auto min_bnd_bx =
m_ls_mf[iline].boxArray().minimalBox();
367 if (bnd_bx == min_bnd_bx) {
break; }
375 amrex::Vector<int>& level_steps,
376 amrex::Vector<amrex::IntVect>& ref_ratio,
377 amrex::Vector<amrex::Geometry>& geom)
390 constexpr
int datwidth = 14;
391 constexpr
int datprecision = 6;
392 constexpr
int timeprecision = 13;
394 int nline =
static_cast<int>(
m_ls_mf.size());
395 int nvar =
static_cast<int>(
m_varnames.size());
396 for (
int iline(0); iline<nline; ++iline) {
397 int dir =
m_dir[iline];
398 int lev =
m_lev[iline];
402 for (
int ivar(0); ivar<nvar; ++ivar) {
404 amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(
m_ls_mf[iline], ivar, 1, m_dom, dir);
406 if (amrex::ParallelDescriptor::IOProcessor()) {
407 int ifile = iline*nvar + ivar;
409 fs << std::setw(datwidth) << std::setprecision(timeprecision) << m_time
410 << std::setw(datwidth) << std::setprecision(datprecision);
411 for (
const auto& val : vec) {
422 amrex::Vector<int>& level_steps,
423 amrex::Vector<amrex::IntVect>& ref_ratio,
424 amrex::Vector<amrex::Geometry>& geom)
426 int nline =
static_cast<int>(
m_ls_mf.size());
427 for (
int iline(0); iline<nline; ++iline) {
429 int dir =
m_dir[iline];
430 int lev =
m_lev[iline];
432 amrex::Vector<int> m_level_steps = {level_steps[lev]};
433 amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
436 auto plo = geom[lev].ProbLo();
437 auto dx = geom[lev].CellSize();
438 amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
439 amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
442 for (
int d(0); d<AMREX_SPACEDIM; ++d) {
450 is_per[d] = geom[lev].isPeriodic(d);
452 m_geom[0].define(m_dom, &m_rb, geom[lev].
Coord(), is_per.data());
455 std::string name_line =
m_name[iline];
456 name_line +=
"_step_";
457 std::string plotfilename = amrex::Concatenate(name_line, m_level_steps[0], 5);
460 amrex::Vector<const amrex::MultiFab*> mf = {&(
m_ls_mf[iline])};
463 WriteMultiLevelPlotfile(plotfilename, 1, mf,
465 m_level_steps, m_ref_ratio);
485 amrex::ParmParse
pp(
"erf");
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) );
495 if (n_plane_lo > 0) {
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]};
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]};
520 for (
int i(0); i < n_plane_hi; i++){
521 amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
526 m_dir.resize(n_plane_dir);
527 pp.queryarr(
"sample_plane_dir",
m_dir,0,n_plane_dir);
530 std::string name_base =
"plt_plane_";
531 m_name.resize(n_plane_lo);
532 int n_names =
pp.countval(
"sample_plane_name");
534 AMREX_ALWAYS_ASSERT( n_names==n_plane_lo );
535 pp.queryarr(
"sample_plane_name",
m_name,0,n_names);
537 for (
int iplane(0); iplane<n_plane_lo; ++iplane) {
538 m_name[iplane] = amrex::Concatenate(name_base, iplane , 5);
543 m_lev.resize(n_plane_dir,0);
553 const amrex::Geometry& geom) {
554 amrex::IntVect slice_lo, slice_hi;
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))););
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))););
564 return amrex::Box(slice_lo, slice_hi) & geom.Domain();
569 amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_new)
571 int nlev =
static_cast<int>(vars_new.size());
572 int nplane =
static_cast<int>(
m_bnd_rbx.size());
574 bool interpolate =
true;
577 for (
int iplane(0); iplane<nplane; ++iplane) {
578 int dir =
m_dir[iplane];
579 amrex::RealBox bnd_rbx =
m_bnd_rbx[iplane];
583 for (
int ilev(nlev-1); ilev>=0; --ilev) {
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],
596 amrex::MultiFab mf_cc_data;
597 mf_cc_data.define(ba, dm, ncomp, 1);
599 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
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
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)) ;
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);
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; }
631 amrex::Vector<int>& level_steps,
632 amrex::Vector<amrex::IntVect>& ref_ratio,
633 amrex::Vector<amrex::Geometry>& geom)
635 amrex::Vector<std::string> varnames = {
"T",
"Wsp"};
638 for (
int iplane(0); iplane<nplane; ++iplane) {
640 int dir =
m_dir[iplane];
641 int lev =
m_lev[iplane];
643 amrex::Vector<int> m_level_steps = {level_steps[lev]};
644 amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
650 amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
651 for (
int d(0); d<AMREX_SPACEDIM; ++d) {
653 m_rb.setLo(d,point-0.5*geom[lev].CellSize(d));
654 m_rb.setHi(d,point+0.5*geom[lev].CellSize(d));
656 is_per[d] = geom[lev].isPeriodic(d);
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());
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);
667 amrex::Vector<const amrex::MultiFab*> mf = {
m_ps_mf[iplane].get()};
670 WriteMultiLevelPlotfile(plotfilename, 1, mf,
671 varnames, m_geom, m_time,
672 m_level_steps, m_ref_ratio);
679 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
m_ps_mf;
bool containerHasElement(const V &iterable, const T &query)
Definition: ERF_Container.H:5
Coord
Definition: ERF_DataStruct.H:85
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:233
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ fs
Definition: ERF_AdvanceMorrison.cpp:117
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
Definition: ERF_SampleData.H:18
amrex::Vector< std::string > m_varnames
Definition: ERF_SampleData.H:476
amrex::Vector< amrex::MultiFab > m_ls_mf
Definition: ERF_SampleData.H:472
void get_sample_data(amrex::Vector< amrex::Geometry > &, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
Definition: ERF_SampleData.H:202
amrex::Vector< std::unique_ptr< std::fstream > > m_datastream
Definition: ERF_SampleData.H:477
void write_coords(amrex::Vector< std::unique_ptr< amrex::MultiFab > > &z_phys_cc)
Definition: ERF_SampleData.H:152
void write_line_plotfile(amrex::Vector< amrex::Real > &time, amrex::Vector< int > &level_steps, amrex::Vector< amrex::IntVect > &ref_ratio, amrex::Vector< amrex::Geometry > &geom)
Definition: ERF_SampleData.H:421
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:469
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:473
void write_line_ascii(amrex::Vector< amrex::Real > &time)
Definition: ERF_SampleData.H:387
LineSampler()
Definition: ERF_SampleData.H:19
amrex::Vector< amrex::Box > m_bnd_bx
Definition: ERF_SampleData.H:471
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)
Definition: ERF_SampleData.H:374
bool m_write_ascii
Definition: ERF_SampleData.H:475
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:470
Definition: ERF_SampleData.H:482
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)
Definition: ERF_SampleData.H:630
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:676
amrex::Box getIndexBox(const amrex::RealBox &real_box, const amrex::Geometry &geom)
Definition: ERF_SampleData.H:552
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
PlaneSampler()
Definition: ERF_SampleData.H:483
void get_sample_data(amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
Definition: ERF_SampleData.H:568