ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
LineSampler Struct Reference

#include <ERF_SampleData.H>

Collaboration diagram for LineSampler:

Public Member Functions

 LineSampler ()
 
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, "x_velocity")) {
91  m_varnames.push_back("x_velocity");
92  amrex::Print() << " " << "x_velocity";
93  }
94  if (containerHasElement(requested_vars, "y_velocity")) {
95  m_varnames.push_back("y_velocity");
96  amrex::Print() << " " << "y_velocity";
97  }
98  if (containerHasElement(requested_vars, "z_velocity")) {
99  m_varnames.push_back("z_velocity");
100  amrex::Print() << " " << "z_velocity";
101  }
102  if (containerHasElement(requested_vars, "magvel")) {
103  m_varnames.push_back("magvel");
104  amrex::Print() << " " << "magvel";
105  }
106  if (containerHasElement(requested_vars, "theta")) {
107  m_varnames.push_back("theta");
108  amrex::Print() << " " << "theta";
109  }
110  if (containerHasElement(requested_vars, "qv")) {
111  m_varnames.push_back("qv");
112  amrex::Print() << " " << "qv";
113  }
114  if (containerHasElement(requested_vars, "pressure")) {
115  m_varnames.push_back("pressure");
116  amrex::Print() << " " << "pressure";
117  }
118  amrex::Print() << std::endl;
119  }
120 
121  // Write outputs to text files, one file per variable, with all
122  // times appended to the same file
123  pp.query("line_sampling_text_output",m_write_ascii);
124  if (m_write_ascii && amrex::ParallelDescriptor::IOProcessor()) {
125  int nvar = m_varnames.size();
126  m_datastream.resize(n_line_lo * nvar);
127  int i = 0;
128  for (int iline(0); iline<n_line_lo; ++iline) {
129  for (int ivar(0); ivar<nvar; ++ivar) {
130  std::string filename = m_name[iline] + "." + m_varnames[ivar];
131  m_datastream[i] = std::make_unique<std::fstream>();
132  m_datastream[i]->open(filename.c_str(),std::ios::out|std::ios::app);
133  if (!m_datastream[i]->good()) {
134  amrex::FileOpenFailed(filename);
135  }
136  i++;
137  }
138  }
139  }
140  }
141  }
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:219
amrex::Vector< std::string > m_varnames
Definition: ERF_SampleData.H:390
amrex::Vector< amrex::MultiFab > m_ls_mf
Definition: ERF_SampleData.H:386
amrex::Vector< std::unique_ptr< std::fstream > > m_datastream
Definition: ERF_SampleData.H:391
amrex::Vector< int > m_dir
Definition: ERF_SampleData.H:383
amrex::Vector< std::string > m_name
Definition: ERF_SampleData.H:387
amrex::Vector< amrex::Box > m_bnd_bx
Definition: ERF_SampleData.H:385
bool m_write_ascii
Definition: ERF_SampleData.H:389
amrex::Vector< int > m_lev
Definition: ERF_SampleData.H:384
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
145  {
146  int nlev = vars_new.size();
147  int nline = m_bnd_bx.size();
148  int ncomp = m_varnames.size();
149 
150  // Loop over each line
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();
155 
156  // Search each level to get the finest data possible
157  for (int ilev(nlev-1); ilev>=0; --ilev) {
158 
159  // Construct CC velocities
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],
166  &vars_new[ilev][Vars::yvel],
167  &vars_new[ilev][Vars::zvel]});
168 
169  // Construct vector of MFs holding T and WSP
170  amrex::MultiFab mf_cc_data;
171  mf_cc_data.define(ba, dm, ncomp, 1);
172 
173  int mf_comp = 0;
174 
175  if (containerHasElement(m_varnames, "x_velocity")) {
176  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 0, mf_comp, 1, 0);
177  mf_comp += 1;
178  }
179  if (containerHasElement(m_varnames, "y_velocity")) {
180  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 1, mf_comp, 1, 0);
181  mf_comp += 1;
182  }
183  if (containerHasElement(m_varnames, "z_velocity")) {
184  amrex::MultiFab::Copy(mf_cc_data, mf_cc_vel, 2, mf_comp, 1, 0);
185  mf_comp += 1;
186  }
187 
188  if (containerHasElement(m_varnames, "magvel")) {
189 #ifdef _OPENMP
190 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
191 #endif
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);
196 
197  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
198  {
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)) ;
202  });
203  }
204  mf_comp += 1;
205  }
206 
207  if (containerHasElement(m_varnames, "theta")) {
208 #ifdef _OPENMP
209 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
210 #endif
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);
215 
216  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
217  {
218  dfab(i,j,k,mf_comp) = cfab(i,j,k,RhoTheta_comp) / cfab(i,j,k,Rho_comp);
219  });
220  }
221  mf_comp += 1;
222  }
223 
224  if (containerHasElement(m_varnames, "qv") && containerHasElement(m_varnames, "pressure")) {
225  // if qv is requested, assume that we have moisture
226 #ifdef _OPENMP
227 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
228 #endif
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);
233 
234  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
235  {
236  dfab(i,j,k,mf_comp ) = cfab(i,j,k,RhoQ1_comp) / cfab(i,j,k,Rho_comp);
237  dfab(i,j,k,mf_comp+1) = getPgivenRTh(cfab(i,j,k,RhoTheta_comp),
238  dfab(i,j,k,mf_comp));
239  });
240  }
241  mf_comp += 2;
242  } else if (containerHasElement(m_varnames, "qv")) {
243  // if qv is requested, assume that we have moisture
244 #ifdef _OPENMP
245 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
246 #endif
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);
251 
252  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
253  {
254  dfab(i,j,k,mf_comp ) = cfab(i,j,k,RhoQ1_comp) / cfab(i,j,k,Rho_comp);
255  });
256  }
257  mf_comp += 1;
258  } else if (containerHasElement(m_varnames, "pressure")) {
259  // pressure requested w/o qv, assume dry
260 #ifdef _OPENMP
261 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
262 #endif
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);
267 
268  amrex::ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
269  {
270  dfab(i,j,k,mf_comp) = getPgivenRTh(cfab(i,j,k,RhoTheta_comp));
271  });
272  }
273  mf_comp += 1;
274  }
275 
276  m_lev[iline] = ilev;
277  m_ls_mf[iline] = get_line_data(mf_cc_data, dir, cell, bnd_bx);
278 
279  // We can stop if we got the entire line
280  auto min_bnd_bx = m_ls_mf[iline].boxArray().minimalBox();
281  if (bnd_bx == min_bnd_bx) { continue; }
282 
283  } // ilev
284  }// iline
285  }
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
@ 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
302  {
303  // same as definitions in ERF.H
304  constexpr int datwidth = 14;
305  constexpr int datprecision = 6;
306  constexpr int timeprecision = 13;
307 
308  int nline = m_ls_mf.size();
309  int nvar = m_varnames.size();
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];
314  amrex::Box m_dom = m_bnd_bx[iline];
315 
316  for (int ivar(0); ivar<nvar; ++ivar) {
317  // Convert multifab to vector
318  amrex::Gpu::HostVector<amrex::Real> vec = sumToLine(m_ls_mf[iline], ivar, 1, m_dom, dir);
319 
320  if (amrex::ParallelDescriptor::IOProcessor()) {
321  int ifile = iline*nvar + ivar;
322  std::ostream& fs = *m_datastream[ifile];
323  fs << std::setw(datwidth) << std::setprecision(timeprecision) << m_time
324  << std::setw(datwidth) << std::setprecision(datprecision);
325  for (const auto& val : vec) {
326  fs << " " << val;
327  }
328  fs << std::endl;
329  }
330  }
331  }
332  }

Referenced by write_line_mfs().

Here is the caller graph for this function:

◆ 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
292  {
293  if (m_write_ascii) {
294  write_line_ascii(time);
295  } else {
296  write_line_plotfile(time, level_steps, ref_ratio, geom);
297  }
298  }
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
void write_line_ascii(amrex::Vector< amrex::Real > &time)
Definition: ERF_SampleData.H:301
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
339  {
340  int nline = m_ls_mf.size();
341  for (int iline(0); iline<nline; ++iline) {
342  // Data members that can be used as-is
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]};
348 
349  // Create modified geometry object corresponding to the line
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);
354  amrex::Box m_dom = m_bnd_bx[iline];
355  amrex::RealBox m_rb;
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];
360 
361  m_rb.setLo(d,lo);
362  m_rb.setHi(d,hi);
363 
364  is_per[d] = geom[lev].isPeriodic(d);
365  }
366  m_geom[0].define(m_dom, &m_rb, geom[lev].Coord(), is_per.data());
367 
368  // Create plotfile name
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);
372 
373  // Get the data
374  amrex::Vector<const amrex::MultiFab*> mf = {&(m_ls_mf[iline])};
375 
376  // Write each line
377  WriteMultiLevelPlotfile(plotfilename, 1, mf,
378  m_varnames, m_geom, m_time,
379  m_level_steps, m_ref_ratio);
380  }
381  }
Coord
Definition: ERF_DataStruct.H:68
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}

Referenced by LineSampler(), and write_line_mfs().


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