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() <<
" " <<
"x_velocity";
96 amrex::Print() <<
" " <<
"y_velocity";
100 amrex::Print() <<
" " <<
"z_velocity";
104 amrex::Print() <<
" " <<
"magvel";
108 amrex::Print() <<
" " <<
"theta";
112 amrex::Print() <<
" " <<
"qv";
116 amrex::Print() <<
" " <<
"pressure";
118 amrex::Print() << std::endl;
124 if (
m_write_ascii && amrex::ParallelDescriptor::IOProcessor()) {
128 for (
int iline(0); iline<n_line_lo; ++iline) {
129 for (
int ivar(0); ivar<nvar; ++ivar) {
132 m_datastream[i]->open(filename.c_str(),std::ios::out|std::ios::app);
134 amrex::FileOpenFailed(filename);
146 int nlev = vars_new.size();
151 for (
int iline(0); iline<nline; ++iline) {
152 int dir =
m_dir[iline];
153 amrex::Box bnd_bx =
m_bnd_bx[iline];
154 amrex::IntVect cell = bnd_bx.smallEnd();
157 for (
int ilev(nlev-1); ilev>=0; --ilev) {
160 amrex::MultiFab mf_cc_vel;
161 auto ba = vars_new[ilev][
Vars::cons].boxArray();
162 auto dm = vars_new[ilev][
Vars::cons].DistributionMap();
163 mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
164 average_face_to_cellcenter(mf_cc_vel,0,
165 amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][
Vars::xvel],
170 amrex::MultiFab mf_cc_data;
171 mf_cc_data.define(ba, dm, ncomp, 1);
176 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 0, mf_comp, 1, 0);
180 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 1, mf_comp, 1, 0);
184 amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 2, mf_comp, 1, 0);
190 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
192 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
193 const amrex::Box& tbx = mfi.tilebox();
194 auto const& dfab = mf_cc_data.array(mfi);
195 auto const& vfab = mf_cc_vel.array(mfi);
197 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
199 dfab(i,j,k,mf_comp) = std::sqrt(vfab(i,j,k,0)*vfab(i,j,k,0)
200 + vfab(i,j,k,1)*vfab(i,j,k,1)
201 + vfab(i,j,k,2)*vfab(i,j,k,2)) ;
209 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
211 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
212 const amrex::Box& tbx = mfi.tilebox();
213 auto const& dfab = mf_cc_data.array(mfi);
214 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
216 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
227 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
229 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
230 const amrex::Box& tbx = mfi.tilebox();
231 auto const& dfab = mf_cc_data.array(mfi);
232 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
234 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
238 dfab(i,j,k,mf_comp));
245 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
247 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
248 const amrex::Box& tbx = mfi.tilebox();
249 auto const& dfab = mf_cc_data.array(mfi);
250 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
252 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
261 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
263 for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
264 const amrex::Box& tbx = mfi.tilebox();
265 auto const& dfab = mf_cc_data.array(mfi);
266 auto const& cfab = vars_new[ilev][
Vars::cons].array(mfi);
268 amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
277 m_ls_mf[iline] = get_line_data(mf_cc_data, dir, cell, bnd_bx);
280 auto min_bnd_bx =
m_ls_mf[iline].boxArray().minimalBox();
281 if (bnd_bx == min_bnd_bx) {
continue; }
289 amrex::Vector<int>& level_steps,
290 amrex::Vector<amrex::IntVect>& ref_ratio,
291 amrex::Vector<amrex::Geometry>& geom)
304 constexpr
int datwidth = 14;
305 constexpr
int datprecision = 6;
306 constexpr
int timeprecision = 13;
310 for (
int iline(0); iline<nline; ++iline) {
311 int dir =
m_dir[iline];
312 int lev =
m_lev[iline];
313 amrex::Real m_time = time[lev];
316 for (
int ivar(0); ivar<nvar; ++ivar) {
318 amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(
m_ls_mf[iline], ivar, 1, m_dom, dir);
320 if (amrex::ParallelDescriptor::IOProcessor()) {
321 int ifile = iline*nvar + ivar;
323 fs << std::setw(datwidth) << std::setprecision(timeprecision) << m_time
324 << std::setw(datwidth) << std::setprecision(datprecision);
325 for (
const auto& val : vec) {
336 amrex::Vector<int>& level_steps,
337 amrex::Vector<amrex::IntVect>& ref_ratio,
338 amrex::Vector<amrex::Geometry>& geom)
341 for (
int iline(0); iline<nline; ++iline) {
343 int dir =
m_dir[iline];
344 int lev =
m_lev[iline];
345 amrex::Real m_time = time[lev];
346 amrex::Vector<int> m_level_steps = {level_steps[lev]};
347 amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
350 auto plo = geom[lev].ProbLo();
351 auto dx = geom[lev].CellSize();
352 amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
353 amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
356 for (
int d(0); d<AMREX_SPACEDIM; ++d) {
357 amrex::Real
offset = (d==dir) ? 0 : 0.5;
358 amrex::Real lo = plo[d] + ( m_dom.smallEnd(d) -
offset ) * dx[d];
359 amrex::Real hi = plo[d] + ( m_dom.bigEnd(d) +
offset ) * dx[d];
364 is_per[d] = geom[lev].isPeriodic(d);
366 m_geom[0].define(m_dom, &m_rb, geom[lev].
Coord(), is_per.data());
369 std::string name_line =
m_name[iline];
370 name_line +=
"_step_";
371 std::string plotfilename = amrex::Concatenate(name_line, m_level_steps[0], 5);
374 amrex::Vector<const amrex::MultiFab*> mf = {&(
m_ls_mf[iline])};
377 WriteMultiLevelPlotfile(plotfilename, 1, mf,
379 m_level_steps, m_ref_ratio);
399 amrex::ParmParse
pp(
"erf");
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) );
409 if (n_plane_lo > 0) {
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]};
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]};
434 for (
int i(0); i < n_plane_hi; i++){
435 amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
440 m_dir.resize(n_plane_dir);
441 pp.queryarr(
"sample_plane_dir",
m_dir,0,n_plane_dir);
444 std::string name_base =
"plt_plane_";
445 m_name.resize(n_plane_lo);
446 int n_names =
pp.countval(
"sample_plane_name");
448 AMREX_ALWAYS_ASSERT( n_names==n_plane_lo );
449 pp.queryarr(
"sample_plane_name",
m_name,0,n_names);
451 for (
int iplane(0); iplane<n_plane_lo; ++iplane) {
452 m_name[iplane] = amrex::Concatenate(name_base, iplane , 5);
457 m_lev.resize(n_plane_dir,0);
467 const amrex::Geometry& geom) {
468 amrex::IntVect slice_lo, slice_hi;
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))););
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))););
478 return amrex::Box(slice_lo, slice_hi) & geom.Domain();
483 amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_new)
485 int nlev = vars_new.size();
488 bool interpolate =
true;
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);
497 for (
int ilev(nlev-1); ilev>=0; --ilev) {
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],
510 amrex::MultiFab mf_cc_data;
511 mf_cc_data.define(ba, dm, ncomp, 1);
513 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
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
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)) ;
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);
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; }
545 amrex::Vector<int>& level_steps,
546 amrex::Vector<amrex::IntVect>& ref_ratio,
547 amrex::Vector<amrex::Geometry>& geom)
549 amrex::Vector<std::string> varnames = {
"T",
"Wsp"};
552 for (
int iplane(0); iplane<nplane; ++iplane) {
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]};
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) {
567 m_rb.setLo(d,point-0.5*geom[lev].CellSize(d));
568 m_rb.setHi(d,point+0.5*geom[lev].CellSize(d));
570 is_per[d] = geom[lev].isPeriodic(d);
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());
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);
581 amrex::Vector<const amrex::MultiFab*> mf = {
m_ps_mf[iplane].get()};
584 WriteMultiLevelPlotfile(plotfilename, 1, mf,
585 varnames, m_geom, m_time,
586 m_level_steps, m_ref_ratio);
593 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
m_ps_mf;
604 if(do_line)
m_ls = std::make_unique<LineSampler >();
605 if(do_plane)
m_ps = std::make_unique<PlaneSampler>();
610 amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_new)
612 if (
m_ls)
m_ls->get_line_mfs(vars_new);
613 if (
m_ps)
m_ps->get_plane_mfs(geom, vars_new);
618 amrex::Vector<int>& level_steps,
619 amrex::Vector<amrex::IntVect>& ref_ratio,
620 amrex::Vector<amrex::Geometry>& geom)
622 if (
m_ls)
m_ls->write_line_mfs(time, level_steps, ref_ratio, geom);
623 if (
m_ps)
m_ps->write_plane_mfs(time, level_steps, ref_ratio, geom);
629 std::unique_ptr<LineSampler>
m_ls =
nullptr;
632 std::unique_ptr<PlaneSampler>
m_ps =
nullptr;
bool containerHasElement(const V &iterable, const T &query)
Definition: ERF_Container.H:5
Coord
Definition: ERF_DataStruct.H:68
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:84
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#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:219
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
Definition: ERF_SampleData.H:599
std::unique_ptr< LineSampler > m_ls
Definition: ERF_SampleData.H:629
void get_sample_data(amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
Definition: ERF_SampleData.H:609
SampleData(bool do_line=false, bool do_plane=false)
Definition: ERF_SampleData.H:601
std::unique_ptr< PlaneSampler > m_ps
Definition: ERF_SampleData.H:632
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:617
@ 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:390
amrex::Vector< amrex::MultiFab > m_ls_mf
Definition: ERF_SampleData.H:386
void get_line_mfs(amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
Definition: ERF_SampleData.H:144
amrex::Vector< std::unique_ptr< std::fstream > > m_datastream
Definition: ERF_SampleData.H:391
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:335
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:383
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:387
void write_line_ascii(amrex::Vector< amrex::Real > &time)
Definition: ERF_SampleData.H:301
LineSampler()
Definition: ERF_SampleData.H:19
amrex::Vector< amrex::Box > m_bnd_bx
Definition: ERF_SampleData.H:385
bool m_write_ascii
Definition: ERF_SampleData.H:389
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)
Definition: ERF_SampleData.H:288
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:384
Definition: ERF_SampleData.H:396
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:590
amrex::Box getIndexBox(const amrex::RealBox &real_box, const amrex::Geometry &geom)
Definition: ERF_SampleData.H:466
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:591
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:594
void get_plane_mfs(amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
Definition: ERF_SampleData.H:482
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
PlaneSampler()
Definition: ERF_SampleData.H:397
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)
Definition: ERF_SampleData.H:544