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

#include <ERF_SampleData.H>

Collaboration diagram for LineSampler:

Public Member Functions

 LineSampler ()
 
void write_line_coords (amrex::Vector< std::unique_ptr< amrex::MultiFab > > &z_phys_cc)
 
void get_line_mfs (amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new)
 
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)
 
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::MultiFab > m_ls_mf
 
amrex::Vector< std::string > m_name
 
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  // Count number of lo and hi points define the line
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) );
29 
30  // Parse the data
31  if (n_line_lo > 0) {
32  // Parse lo
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]);
40  iv_lo[i] = iv;
41  }
42 
43  // Parse hi
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]);
51  iv_hi[i] = iv;
52  }
53 
54  // Construct vector of bounding boxes
55  m_bnd_bx.resize(n_line_lo);
56  for (int i = 0; i < n_line_hi; i++){
57  amrex::Box lbx(iv_lo[i],iv_hi[i]);
58  m_bnd_bx[i] = lbx;
59  }
60 
61  // Parse directionality
62  m_dir.resize(n_line_dir);
63  pp.queryarr("sample_line_dir",m_dir,0,n_line_dir);
64 
65  // Parse names
66  std::string name_base = "plt_line_";
67  m_name.resize(n_line_lo);
68  int n_names = pp.countval("sample_line_name");
69  if (n_names > 0) {
70  AMREX_ALWAYS_ASSERT( n_names==n_line_lo );
71  pp.queryarr("sample_line_name",m_name,0,n_names);
72  } else {
73  for (int iline(0); iline<n_line_lo; ++iline) {
74  m_name[iline] = amrex::Concatenate(name_base, iline , 5);
75  }
76  }
77 
78  // Allocate space for level indicator
79  m_lev.resize(n_line_dir,0);
80 
81  // Allocate space for MF pointers
82  m_ls_mf.resize(n_line_lo);
83 
84  // Get requested vars
85  if (pp.countval("line_sampling_vars") > 0) {
86  m_varnames.clear();
87  amrex::Vector<std::string> requested_vars;
88  pp.queryarr("line_sampling_vars",requested_vars);
89  amrex::Print() << "Selected line sampling vars :";
90  if (containerHasElement(requested_vars, "density")) {
91  m_varnames.push_back("density");
92  amrex::Print() << " " << "density";
93  }
94  if (containerHasElement(requested_vars, "x_velocity")) {
95  m_varnames.push_back("x_velocity");
96  amrex::Print() << " " << "x_velocity";
97  }
98  if (containerHasElement(requested_vars, "y_velocity")) {
99  m_varnames.push_back("y_velocity");
100  amrex::Print() << " " << "y_velocity";
101  }
102  if (containerHasElement(requested_vars, "z_velocity")) {
103  m_varnames.push_back("z_velocity");
104  amrex::Print() << " " << "z_velocity";
105  }
106  if (containerHasElement(requested_vars, "magvel")) {
107  m_varnames.push_back("magvel");
108  amrex::Print() << " " << "magvel";
109  }
110  if (containerHasElement(requested_vars, "theta")) {
111  m_varnames.push_back("theta");
112  amrex::Print() << " " << "theta";
113  }
114  if (containerHasElement(requested_vars, "qv")) {
115  m_varnames.push_back("qv");
116  amrex::Print() << " " << "qv";
117  }
118  if (containerHasElement(requested_vars, "qc")) {
119  m_varnames.push_back("qc");
120  amrex::Print() << " " << "qc";
121  }
122  if (containerHasElement(requested_vars, "pressure")) {
123  m_varnames.push_back("pressure");
124  amrex::Print() << " " << "pressure";
125  }
126  amrex::Print() << std::endl;
127  }
128 
129  // Write outputs to text files, one file per variable, with all
130  // times appended to the same file
131  pp.query("line_sampling_text_output",m_write_ascii);
132  if (m_write_ascii && amrex::ParallelDescriptor::IOProcessor()) {
133  int nvar = static_cast<int>(m_varnames.size());
134  m_datastream.resize(n_line_lo * nvar);
135  int i = 0;
136  for (int iline(0); iline<n_line_lo; ++iline) {
137  for (int ivar(0); ivar<nvar; ++ivar) {
138  std::string filename = m_name[iline] + "." + m_varnames[ivar];
139  m_datastream[i] = std::make_unique<std::fstream>();
140  m_datastream[i]->open(filename.c_str(),std::ios::out|std::ios::app);
141  if (!m_datastream[i]->good()) {
142  amrex::FileOpenFailed(filename);
143  }
144  i++;
145  }
146  }
147  }
148  }
149  }
bool containerHasElement(const V &iterable, const T &query)
Definition: ERF_Container.H:5
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:230
amrex::Vector< std::string > m_varnames
Definition: ERF_SampleData.H:475
amrex::Vector< amrex::MultiFab > m_ls_mf
Definition: ERF_SampleData.H:471
amrex::Vector< std::unique_ptr< std::fstream > > m_datastream
Definition: ERF_SampleData.H:476
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:468
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:472
amrex::Vector< amrex::Box > m_bnd_bx
Definition: ERF_SampleData.H:470
bool m_write_ascii
Definition: ERF_SampleData.H:474
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:469
Here is the call graph for this function:

Member Function Documentation

◆ get_line_mfs()

void LineSampler::get_line_mfs ( amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars_new)
inline
203  {
204  int nlev = static_cast<int>(vars_new.size());
205  int nline = static_cast<int>(m_bnd_bx.size());
206  int ncomp = static_cast<int>(m_varnames.size());
207 
208  int qv_comp = -1;
209 
210  // Loop over each line
211  for (int iline(0); iline<nline; ++iline) {
212  int dir = m_dir[iline];
213  amrex::Box bnd_bx = m_bnd_bx[iline];
214  amrex::IntVect cell = bnd_bx.smallEnd();
215 
216  // Search each level to get the finest data possible
217  for (int ilev(nlev-1); ilev>=0; --ilev) {
218 
219  // Construct CC velocities
220  amrex::MultiFab mf_cc_vel;
221  auto ba = vars_new[ilev][Vars::cons].boxArray();
222  auto dm = vars_new[ilev][Vars::cons].DistributionMap();
223  mf_cc_vel.define(ba, dm, AMREX_SPACEDIM, amrex::IntVect(1,1,1));
224  average_face_to_cellcenter(mf_cc_vel,0,
225  amrex::Array<const amrex::MultiFab*,3>{&vars_new[ilev][Vars::xvel],
226  &vars_new[ilev][Vars::yvel],
227  &vars_new[ilev][Vars::zvel]});
228 
229  // Construct vector of MFs holding T and WSP
230  amrex::MultiFab mf_cc_data;
231  mf_cc_data.define(ba, dm, ncomp, 1);
232 
233  int mf_comp = 0;
234 
235  if (containerHasElement(m_varnames, "density")) {
236  amrex::MultiFab::Copy(mf_cc_data, vars_new[ilev][Vars::cons], Rho_comp, mf_comp, 1, 0);
237  mf_comp += 1;
238  }
239 
240  if (containerHasElement(m_varnames, "x_velocity")) {
241  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 0, mf_comp, 1, 0);
242  mf_comp += 1;
243  }
244  if (containerHasElement(m_varnames, "y_velocity")) {
245  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 1, mf_comp, 1, 0);
246  mf_comp += 1;
247  }
248  if (containerHasElement(m_varnames, "z_velocity")) {
249  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 2, mf_comp, 1, 0);
250  mf_comp += 1;
251  }
252 
253  if (containerHasElement(m_varnames, "magvel")) {
254 #ifdef _OPENMP
255 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
256 #endif
257  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
258  const amrex::Box& tbx = mfi.tilebox();
259  auto const& dfab = mf_cc_data.array(mfi);
260  auto const& vfab = mf_cc_vel.array(mfi);
261 
262  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
263  {
264  dfab(i,j,k,mf_comp) = std::sqrt(vfab(i,j,k,0)*vfab(i,j,k,0)
265  + vfab(i,j,k,1)*vfab(i,j,k,1)
266  + vfab(i,j,k,2)*vfab(i,j,k,2)) ;
267  });
268  }
269  mf_comp += 1;
270  }
271 
272  if (containerHasElement(m_varnames, "theta")) {
273 #ifdef _OPENMP
274 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
275 #endif
276  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
277  const amrex::Box& tbx = mfi.tilebox();
278  auto const& dfab = mf_cc_data.array(mfi);
279  auto const& cfab = vars_new[ilev][Vars::cons].array(mfi);
280 
281  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
282  {
283  dfab(i,j,k,mf_comp) = cfab(i,j,k,RhoTheta_comp) / cfab(i,j,k,Rho_comp);
284  });
285  }
286  mf_comp += 1;
287  }
288 
289  if (containerHasElement(m_varnames, "qv")) {
290  // if qv is requested, assume that we have moisture
291  if (qv_comp >= 0) AMREX_ALWAYS_ASSERT(qv_comp == mf_comp);
292  qv_comp = mf_comp;
293 #ifdef _OPENMP
294 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
295 #endif
296  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
297  const amrex::Box& tbx = mfi.tilebox();
298  auto const& dfab = mf_cc_data.array(mfi);
299  auto const& cfab = vars_new[ilev][Vars::cons].array(mfi);
300 
301  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
302  {
303  dfab(i,j,k,mf_comp) = cfab(i,j,k,RhoQ1_comp) / cfab(i,j,k,Rho_comp);
304  });
305  }
306  mf_comp += 1;
307  }
308  if (containerHasElement(m_varnames, "qc")) {
309 #ifdef _OPENMP
310 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
311 #endif
312  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
313  const amrex::Box& tbx = mfi.tilebox();
314  auto const& dfab = mf_cc_data.array(mfi);
315  auto const& cfab = vars_new[ilev][Vars::cons].array(mfi);
316 
317  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
318  {
319  dfab(i,j,k,mf_comp) = cfab(i,j,k,RhoQ2_comp) / cfab(i,j,k,Rho_comp);
320  });
321  }
322  mf_comp += 1;
323  }
324 
325  if (containerHasElement(m_varnames, "pressure") && (qv_comp >= 0)) {
326  // if qv is requested, assume that we have moisture
327 #ifdef _OPENMP
328 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
329 #endif
330  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
331  const amrex::Box& tbx = mfi.tilebox();
332  auto const& dfab = mf_cc_data.array(mfi);
333  auto const& cfab = vars_new[ilev][Vars::cons].array(mfi);
334 
335  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
336  {
337  dfab(i,j,k,mf_comp) = getPgivenRTh(cfab(i,j,k,RhoTheta_comp),
338  dfab(i,j,k,qv_comp));
339  });
340  }
341  mf_comp += 1;
342  }
343  else if (containerHasElement(m_varnames, "pressure")) {
344  // pressure requested w/o qv, assume dry
345 #ifdef _OPENMP
346 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
347 #endif
348  for (amrex::MFIter mfi(mf_cc_data, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
349  const amrex::Box& tbx = mfi.tilebox();
350  auto const& dfab = mf_cc_data.array(mfi);
351  auto const& cfab = vars_new[ilev][Vars::cons].array(mfi);
352 
353  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
354  {
355  dfab(i,j,k,mf_comp) = getPgivenRTh(cfab(i,j,k,RhoTheta_comp));
356  });
357  }
358  mf_comp += 1;
359  }
360 
361  m_lev[iline] = ilev;
362  m_ls_mf[iline] = get_line_data(mf_cc_data, dir, cell, bnd_bx);
363 
364  // We can stop if we got the entire line
365  auto min_bnd_bx = m_ls_mf[iline].boxArray().minimalBox();
366  if (bnd_bx == min_bnd_bx) { break; }
367 
368  } // ilev
369  }// iline
370  }
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
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
Here is the call graph for this function:

◆ write_line_ascii()

void LineSampler::write_line_ascii ( amrex::Vector< amrex::Real > &  time)
inline
387  {
388  // same as definitions in ERF.H
389  constexpr int datwidth = 14;
390  constexpr int datprecision = 6;
391  constexpr int timeprecision = 13;
392 
393  int nline = static_cast<int>(m_ls_mf.size());
394  int nvar = static_cast<int>(m_varnames.size());
395  for (int iline(0); iline<nline; ++iline) {
396  int dir = m_dir[iline];
397  int lev = m_lev[iline];
398  amrex::Real m_time = time[lev];
399  amrex::Box m_dom = m_bnd_bx[iline];
400 
401  for (int ivar(0); ivar<nvar; ++ivar) {
402  // Convert multifab to vector
403  amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(m_ls_mf[iline], ivar, 1, m_dom, dir);
404 
405  if (amrex::ParallelDescriptor::IOProcessor()) {
406  int ifile = iline*nvar + ivar;
407  std::ostream& fs = *m_datastream[ifile];
408  fs << std::setw(datwidth) << std::setprecision(timeprecision) << m_time
409  << std::setw(datwidth) << std::setprecision(datprecision);
410  for (const auto& val : vec) {
411  fs << " " << val;
412  }
413  fs << std::endl;
414  }
415  }
416  }
417  }
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ fs
Definition: ERF_AdvanceMorrison.cpp:117

Referenced by write_line_mfs().

Here is the caller graph for this function:

◆ write_line_coords()

void LineSampler::write_line_coords ( amrex::Vector< std::unique_ptr< amrex::MultiFab > > &  z_phys_cc)
inline
153  {
154  if (!m_write_ascii) return;
155 
156  amrex::Print() << "Writing out line coordinates to text" << std::endl;
157 
158  for (int lev(0); lev < z_phys_cc.size(); ++lev) {
159  // Write one text file per level
160  std::ofstream outfile;
161  if (amrex::ParallelDescriptor::IOProcessor()) {
162  std::string fname = amrex::Concatenate("plt_line_lev", lev, 1);
163  fname += ".zcc";
164  outfile.open(fname);
165 
166  if (!outfile.is_open()) {
167  amrex::AllPrint() << "Could not open " << fname << std::endl;
168  }
169  }
170 
171  // Loop over each line
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();
177 
178  // Create multifab with "sampled" z_phys values
179  amrex::MultiFab line_coords_mf = get_line_data(
180  *z_phys_cc[lev], dir, first_cell, bnd_bx
181  );
182 
183  // Convert multifab to vector
184  amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(
185  line_coords_mf, 0, 1, bnd_bx, dir
186  );
187 
188  // Append to file
189  if (amrex::ParallelDescriptor::IOProcessor()) {
190  for (const auto& zval : vec) {
191  outfile << " " << zval;
192  }
193  outfile << std::endl;
194  }
195  } // line loop
196 
197  outfile.close();
198  } // level loop
199  }

◆ write_line_mfs()

void LineSampler::write_line_mfs ( amrex::Vector< amrex::Real > &  time,
amrex::Vector< int > &  level_steps,
amrex::Vector< amrex::IntVect > &  ref_ratio,
amrex::Vector< amrex::Geometry > &  geom 
)
inline
377  {
378  if (m_write_ascii) {
379  write_line_ascii(time);
380  } else {
381  write_line_plotfile(time, level_steps, ref_ratio, geom);
382  }
383  }
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:420
void write_line_ascii(amrex::Vector< amrex::Real > &time)
Definition: ERF_SampleData.H:386
Here is the call 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
424  {
425  int nline = static_cast<int>(m_ls_mf.size());
426  for (int iline(0); iline<nline; ++iline) {
427  // Data members that can be used as-is
428  int dir = m_dir[iline];
429  int lev = m_lev[iline];
430  amrex::Real m_time = time[lev];
431  amrex::Vector<int> m_level_steps = {level_steps[lev]};
432  amrex::Vector<amrex::IntVect> m_ref_ratio = {ref_ratio[lev]};
433 
434  // Create modified geometry object corresponding to the line
435  auto plo = geom[lev].ProbLo();
436  auto dx = geom[lev].CellSize();
437  amrex::Vector<amrex::Geometry> m_geom; m_geom.resize(1);
438  amrex::Vector<int> is_per(AMREX_SPACEDIM,0);
439  amrex::Box m_dom = m_bnd_bx[iline];
440  amrex::RealBox m_rb;
441  for (int d(0); d<AMREX_SPACEDIM; ++d) {
442  amrex::Real offset = (d==dir) ? 0 : 0.5;
443  amrex::Real lo = plo[d] + ( m_dom.smallEnd(d) - offset ) * dx[d];
444  amrex::Real hi = plo[d] + ( m_dom.bigEnd(d) + offset ) * dx[d];
445 
446  m_rb.setLo(d,lo);
447  m_rb.setHi(d,hi);
448 
449  is_per[d] = geom[lev].isPeriodic(d);
450  }
451  m_geom[0].define(m_dom, &m_rb, geom[lev].Coord(), is_per.data());
452 
453  // Create plotfile name
454  std::string name_line = m_name[iline];
455  name_line += "_step_";
456  std::string plotfilename = amrex::Concatenate(name_line, m_level_steps[0], 5);
457 
458  // Get the data
459  amrex::Vector<const amrex::MultiFab*> mf = {&(m_ls_mf[iline])};
460 
461  // Write each line
462  WriteMultiLevelPlotfile(plotfilename, 1, mf,
463  m_varnames, m_geom, m_time,
464  m_level_steps, m_ref_ratio);
465  }
466  }
Coord
Definition: ERF_DataStruct.H:85
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28

Referenced by write_line_mfs().

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

Member Data Documentation

◆ m_bnd_bx

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

◆ 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_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: