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");
23 bool has_line_idx =
pp.contains(
"sample_line_lo") ||
pp.contains(
"sample_line_hi");
24 bool has_line_real =
pp.contains(
"sample_line_lo_real") ||
pp.contains(
"sample_line_hi_real");
25 if (has_line_idx && has_line_real) {
26 amrex::Abort(
"Specify only one of erf.sample_line_lo/hi or erf.sample_line_lo_real/hi_real");
30 int n_line_lo =
pp.countval(
"sample_line_lo") / AMREX_SPACEDIM;
31 int n_line_hi =
pp.countval(
"sample_line_hi") / AMREX_SPACEDIM;
32 int n_line_lo_real =
pp.countval(
"sample_line_lo_real") / AMREX_SPACEDIM;
33 int n_line_hi_real =
pp.countval(
"sample_line_hi_real") / AMREX_SPACEDIM;
34 int n_line_dir =
pp.countval(
"sample_line_dir");
37 (n_line_lo==n_line_dir) );
38 }
else if (has_line_real) {
40 (n_line_lo_real==n_line_dir) );
49 amrex::Vector<int> idx_lo; idx_lo.resize(nline*AMREX_SPACEDIM);
50 amrex::Vector<amrex::IntVect> iv_lo; iv_lo.resize(nline);
51 pp.queryarr(
"sample_line_lo",idx_lo,0,nline*AMREX_SPACEDIM);
52 for (
int i(0); i < nline; i++) {
53 amrex::IntVect iv(idx_lo[AMREX_SPACEDIM*i+0],
54 idx_lo[AMREX_SPACEDIM*i+1],
55 idx_lo[AMREX_SPACEDIM*i+2]);
60 amrex::Vector<int> idx_hi; idx_hi.resize(nline*AMREX_SPACEDIM);
61 amrex::Vector<amrex::IntVect> iv_hi; iv_hi.resize(nline);
62 pp.queryarr(
"sample_line_hi",idx_hi,0,nline*AMREX_SPACEDIM);
63 for (
int i(0); i < nline; i++) {
64 amrex::IntVect iv(idx_hi[AMREX_SPACEDIM*i+0],
65 idx_hi[AMREX_SPACEDIM*i+1],
66 idx_hi[AMREX_SPACEDIM*i+2]);
72 for (
int i = 0; i < nline; i++){
73 amrex::Box lbx(iv_lo[i],iv_hi[i]);
77 }
else if (has_line_real) {
79 nline = n_line_lo_real;
82 amrex::Vector<amrex::Real> real_lo; real_lo.resize(nline*AMREX_SPACEDIM);
83 amrex::Vector<amrex::Vector<amrex::Real>> rv_lo;
84 pp.queryarr(
"sample_line_lo_real",real_lo,0,nline*AMREX_SPACEDIM);
85 for (
int i(0); i < nline; i++) {
86 amrex::Vector<amrex::Real> rv = {real_lo[AMREX_SPACEDIM*i+0],
87 real_lo[AMREX_SPACEDIM*i+1],
88 real_lo[AMREX_SPACEDIM*i+2]};
93 amrex::Vector<amrex::Real> real_hi; real_hi.resize(nline*AMREX_SPACEDIM);
94 amrex::Vector<amrex::Vector<amrex::Real>> rv_hi;
95 pp.queryarr(
"sample_line_hi_real",real_hi,0,nline*AMREX_SPACEDIM);
96 for (
int i(0); i < nline; i++) {
97 amrex::Vector<amrex::Real> rv = {real_hi[AMREX_SPACEDIM*i+0],
98 real_hi[AMREX_SPACEDIM*i+1],
99 real_hi[AMREX_SPACEDIM*i+2]};
105 for (
int i = 0; i < nline; i++){
106 amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
117 m_dir.resize(n_line_dir);
118 pp.queryarr(
"sample_line_dir",
m_dir,0,n_line_dir);
121 std::string name_base =
"plt_line_";
123 int n_names =
pp.countval(
"sample_line_name");
126 pp.queryarr(
"sample_line_name",
m_name,0,n_names);
128 for (
int iline(0); iline<nline; ++iline) {
129 m_name[iline] = amrex::Concatenate(name_base, iline , 5);
134 m_lev.resize(n_line_dir,0);
140 if (
pp.countval(
"line_sampling_vars") > 0) {
142 amrex::Vector<std::string> requested_vars;
143 pp.queryarr(
"line_sampling_vars",requested_vars);
144 amrex::Print() <<
"Selected line sampling vars :";
147 amrex::Print() <<
" " <<
"density";
151 amrex::Print() <<
" " <<
"x_velocity";
155 amrex::Print() <<
" " <<
"y_velocity";
159 amrex::Print() <<
" " <<
"z_velocity";
163 amrex::Print() <<
" " <<
"magvel";
167 amrex::Print() <<
" " <<
"theta";
171 amrex::Print() <<
" " <<
"qv";
175 amrex::Print() <<
" " <<
"qc";
179 amrex::Print() <<
" " <<
"pressure";
181 amrex::Print() << std::endl;
187 if (
m_write_ascii && amrex::ParallelDescriptor::IOProcessor()) {
188 int nvar =
static_cast<int>(
m_varnames.size());
191 for (
int iline(0); iline<nline; ++iline) {
192 for (
int ivar(0); ivar<nvar; ++ivar) {
197 amrex::FileOpenFailed(filename);
208 const amrex::Geometry& geom) {
209 amrex::IntVect slice_lo, slice_hi;
211 AMREX_D_TERM(slice_lo[0]=
static_cast<int>(std::floor((real_box.lo(0) - geom.ProbLo(0))/geom.CellSize(0)));,
212 slice_lo[1]=
static_cast<int>(std::floor((real_box.lo(1) - geom.ProbLo(1))/geom.CellSize(1)));,
213 slice_lo[2]=
static_cast<int>(std::floor((real_box.lo(2) - geom.ProbLo(2))/geom.CellSize(2))););
215 AMREX_D_TERM(slice_hi[0]=
static_cast<int>(std::floor((real_box.hi(0) - geom.ProbLo(0))/geom.CellSize(0)));,
216 slice_hi[1]=
static_cast<int>(std::floor((real_box.hi(1) - geom.ProbLo(1))/geom.CellSize(1)));,
217 slice_hi[2]=
static_cast<int>(std::floor((real_box.hi(2) - geom.ProbLo(2))/geom.CellSize(2))););
219 return amrex::Box(slice_lo, slice_hi) & geom.Domain();
223 write_coords (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& z_phys_cc,
224 amrex::Vector<amrex::Geometry>& geom)
228 amrex::Print() <<
"Writing out line coordinates to text" << std::endl;
230 for (
int lev(0); lev < z_phys_cc.size(); ++lev) {
232 std::ofstream outfile;
233 if (amrex::ParallelDescriptor::IOProcessor()) {
234 std::string fname = amrex::Concatenate(
"plt_line_lev", lev, 1);
238 if (!outfile.is_open()) {
239 amrex::AllPrint() <<
"Could not open " << fname << std::endl;
244 int nline =
static_cast<int>(
m_ls_mf.size());
245 for (
int iline(0); iline<nline; ++iline) {
246 int dir =
m_dir[iline];
250 amrex::IntVect first_cell = bnd_bx.smallEnd();
253 amrex::MultiFab line_coords_mf = get_line_data(
254 *z_phys_cc[lev], dir, first_cell, bnd_bx
258 amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(
259 line_coords_mf, 0, 1, bnd_bx, dir
263 if (amrex::ParallelDescriptor::IOProcessor()) {
264 for (
const auto& zval : vec) {
265 outfile <<
" " << zval;
267 outfile << std::endl;
277 amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_new)
279 int nlev =
static_cast<int>(vars_new.size());
280 int nline =
static_cast<int>(
m_bnd_bx.size());
281 int ncomp =
static_cast<int>(
m_varnames.size());
286 for (
int iline(0); iline<nline; ++iline) {
287 int dir =
m_dir[iline];
290 for (
int ilev(nlev-1); ilev>=0; --ilev) {
294 amrex::IntVect cell = bnd_bx.smallEnd();
297 amrex::MultiFab mf_cc_vel;
298 auto ba = vars_new[ilev][
Vars::cons].boxArray();
299 auto dm = vars_new[ilev][
Vars::cons].DistributionMap();
300 mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
301 average_face_to_cellcenter(mf_cc_vel,0,
302 amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][
Vars::xvel],
307 amrex::MultiFab mf_cc_data;
308 mf_cc_data.define(ba, dm, ncomp, 1);
313 amrex::MultiFab::Copy(mf_cc_data, vars_new[ilev][
Vars::cons],
Rho_comp, mf_comp, 1, 0);
318 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 0, mf_comp, 1, 0);
322 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 1, mf_comp, 1, 0);
326 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 2, mf_comp, 1, 0);
332 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
334 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
335 const amrex::Box& tbx = mfi.tilebox();
336 auto const& dfab = mf_cc_data.array(mfi);
337 auto const& vfab = mf_cc_vel.array(mfi);
341 dfab(i,j,k,mf_comp) = std::sqrt(vfab(i,j,k,0)*vfab(i,j,k,0)
342 + vfab(i,j,k,1)*vfab(i,j,k,1)
343 + vfab(i,j,k,2)*vfab(i,j,k,2)) ;
351 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
353 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
354 const amrex::Box& tbx = mfi.tilebox();
355 auto const& dfab = mf_cc_data.array(mfi);
356 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
371 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
373 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
374 const amrex::Box& tbx = mfi.tilebox();
375 auto const& dfab = mf_cc_data.array(mfi);
376 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
387 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
389 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
390 const amrex::Box& tbx = mfi.tilebox();
391 auto const& dfab = mf_cc_data.array(mfi);
392 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
405 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
407 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
408 const amrex::Box& tbx = mfi.tilebox();
409 auto const& dfab = mf_cc_data.array(mfi);
410 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
415 dfab(i,j,k,qv_comp));
423 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
425 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
426 const amrex::Box& tbx = mfi.tilebox();
427 auto const& dfab = mf_cc_data.array(mfi);
428 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
439 m_ls_mf[iline] = get_line_data(mf_cc_data, dir, cell, bnd_bx);
442 auto min_bnd_bx =
m_ls_mf[iline].boxArray().minimalBox();
443 if (bnd_bx == min_bnd_bx) {
break; }
457 amrex::Vector<int>& level_steps,
458 amrex::Vector<amrex::IntVect>& ref_ratio,
459 amrex::Vector<amrex::Geometry>& geom)
472 constexpr
int datwidth = 14;
473 constexpr
int datprecision = 6;
474 constexpr
int timeprecision = 13;
476 int nline =
static_cast<int>(
m_ls_mf.size());
477 int nvar =
static_cast<int>(
m_varnames.size());
478 for (
int iline(0); iline<nline; ++iline) {
479 int dir =
m_dir[iline];
480 int lev =
m_lev[iline];
484 for (
int ivar(0); ivar<nvar; ++ivar) {
486 amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(
m_ls_mf[iline], ivar, 1, m_dom, dir);
488 if (amrex::ParallelDescriptor::IOProcessor()) {
489 int ifile = iline*nvar + ivar;
491 fs << std::setw(datwidth) << std::setprecision(timeprecision) << m_time
492 << std::setw(datwidth) << std::setprecision(datprecision);
493 for (
const auto& val : vec) {
504 amrex::Vector<int>& level_steps,
505 amrex::Vector<amrex::IntVect>& ref_ratio,
506 amrex::Vector<amrex::Geometry>& geom)
508 int nline =
static_cast<int>(
m_ls_mf.size());
509 for (
int iline(0); iline<nline; ++iline) {
511 int dir =
m_dir[iline];
512 int lev =
m_lev[iline];
514 amrex::Vector<int> m_level_steps = {level_steps[lev]};
515 amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
518 auto plo = geom[lev].ProbLo();
519 auto dx = geom[lev].CellSize();
520 amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
521 amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
524 for (
int d(0); d<AMREX_SPACEDIM; ++d) {
532 is_per[d] = geom[lev].isPeriodic(d);
534 m_geom[0].define(m_dom, &m_rb, geom[lev].
Coord(), is_per.data());
537 std::string name_line =
m_name[iline];
538 name_line +=
"_step_";
539 std::string plotfilename = amrex::Concatenate(name_line, m_level_steps[0], 5);
542 amrex::Vector<const amrex::MultiFab*> mf = {&(
m_ls_mf[iline])};
545 WriteMultiLevelPlotfile(plotfilename, 1, mf,
547 m_level_steps, m_ref_ratio);
569 amrex::ParmParse
pp(
"erf");
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");
576 (n_plane_lo==n_plane_dir) );
579 if (n_plane_lo > 0) {
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]};
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]};
604 for (
int i(0); i < n_plane_hi; i++){
605 amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
610 m_dir.resize(n_plane_dir);
611 pp.queryarr(
"sample_plane_dir",
m_dir,0,n_plane_dir);
614 std::string name_base =
"plt_plane_";
615 m_name.resize(n_plane_lo);
616 int n_names =
pp.countval(
"sample_plane_name");
619 pp.queryarr(
"sample_plane_name",
m_name,0,n_names);
621 for (
int iplane(0); iplane<n_plane_lo; ++iplane) {
622 m_name[iplane] = amrex::Concatenate(name_base, iplane , 5);
627 m_lev.resize(n_plane_dir,0);
637 const amrex::Geometry& geom) {
638 amrex::IntVect slice_lo, slice_hi;
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))););
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))););
648 return amrex::Box(slice_lo, slice_hi) & geom.Domain();
653 amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_new)
655 int nlev =
static_cast<int>(vars_new.size());
656 int nplane =
static_cast<int>(
m_bnd_rbx.size());
658 bool interpolate =
true;
661 for (
int iplane(0); iplane<nplane; ++iplane) {
662 int dir =
m_dir[iplane];
663 amrex::RealBox bnd_rbx =
m_bnd_rbx[iplane];
667 for (
int ilev(nlev-1); ilev>=0; --ilev) {
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],
680 amrex::MultiFab mf_cc_data;
681 mf_cc_data.define(ba, dm, ncomp, 1);
683 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
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);
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)) ;
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);
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; }
715 amrex::Vector<int>& level_steps,
716 amrex::Vector<amrex::IntVect>& ref_ratio,
717 amrex::Vector<amrex::Geometry>& geom)
719 amrex::Vector<std::string> varnames = {
"T",
"Wsp"};
722 for (
int iplane(0); iplane<nplane; ++iplane) {
724 int dir =
m_dir[iplane];
725 int lev =
m_lev[iplane];
727 amrex::Vector<int> m_level_steps = {level_steps[lev]};
728 amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
734 amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
735 for (
int d(0); d<AMREX_SPACEDIM; ++d) {
737 m_rb.setLo(d,point-
myhalf*geom[lev].CellSize(d));
738 m_rb.setHi(d,point+
myhalf*geom[lev].CellSize(d));
740 is_per[d] = geom[lev].isPeriodic(d);
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());
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);
751 amrex::Vector<const amrex::MultiFab*> mf = {
m_ps_mf[iplane].get()};
754 WriteMultiLevelPlotfile(plotfilename, 1, mf,
755 varnames, m_geom, m_time,
756 m_level_steps, m_ref_ratio);
763 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
m_ps_mf;
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
bool containerHasElement(const V &iterable, const T &query)
Definition: ERF_Container.H:5
Coord
Definition: ERF_DataStruct.H:92
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(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
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
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_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:119
@ xvel
Definition: ERF_IndexDefines.H:159
@ cons
Definition: ERF_IndexDefines.H:158
@ zvel
Definition: ERF_IndexDefines.H:161
@ yvel
Definition: ERF_IndexDefines.H:160
Definition: ERF_SampleData.H:18
amrex::Vector< std::string > m_varnames
Definition: ERF_SampleData.H:560
amrex::Vector< amrex::MultiFab > m_ls_mf
Definition: ERF_SampleData.H:555
void write_coords(amrex::Vector< std::unique_ptr< amrex::MultiFab > > &z_phys_cc, amrex::Vector< amrex::Geometry > &geom)
Definition: ERF_SampleData.H:223
amrex::Vector< std::unique_ptr< std::fstream > > m_datastream
Definition: ERF_SampleData.H:561
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:503
bool m_use_real_bx
Definition: ERF_SampleData.H:558
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:551
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:556
void write_line_ascii(amrex::Vector< amrex::Real > &time)
Definition: ERF_SampleData.H:469
LineSampler()
Definition: ERF_SampleData.H:19
amrex::Vector< amrex::Box > m_bnd_bx
Definition: ERF_SampleData.H:553
void get_sample_data(amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
Definition: ERF_SampleData.H:276
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:456
amrex::Box getIndexBox(const amrex::RealBox &real_box, const amrex::Geometry &geom)
Definition: ERF_SampleData.H:207
bool m_write_ascii
Definition: ERF_SampleData.H:559
amrex::Vector< amrex::RealBox > m_bnd_rbx
Definition: ERF_SampleData.H:554
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:552
Definition: ERF_SampleData.H:566
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:714
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:760
amrex::Box getIndexBox(const amrex::RealBox &real_box, const amrex::Geometry &geom)
Definition: ERF_SampleData.H:636
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
PlaneSampler()
Definition: ERF_SampleData.H:567
void get_sample_data(amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
Definition: ERF_SampleData.H:652