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);
368 "qv sampling requested but moisture components not present in state");
372 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
374 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
375 const amrex::Box& tbx = mfi.tilebox();
376 auto const& dfab = mf_cc_data.array(mfi);
377 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
388 "qc sampling requested but moisture components not present in state");
390 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
392 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
393 const amrex::Box& tbx = mfi.tilebox();
394 auto const& dfab = mf_cc_data.array(mfi);
395 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
409 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
411 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
412 const amrex::Box& tbx = mfi.tilebox();
413 auto const& dfab = mf_cc_data.array(mfi);
414 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
418 amrex::Real qv_val = (qv_comp >= 0) ? dfab(i,j,k,qv_comp)
426 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
428 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
429 const amrex::Box& tbx = mfi.tilebox();
430 auto const& dfab = mf_cc_data.array(mfi);
431 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
443 m_ls_mf[iline] = get_line_data(mf_cc_data, dir, cell, bnd_bx);
446 auto min_bnd_bx =
m_ls_mf[iline].boxArray().minimalBox();
447 if (bnd_bx == min_bnd_bx) {
break; }
461 amrex::Vector<int>& level_steps,
462 amrex::Vector<amrex::IntVect>& ref_ratio,
463 amrex::Vector<amrex::Geometry>& geom)
476 constexpr
int datwidth = 14;
477 constexpr
int datprecision = 6;
478 constexpr
int timeprecision = 13;
480 int nline =
static_cast<int>(
m_ls_mf.size());
481 int nvar =
static_cast<int>(
m_varnames.size());
482 for (
int iline(0); iline<nline; ++iline) {
483 int dir =
m_dir[iline];
484 int lev =
m_lev[iline];
488 for (
int ivar(0); ivar<nvar; ++ivar) {
490 amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(
m_ls_mf[iline], ivar, 1, m_dom, dir);
492 if (amrex::ParallelDescriptor::IOProcessor()) {
493 int ifile = iline*nvar + ivar;
495 fs << std::setw(datwidth) << std::setprecision(timeprecision) << m_time
496 << std::setw(datwidth) << std::setprecision(datprecision);
497 for (
const auto& val : vec) {
508 amrex::Vector<int>& level_steps,
509 amrex::Vector<amrex::IntVect>& ref_ratio,
510 amrex::Vector<amrex::Geometry>& geom)
512 int nline =
static_cast<int>(
m_ls_mf.size());
513 for (
int iline(0); iline<nline; ++iline) {
515 int dir =
m_dir[iline];
516 int lev =
m_lev[iline];
518 amrex::Vector<int> m_level_steps = {level_steps[lev]};
519 amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
522 auto plo = geom[lev].ProbLo();
523 auto dx = geom[lev].CellSize();
524 amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
525 amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
528 for (
int d(0); d<AMREX_SPACEDIM; ++d) {
536 is_per[d] = geom[lev].isPeriodic(d);
538 m_geom[0].define(m_dom, &m_rb, geom[lev].
Coord(), is_per.data());
541 std::string name_line =
m_name[iline];
542 name_line +=
"_step_";
543 std::string plotfilename = amrex::Concatenate(name_line, m_level_steps[0], 5);
546 amrex::Vector<const amrex::MultiFab*> mf = {&(
m_ls_mf[iline])};
549 WriteMultiLevelPlotfile(plotfilename, 1, mf,
551 m_level_steps, m_ref_ratio);
573 amrex::ParmParse
pp(
"erf");
576 int n_plane_lo =
pp.countval(
"sample_plane_lo") / AMREX_SPACEDIM;
577 int n_plane_hi =
pp.countval(
"sample_plane_hi") / AMREX_SPACEDIM;
578 int n_plane_dir =
pp.countval(
"sample_plane_dir");
580 (n_plane_lo==n_plane_dir) );
583 if (n_plane_lo > 0) {
585 amrex::Vector<amrex::Real> r_lo; r_lo.resize(n_plane_lo*AMREX_SPACEDIM);
586 amrex::Vector<amrex::Vector<amrex::Real>> rv_lo;
587 pp.queryarr(
"sample_plane_lo",r_lo,0,n_plane_lo*AMREX_SPACEDIM);
588 for (
int i(0); i < n_plane_lo; i++) {
589 amrex::Vector<amrex::Real> rv = {r_lo[AMREX_SPACEDIM*i+0],
590 r_lo[AMREX_SPACEDIM*i+1],
591 r_lo[AMREX_SPACEDIM*i+2]};
596 amrex::Vector<amrex::Real> r_hi; r_hi.resize(n_plane_hi*AMREX_SPACEDIM);
597 amrex::Vector<amrex::Vector<amrex::Real>> rv_hi;
598 pp.queryarr(
"sample_plane_hi",r_hi,0,n_plane_hi*AMREX_SPACEDIM);
599 for (
int i(0); i < n_plane_hi; i++) {
600 amrex::Vector<amrex::Real> rv = {r_hi[AMREX_SPACEDIM*i+0],
601 r_hi[AMREX_SPACEDIM*i+1],
602 r_hi[AMREX_SPACEDIM*i+2]};
608 for (
int i(0); i < n_plane_hi; i++){
609 amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
614 m_dir.resize(n_plane_dir);
615 pp.queryarr(
"sample_plane_dir",
m_dir,0,n_plane_dir);
618 std::string name_base =
"plt_plane_";
619 m_name.resize(n_plane_lo);
620 int n_names =
pp.countval(
"sample_plane_name");
623 pp.queryarr(
"sample_plane_name",
m_name,0,n_names);
625 for (
int iplane(0); iplane<n_plane_lo; ++iplane) {
626 m_name[iplane] = amrex::Concatenate(name_base, iplane , 5);
631 m_lev.resize(n_plane_dir,0);
637 if (
pp.countval(
"plane_sampling_vars") > 0) {
639 amrex::Vector<std::string> requested_vars;
640 pp.queryarr(
"plane_sampling_vars",requested_vars);
641 amrex::Print() <<
"Selected plane sampling vars :";
644 amrex::Print() <<
" " <<
"density";
648 amrex::Print() <<
" " <<
"x_velocity";
652 amrex::Print() <<
" " <<
"y_velocity";
656 amrex::Print() <<
" " <<
"z_velocity";
660 amrex::Print() <<
" " <<
"magvel";
664 amrex::Print() <<
" " <<
"theta";
668 amrex::Print() <<
" " <<
"qv";
672 amrex::Print() <<
" " <<
"qc";
676 amrex::Print() <<
" " <<
"pressure";
678 amrex::Print() << std::endl;
686 const amrex::Geometry& geom) {
687 amrex::IntVect slice_lo, slice_hi;
689 AMREX_D_TERM(slice_lo[0]=
static_cast<int>(std::floor((real_box.lo(0) - geom.ProbLo(0))/geom.CellSize(0)));,
690 slice_lo[1]=
static_cast<int>(std::floor((real_box.lo(1) - geom.ProbLo(1))/geom.CellSize(1)));,
691 slice_lo[2]=
static_cast<int>(std::floor((real_box.lo(2) - geom.ProbLo(2))/geom.CellSize(2))););
693 AMREX_D_TERM(slice_hi[0]=
static_cast<int>(std::floor((real_box.hi(0) - geom.ProbLo(0))/geom.CellSize(0)));,
694 slice_hi[1]=
static_cast<int>(std::floor((real_box.hi(1) - geom.ProbLo(1))/geom.CellSize(1)));,
695 slice_hi[2]=
static_cast<int>(std::floor((real_box.hi(2) - geom.ProbLo(2))/geom.CellSize(2))););
697 return amrex::Box(slice_lo, slice_hi) & geom.Domain();
702 amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_new)
704 int nlev =
static_cast<int>(vars_new.size());
705 int nplane =
static_cast<int>(
m_bnd_rbx.size());
706 int ncomp =
static_cast<int>(
m_varnames.size());
707 bool interpolate =
true;
712 for (
int iplane(0); iplane<nplane; ++iplane) {
713 int dir =
m_dir[iplane];
714 amrex::RealBox bnd_rbx =
m_bnd_rbx[iplane];
718 for (
int ilev(nlev-1); ilev>=0; --ilev) {
721 amrex::MultiFab mf_cc_vel;
722 auto ba = vars_new[ilev][
Vars::cons].boxArray();
723 auto dm = vars_new[ilev][
Vars::cons].DistributionMap();
724 mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
725 average_face_to_cellcenter(mf_cc_vel,0,
726 amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][
Vars::xvel],
731 amrex::MultiFab mf_cc_data;
732 mf_cc_data.define(ba, dm, ncomp, 1);
737 amrex::MultiFab::Copy(mf_cc_data, vars_new[ilev][
Vars::cons],
Rho_comp, mf_comp, 1, 0);
742 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 0, mf_comp, 1, 0);
746 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 1, mf_comp, 1, 0);
750 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 2, mf_comp, 1, 0);
756 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
758 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
759 const amrex::Box& tbx = mfi.tilebox();
760 auto const& dfab = mf_cc_data.array(mfi);
761 auto const& vfab = mf_cc_vel.array(mfi);
765 dfab(i,j,k,mf_comp) = std::sqrt(vfab(i,j,k,0)*vfab(i,j,k,0)
766 + vfab(i,j,k,1)*vfab(i,j,k,1)
767 + vfab(i,j,k,2)*vfab(i,j,k,2));
775 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
777 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
778 const amrex::Box& tbx = mfi.tilebox();
779 auto const& dfab = mf_cc_data.array(mfi);
780 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
792 "qv sampling requested but moisture components not present in state");
796 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
798 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
799 const amrex::Box& tbx = mfi.tilebox();
800 auto const& dfab = mf_cc_data.array(mfi);
801 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
812 "qc sampling requested but moisture components not present in state");
814 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
816 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
817 const amrex::Box& tbx = mfi.tilebox();
818 auto const& dfab = mf_cc_data.array(mfi);
819 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
833 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
835 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
836 const amrex::Box& tbx = mfi.tilebox();
837 auto const& dfab = mf_cc_data.array(mfi);
838 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
842 amrex::Real qv_val = (qv_comp >= 0) ? dfab(i,j,k,qv_comp)
850 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
852 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
853 const amrex::Box& tbx = mfi.tilebox();
854 auto const& dfab = mf_cc_data.array(mfi);
855 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
866 m_lev[iplane] = ilev;
867 m_ps_mf[iplane] = get_slice_data(dir, point, mf_cc_data, geom[ilev],
868 0, ncomp, interpolate, bnd_rbx);
871 auto min_bnd_bx =
m_ps_mf[iplane]->boxArray().minimalBox();
872 amrex::Box bnd_bx =
getIndexBox(bnd_rbx, geom[ilev]);
873 if (bnd_bx == min_bnd_bx) {
break; }
881 amrex::Vector<int>& level_steps,
882 amrex::Vector<amrex::IntVect>& ref_ratio,
883 amrex::Vector<amrex::Geometry>& geom)
886 for (
int iplane(0); iplane<nplane; ++iplane) {
888 int dir =
m_dir[iplane];
889 int lev =
m_lev[iplane];
891 amrex::Vector<int> m_level_steps = {level_steps[lev]};
892 amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
898 amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
899 for (
int d(0); d<AMREX_SPACEDIM; ++d) {
901 m_rb.setLo(d,point-
myhalf*geom[lev].CellSize(d));
902 m_rb.setHi(d,point+
myhalf*geom[lev].CellSize(d));
904 is_per[d] = geom[lev].isPeriodic(d);
906 amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
907 m_geom[0].define(m_dom, &m_rb, geom[lev].
Coord(), is_per.data());
910 std::string name_plane =
m_name[iplane];
911 name_plane +=
"_step_";
912 std::string plotfilename = amrex::Concatenate(name_plane, m_level_steps[0], 5);
915 amrex::Vector<const amrex::MultiFab*> mf = {
m_ps_mf[iplane].get()};
918 WriteMultiLevelPlotfile(plotfilename, 1, mf,
920 m_level_steps, m_ref_ratio);
927 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(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
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:122
@ xvel
Definition: ERF_IndexDefines.H:175
@ cons
Definition: ERF_IndexDefines.H:174
@ zvel
Definition: ERF_IndexDefines.H:177
@ yvel
Definition: ERF_IndexDefines.H:176
Definition: ERF_SampleData.H:18
amrex::Vector< std::string > m_varnames
Definition: ERF_SampleData.H:564
amrex::Vector< amrex::MultiFab > m_ls_mf
Definition: ERF_SampleData.H:559
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:565
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:507
bool m_use_real_bx
Definition: ERF_SampleData.H:562
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:555
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:560
void write_line_ascii(amrex::Vector< amrex::Real > &time)
Definition: ERF_SampleData.H:473
LineSampler()
Definition: ERF_SampleData.H:19
amrex::Vector< amrex::Box > m_bnd_bx
Definition: ERF_SampleData.H:557
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:460
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:563
amrex::Vector< amrex::RealBox > m_bnd_rbx
Definition: ERF_SampleData.H:558
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:556
Definition: ERF_SampleData.H:570
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:880
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:924
amrex::Box getIndexBox(const amrex::RealBox &real_box, const amrex::Geometry &geom)
Definition: ERF_SampleData.H:685
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:925
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:928
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_ps_mf
Definition: ERF_SampleData.H:927
amrex::Vector< amrex::RealBox > m_bnd_rbx
Definition: ERF_SampleData.H:926
PlaneSampler()
Definition: ERF_SampleData.H:571
amrex::Vector< std::string > m_varnames
Definition: ERF_SampleData.H:930
void get_sample_data(amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
Definition: ERF_SampleData.H:701