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

#include <ERF_SampleData.H>

Collaboration diagram for LineSampler:

Public Member Functions

 LineSampler ()
 
amrex::Box getIndexBox (const amrex::RealBox &real_box, const amrex::Geometry &geom)
 
void write_coords (amrex::Vector< std::unique_ptr< amrex::MultiFab > > &z_phys_cc, amrex::Vector< amrex::Geometry > &geom)
 
void get_sample_data (amrex::Vector< amrex::Geometry > &geom, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
 
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)
 
void write_line_ascii (amrex::Vector< amrex::Real > &time)
 
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)
 

Public Attributes

amrex::Vector< int > m_dir
 
amrex::Vector< int > m_lev
 
amrex::Vector< amrex::Box > m_bnd_bx
 
amrex::Vector< amrex::RealBox > m_bnd_rbx
 
amrex::Vector< amrex::MultiFab > m_ls_mf
 
amrex::Vector< std::string > m_name
 
bool m_use_real_bx {false}
 
bool m_write_ascii {false}
 
amrex::Vector< std::string > m_varnames {"magvel","theta"}
 
amrex::Vector< std::unique_ptr< std::fstream > > m_datastream
 

Constructor & Destructor Documentation

◆ LineSampler()

LineSampler::LineSampler ( )
inline
20  {
21  amrex::ParmParse pp("erf");
22 
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");
27  }
28 
29  // Count number of lo and hi points define the line
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");
35  if (has_line_idx) {
36  AMREX_ALWAYS_ASSERT( (n_line_lo==n_line_hi ) &&
37  (n_line_lo==n_line_dir) );
38  } else if (has_line_real) {
39  AMREX_ALWAYS_ASSERT( (n_line_lo_real==n_line_hi_real ) &&
40  (n_line_lo_real==n_line_dir) );
41  }
42 
43  // Parse the data
44  int nline = 0;
45  if (has_line_idx) {
46  nline = n_line_lo;
47  if (nline > 0) {
48  // Parse lo
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]);
56  iv_lo[i] = iv;
57  }
58 
59  // Parse hi
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]);
67  iv_hi[i] = iv;
68  }
69 
70  // Construct vector of bounding boxes
71  m_bnd_bx.resize(nline);
72  for (int i = 0; i < nline; i++){
73  amrex::Box lbx(iv_lo[i],iv_hi[i]);
74  m_bnd_bx[i] = lbx;
75  }
76  }
77  } else if (has_line_real) {
78  m_use_real_bx = true;
79  nline = n_line_lo_real;
80  if (nline > 0) {
81  // Parse lo
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]};
89  rv_lo.push_back(rv);
90  }
91 
92  // Parse hi
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]};
100  rv_hi.push_back(rv);
101  }
102 
103  // Construct vector of bounding real boxes
104  m_bnd_rbx.resize(nline);
105  for (int i = 0; i < nline; i++){
106  amrex::RealBox rbx(rv_lo[i].data(),rv_hi[i].data());
107  m_bnd_rbx[i] = rbx;
108  }
109 
110  // Construct vector of bounding boxes
111  m_bnd_bx.resize(nline);
112  }
113  }
114 
115  if (nline > 0) {
116  // Parse directionality
117  m_dir.resize(n_line_dir);
118  pp.queryarr("sample_line_dir",m_dir,0,n_line_dir);
119 
120  // Parse names
121  std::string name_base = "plt_line_";
122  m_name.resize(nline);
123  int n_names = pp.countval("sample_line_name");
124  if (n_names > 0) {
125  AMREX_ALWAYS_ASSERT( n_names==nline );
126  pp.queryarr("sample_line_name",m_name,0,n_names);
127  } else {
128  for (int iline(0); iline<nline; ++iline) {
129  m_name[iline] = amrex::Concatenate(name_base, iline , 5);
130  }
131  }
132 
133  // Allocate space for level indicator
134  m_lev.resize(n_line_dir,0);
135 
136  // Allocate space for MF pointers
137  m_ls_mf.resize(nline);
138 
139  // Get requested vars
140  if (pp.countval("line_sampling_vars") > 0) {
141  m_varnames.clear();
142  amrex::Vector<std::string> requested_vars;
143  pp.queryarr("line_sampling_vars",requested_vars);
144  amrex::Print() << "Selected line sampling vars :";
145  if (containerHasElement(requested_vars, "density")) {
146  m_varnames.push_back("density");
147  amrex::Print() << " " << "density";
148  }
149  if (containerHasElement(requested_vars, "x_velocity")) {
150  m_varnames.push_back("x_velocity");
151  amrex::Print() << " " << "x_velocity";
152  }
153  if (containerHasElement(requested_vars, "y_velocity")) {
154  m_varnames.push_back("y_velocity");
155  amrex::Print() << " " << "y_velocity";
156  }
157  if (containerHasElement(requested_vars, "z_velocity")) {
158  m_varnames.push_back("z_velocity");
159  amrex::Print() << " " << "z_velocity";
160  }
161  if (containerHasElement(requested_vars, "magvel")) {
162  m_varnames.push_back("magvel");
163  amrex::Print() << " " << "magvel";
164  }
165  if (containerHasElement(requested_vars, "theta")) {
166  m_varnames.push_back("theta");
167  amrex::Print() << " " << "theta";
168  }
169  if (containerHasElement(requested_vars, "qv")) {
170  m_varnames.push_back("qv");
171  amrex::Print() << " " << "qv";
172  }
173  if (containerHasElement(requested_vars, "qc")) {
174  m_varnames.push_back("qc");
175  amrex::Print() << " " << "qc";
176  }
177  if (containerHasElement(requested_vars, "pressure")) {
178  m_varnames.push_back("pressure");
179  amrex::Print() << " " << "pressure";
180  }
181  amrex::Print() << std::endl;
182  }
183 
184  // Write outputs to text files, one file per variable, with all
185  // times appended to the same file
186  pp.query("line_sampling_text_output",m_write_ascii);
187  if (m_write_ascii && amrex::ParallelDescriptor::IOProcessor()) {
188  int nvar = static_cast<int>(m_varnames.size());
189  m_datastream.resize(nline * nvar);
190  int i = 0;
191  for (int iline(0); iline<nline; ++iline) {
192  for (int ivar(0); ivar<nvar; ++ivar) {
193  std::string filename = m_name[iline] + "." + m_varnames[ivar];
194  m_datastream[i] = std::make_unique<std::fstream>();
195  m_datastream[i]->open(filename.c_str(),std::ios::out|std::ios::app);
196  if (!m_datastream[i]->good()) {
197  amrex::FileOpenFailed(filename);
198  }
199  i++;
200  }
201  }
202  }
203  }
204  }
bool containerHasElement(const V &iterable, const T &query)
Definition: ERF_Container.H:5
struct @22 out
ParmParse pp("prob")
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
amrex::Vector< std::string > m_varnames
Definition: ERF_SampleData.H:560
amrex::Vector< amrex::MultiFab > m_ls_mf
Definition: ERF_SampleData.H:555
amrex::Vector< std::unique_ptr< std::fstream > > m_datastream
Definition: ERF_SampleData.H:561
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
amrex::Vector< amrex::Box > m_bnd_bx
Definition: ERF_SampleData.H:553
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
Here is the call graph for this function:

Member Function Documentation

◆ get_sample_data()

void LineSampler::get_sample_data ( amrex::Vector< amrex::Geometry > &  geom,
amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars_new 
)
inline
278  {
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());
282 
283  int qv_comp = -1;
284 
285  // Loop over each line
286  for (int iline(0); iline<nline; ++iline) {
287  int dir = m_dir[iline];
288 
289  // Search each level to get the finest data possible
290  for (int ilev(nlev-1); ilev>=0; --ilev) {
291  amrex::Box bnd_bx = (m_use_real_bx)
292  ? getIndexBox(m_bnd_rbx[iline], geom[ilev])
293  : m_bnd_bx[iline];
294  amrex::IntVect cell = bnd_bx.smallEnd();
295 
296  // Construct CC velocities
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],
303  &vars_new[ilev][Vars::yvel],
304  &vars_new[ilev][Vars::zvel]});
305 
306  // Construct vector of MFs holding T and WSP
307  amrex::MultiFab mf_cc_data;
308  mf_cc_data.define(ba, dm, ncomp, 1);
309 
310  int mf_comp = 0;
311 
312  if (containerHasElement(m_varnames, "density")) {
313  amrex::MultiFab::Copy(mf_cc_data, vars_new[ilev][Vars::cons], Rho_comp, mf_comp, 1, 0);
314  mf_comp += 1;
315  }
316 
317  if (containerHasElement(m_varnames, "x_velocity")) {
318  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 0, mf_comp, 1, 0);
319  mf_comp += 1;
320  }
321  if (containerHasElement(m_varnames, "y_velocity")) {
322  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 1, mf_comp, 1, 0);
323  mf_comp += 1;
324  }
325  if (containerHasElement(m_varnames, "z_velocity")) {
326  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 2, mf_comp, 1, 0);
327  mf_comp += 1;
328  }
329 
330  if (containerHasElement(m_varnames, "magvel")) {
331 #ifdef _OPENMP
332 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
333 #endif
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);
338 
339  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
340  {
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)) ;
344  });
345  }
346  mf_comp += 1;
347  }
348 
349  if (containerHasElement(m_varnames, "theta")) {
350 #ifdef _OPENMP
351 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
352 #endif
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);
357 
358  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
359  {
360  dfab(i,j,k,mf_comp) = cfab(i,j,k,RhoTheta_comp) / cfab(i,j,k,Rho_comp);
361  });
362  }
363  mf_comp += 1;
364  }
365 
366  if (containerHasElement(m_varnames, "qv")) {
367  // if qv is requested, assume that we have moisture
368  if (qv_comp >= 0) AMREX_ALWAYS_ASSERT(qv_comp == mf_comp);
369  qv_comp = mf_comp;
370 #ifdef _OPENMP
371 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
372 #endif
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);
377 
378  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
379  {
380  dfab(i,j,k,mf_comp) = cfab(i,j,k,RhoQ1_comp) / cfab(i,j,k,Rho_comp);
381  });
382  }
383  mf_comp += 1;
384  }
385  if (containerHasElement(m_varnames, "qc")) {
386 #ifdef _OPENMP
387 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
388 #endif
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);
393 
394  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
395  {
396  dfab(i,j,k,mf_comp) = cfab(i,j,k,RhoQ2_comp) / cfab(i,j,k,Rho_comp);
397  });
398  }
399  mf_comp += 1;
400  }
401 
402  if (containerHasElement(m_varnames, "pressure") && (qv_comp >= 0)) {
403  // if qv is requested, assume that we have moisture
404 #ifdef _OPENMP
405 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
406 #endif
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);
411 
412  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
413  {
414  dfab(i,j,k,mf_comp) = getPgivenRTh(cfab(i,j,k,RhoTheta_comp),
415  dfab(i,j,k,qv_comp));
416  });
417  }
418  mf_comp += 1;
419  }
420  else if (containerHasElement(m_varnames, "pressure")) {
421  // pressure requested w/o qv, assume dry
422 #ifdef _OPENMP
423 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
424 #endif
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);
429 
430  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
431  {
432  dfab(i,j,k,mf_comp) = getPgivenRTh(cfab(i,j,k,RhoTheta_comp));
433  });
434  }
435  mf_comp += 1;
436  }
437 
438  m_lev[iline] = ilev;
439  m_ls_mf[iline] = get_line_data(mf_cc_data, dir, cell, bnd_bx);
440 
441  // We can stop if we got the entire line
442  auto min_bnd_bx = m_ls_mf[iline].boxArray().minimalBox();
443  if (bnd_bx == min_bnd_bx) { break; }
444 
445  } // ilev
446 
447  // Fill bnd_bx if use_rbx
448  if (m_use_real_bx) {
449  m_bnd_bx[iline] = m_ls_mf[iline].boxArray().minimalBox();
450  }
451 
452  }// iline
453  }
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
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));})
@ xvel
Definition: ERF_IndexDefines.H:159
@ cons
Definition: ERF_IndexDefines.H:158
@ zvel
Definition: ERF_IndexDefines.H:161
@ yvel
Definition: ERF_IndexDefines.H:160
amrex::Box getIndexBox(const amrex::RealBox &real_box, const amrex::Geometry &geom)
Definition: ERF_SampleData.H:207
Here is the call graph for this function:

◆ getIndexBox()

amrex::Box LineSampler::getIndexBox ( const amrex::RealBox &  real_box,
const amrex::Geometry &  geom 
)
inline
208  {
209  amrex::IntVect slice_lo, slice_hi;
210 
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))););
214 
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))););
218 
219  return amrex::Box(slice_lo, slice_hi) & geom.Domain();
220  }

Referenced by get_sample_data(), and write_coords().

Here is the caller graph for this function:

◆ write_coords()

void LineSampler::write_coords ( amrex::Vector< std::unique_ptr< amrex::MultiFab > > &  z_phys_cc,
amrex::Vector< amrex::Geometry > &  geom 
)
inline
225  {
226  if (!m_write_ascii) return;
227 
228  amrex::Print() << "Writing out line coordinates to text" << std::endl;
229 
230  for (int lev(0); lev < z_phys_cc.size(); ++lev) {
231  // Write one text file per level
232  std::ofstream outfile;
233  if (amrex::ParallelDescriptor::IOProcessor()) {
234  std::string fname = amrex::Concatenate("plt_line_lev", lev, 1);
235  fname += ".zcc";
236  outfile.open(fname);
237 
238  if (!outfile.is_open()) {
239  amrex::AllPrint() << "Could not open " << fname << std::endl;
240  }
241  }
242 
243  // Loop over each line
244  int nline = static_cast<int>(m_ls_mf.size());
245  for (int iline(0); iline<nline; ++iline) {
246  int dir = m_dir[iline];
247  amrex::Box bnd_bx = (m_use_real_bx)
248  ? getIndexBox(m_bnd_rbx[iline], geom[lev])
249  : m_bnd_bx[iline];
250  amrex::IntVect first_cell = bnd_bx.smallEnd();
251 
252  // Create multifab with "sampled" z_phys values
253  amrex::MultiFab line_coords_mf = get_line_data(
254  *z_phys_cc[lev], dir, first_cell, bnd_bx
255  );
256 
257  // Convert multifab to vector
258  amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(
259  line_coords_mf, 0, 1, bnd_bx, dir
260  );
261 
262  // Append to file
263  if (amrex::ParallelDescriptor::IOProcessor()) {
264  for (const auto& zval : vec) {
265  outfile << " " << zval;
266  }
267  outfile << std::endl;
268  }
269  } // line loop
270 
271  outfile.close();
272  } // level loop
273  }
Here is the call graph for this function:

◆ write_line_ascii()

void LineSampler::write_line_ascii ( amrex::Vector< amrex::Real > &  time)
inline
470  {
471  // same as definitions in ERF.H
472  constexpr int datwidth = 14;
473  constexpr int datprecision = 6;
474  constexpr int timeprecision = 13;
475 
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];
481  amrex::Real m_time = time[lev];
482  amrex::Box m_dom = m_bnd_bx[iline];
483 
484  for (int ivar(0); ivar<nvar; ++ivar) {
485  // Convert multifab to vector
486  amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(m_ls_mf[iline], ivar, 1, m_dom, dir);
487 
488  if (amrex::ParallelDescriptor::IOProcessor()) {
489  int ifile = iline*nvar + ivar;
490  std::ostream& fs = *m_datastream[ifile];
491  fs << std::setw(datwidth) << std::setprecision(timeprecision) << m_time
492  << std::setw(datwidth) << std::setprecision(datprecision);
493  for (const auto& val : vec) {
494  fs << " " << val;
495  }
496  fs << std::endl;
497  }
498  }
499  }
500  }
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ fs
Definition: ERF_AdvanceMorrison.cpp:119

Referenced by write_sample_data().

Here is the caller graph for this function:

◆ write_line_plotfile()

void LineSampler::write_line_plotfile ( amrex::Vector< amrex::Real > &  time,
amrex::Vector< int > &  level_steps,
amrex::Vector< amrex::IntVect > &  ref_ratio,
amrex::Vector< amrex::Geometry > &  geom 
)
inline
507  {
508  int nline = static_cast<int>(m_ls_mf.size());
509  for (int iline(0); iline<nline; ++iline) {
510  // Data members that can be used as-is
511  int dir = m_dir[iline];
512  int lev = m_lev[iline];
513  amrex::Real m_time = time[lev];
514  amrex::Vector<int> m_level_steps = {level_steps[lev]};
515  amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
516 
517  // Create modified geometry object corresponding to the line
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);
522  amrex::Box m_dom = m_bnd_bx[iline];
523  amrex::RealBox m_rb;
524  for (int d(0); d<AMREX_SPACEDIM; ++d) {
525  amrex::Real offset = (d==dir) ? 0 : myhalf;
526  amrex::Real lo = plo[d] + ( m_dom.smallEnd(d) - offset ) * dx[d];
527  amrex::Real hi = plo[d] + ( m_dom.bigEnd(d) + offset ) * dx[d];
528 
529  m_rb.setLo(d,lo);
530  m_rb.setHi(d,hi);
531 
532  is_per[d] = geom[lev].isPeriodic(d);
533  }
534  m_geom[0].define(m_dom, &m_rb, geom[lev].Coord(), is_per.data());
535 
536  // Create plotfile name
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);
540 
541  // Get the data
542  amrex::Vector<const amrex::MultiFab*> mf = {&(m_ls_mf[iline])};
543 
544  // Write each line
545  WriteMultiLevelPlotfile(plotfilename, 1, mf,
546  m_varnames, m_geom, m_time,
547  m_level_steps, m_ref_ratio);
548  }
549  }
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
Coord
Definition: ERF_DataStruct.H:92
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28

Referenced by write_sample_data().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ write_sample_data()

void LineSampler::write_sample_data ( amrex::Vector< amrex::Real > &  time,
amrex::Vector< int > &  level_steps,
amrex::Vector< amrex::IntVect > &  ref_ratio,
amrex::Vector< amrex::Geometry > &  geom 
)
inline
460  {
461  if (m_write_ascii) {
462  write_line_ascii(time);
463  } else {
464  write_line_plotfile(time, level_steps, ref_ratio, geom);
465  }
466  }
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
void write_line_ascii(amrex::Vector< amrex::Real > &time)
Definition: ERF_SampleData.H:469
Here is the call graph for this function:

Member Data Documentation

◆ m_bnd_bx

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

◆ m_bnd_rbx

amrex::Vector<amrex::RealBox> LineSampler::m_bnd_rbx

◆ m_datastream

amrex::Vector<std::unique_ptr<std::fstream> > LineSampler::m_datastream

Referenced by LineSampler(), and write_line_ascii().

◆ 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

◆ m_name

amrex::Vector<std::string> LineSampler::m_name

Referenced by LineSampler(), and write_line_plotfile().

◆ m_use_real_bx

bool LineSampler::m_use_real_bx {false}

◆ m_varnames

amrex::Vector<std::string> LineSampler::m_varnames {"magvel","theta"}

◆ m_write_ascii

bool LineSampler::m_write_ascii {false}

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