ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
WriteBndryPlanes Class Reference

#include <ERF_WriteBndryPlanes.H>

Collaboration diagram for WriteBndryPlanes:

Public Member Functions

 WriteBndryPlanes (amrex::Vector< amrex::BoxArray > &grids, amrex::Vector< amrex::Geometry > &geom)
 
void write_planes (int t_step, amrex::Real time, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_new, bool is_moist)
 

Private Attributes

amrex::Box target_box
 IO output box region. More...
 
amrex::Vector< amrex::Geometry > & m_geom
 Geometry objects for all levels. More...
 
std::string m_filename {""}
 File name for IO. More...
 
std::string m_time_file {""}
 File name for Native time file. More...
 
amrex::Vector< std::string > m_var_names
 Variables for IO. More...
 
amrex::Vector< amrex::Realm_in_times
 Timestep and times to be stored in time.dat. More...
 
amrex::Vector< int > m_in_timesteps
 
int m_in_rad = 1
 controls extents on native bndry output More...
 
const int m_out_rad = 1
 
const int m_extent_rad = 0
 

Static Private Attributes

static int bndry_lev = 0
 

Detailed Description

Interface for writing boundary planes

This class performs the necessary file operations to write boundary planes

Constructor & Destructor Documentation

◆ WriteBndryPlanes()

WriteBndryPlanes::WriteBndryPlanes ( amrex::Vector< amrex::BoxArray > &  grids,
amrex::Vector< amrex::Geometry > &  geom 
)
explicit

Constructor for the WriteBndryPlanes class for writing the contents of boundary planes to an output file

Parameters
gridsVector of BoxArrays containing the grids at each level in the AMR
geomVector of Geometry containing the geometry at each level in the AMR
55  : m_geom(geom)
56 {
57  ParmParse pp("erf");
58 
59  // Get the radius inside the domain
60  pp.query("in_rad",m_in_rad);
61 
62  // User-specified region is given in physical coordinates, not index space
63  std::vector<Real> box_lo(3), box_hi(3);
64  pp.getarr("bndry_output_box_lo",box_lo,0,2);
65  pp.getarr("bndry_output_box_hi",box_hi,0,2);
66 
67  // If the target area is contained at a finer level, use the finest data possible
68  for (int ilev = 0; ilev < grids.size(); ilev++) {
69 
70  const Real* xLo = m_geom[ilev].ProbLo();
71  auto const dxi = geom[ilev].InvCellSizeArray();
72  const Box& domain = m_geom[ilev].Domain();
73 
74  // We create the smallest box that contains all of the cell centers
75  // in the physical region specified
76  int ilo = static_cast<int>(Math::floor((box_lo[0] - xLo[0]) * dxi[0])+.5);
77  int jlo = static_cast<int>(Math::floor((box_lo[1] - xLo[1]) * dxi[1])+.5);
78  int ihi = static_cast<int>(Math::floor((box_hi[0] - xLo[0]) * dxi[0])+.5)-1;
79  int jhi = static_cast<int>(Math::floor((box_hi[1] - xLo[1]) * dxi[1])+.5)-1;
80 
81  // Map this to index space -- for now we do no interpolation
82  target_box.setSmall(IntVect(ilo,jlo,0));
83  target_box.setBig(IntVect(ihi,jhi,domain.bigEnd(2)));
84 
85  // Test if the target box at this level fits in the grids at this level
86  Box gbx = target_box; gbx.grow(IntVect(1,1,0));
87 
88  // Ensure that the box is no larger than can fit in the (periodically grown) domain
89  // at level 0
90  if (ilev == 0) {
91  Box per_grown_domain(domain);
92  int growx = (geom[0].isPeriodic(0)) ? 1 : 0;
93  int growy = (geom[0].isPeriodic(1)) ? 1 : 0;
94  per_grown_domain.grow(IntVect(growx,growy,0));
95  /*
96  if (!per_grown_domain.contains(gbx))
97  Error("WriteBndryPlanes: Requested box is too large to fill");
98  */
99  }
100 
101  if (grids[ilev].contains(gbx)) bndry_lev = ilev;
102  }
103 
104  // The folder "m_filename" will contain the time series of data and the time.dat file
105  pp.get("bndry_output_planes_file", m_filename);
106 
107  m_time_file = m_filename + "/time.dat";
108 
109  if (pp.contains("bndry_output_var_names"))
110  {
111  int num_vars = pp.countval("bndry_output_var_names");
112  m_var_names.resize(num_vars);
113  pp.queryarr("bndry_output_var_names",m_var_names,0,num_vars);
114  }
115 }
ParmParse pp("prob")
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Vector< amrex::Geometry > & m_geom
Geometry objects for all levels.
Definition: ERF_WriteBndryPlanes.H:30
int m_in_rad
controls extents on native bndry output
Definition: ERF_WriteBndryPlanes.H:50
static int bndry_lev
Definition: ERF_WriteBndryPlanes.H:47
amrex::Box target_box
IO output box region.
Definition: ERF_WriteBndryPlanes.H:27
amrex::Vector< std::string > m_var_names
Variables for IO.
Definition: ERF_WriteBndryPlanes.H:39
std::string m_time_file
File name for Native time file.
Definition: ERF_WriteBndryPlanes.H:36
std::string m_filename
File name for IO.
Definition: ERF_WriteBndryPlanes.H:33
Here is the call graph for this function:

Member Function Documentation

◆ write_planes()

void WriteBndryPlanes::write_planes ( int  t_step,
amrex::Real  time,
amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars_new,
bool  is_moist 
)

Function to write the specified grid data to an output file

Parameters
t_stepTimestep number
timeCurrent time
vars_newGrid data for all variables across the AMR hierarchy
127 {
128  BL_PROFILE("ERF::WriteBndryPlanes::write_planes");
129 
130  MultiFab& S = vars_new[bndry_lev][Vars::cons];
131  MultiFab& xvel = vars_new[bndry_lev][Vars::xvel];
132  MultiFab& yvel = vars_new[bndry_lev][Vars::yvel];
133  MultiFab& zvel = vars_new[bndry_lev][Vars::zvel];
134 
135  const std::string chkname =
136  m_filename + Concatenate("/bndry_output", t_step);
137 
138  //Print() << "Writing boundary planes at time " << time << std::endl;
139 
140  const std::string level_prefix = "Level_";
141  PreBuildDirectorHierarchy(chkname, level_prefix, 1, true);
142 
143  // note: by using the entire domain box we end up using 1 processor
144  // to hold all boundaries
145  BoxArray ba(target_box);
146  DistributionMapping dm{ba};
147 
148  IntVect new_hi = target_box.bigEnd() - target_box.smallEnd();
149  Box target_box_shifted(IntVect(0,0,0),new_hi);
150  BoxArray ba_shifted(target_box_shifted);
151 
152  for (int i = 0; i < m_var_names.size(); i++)
153  {
154  std::string var_name = m_var_names[i];
155  std::string filename = MultiFabFileFullPrefix(bndry_lev, chkname, level_prefix, var_name);
156 
157  int ncomp;
158  if (var_name == "velocity") {
159  ncomp = AMREX_SPACEDIM;
160  } else {
161  ncomp = 1;
162  }
163 
164  BndryRegister bndry (ba , dm, m_in_rad, m_out_rad, m_extent_rad, ncomp);
165  BndryRegister bndry_shifted(ba_shifted, dm, m_in_rad, m_out_rad, m_extent_rad, ncomp);
166 
167  int nghost = 0;
168  if (var_name == "density")
169  {
170  bndry.copyFrom(S, nghost, Rho_comp, 0, ncomp, m_geom[bndry_lev].periodicity());
171 
172  } else if (var_name == "temperature") {
173  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
174  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
175  {
176  const Box& bx = mfi.tilebox();
177  if (is_moist) {
178  derived::erf_dermoisttemp(bx, Temp[mfi], 0, 1, S[mfi], m_geom[bndry_lev], time, nullptr, bndry_lev);
179  } else {
180  derived::erf_dertemp(bx, Temp[mfi], 0, 1, S[mfi], m_geom[bndry_lev], time, nullptr, bndry_lev);
181  }
182  }
183  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
184  } else if (var_name == "theta") {
185  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
186  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
187  {
188  const Box& bx = mfi.tilebox();
189  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoTheta_comp);
190  }
191  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
192  } else if (var_name == "scalar") {
193  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
194  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
195  {
196  const Box& bx = mfi.tilebox();
197  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoScalar_comp);
198  }
199  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
200  } else if (var_name == "ke") {
201  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
202  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
203  {
204  const Box& bx = mfi.tilebox();
205  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoKE_comp);
206  }
207  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
208  } else if (var_name == "qv") {
209  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
210  if (S.nComp() > RhoQ2_comp) {
211  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
212  {
213  const Box& bx = mfi.tilebox();
214  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoQ1_comp);
215  }
216  } else {
217  Temp.setVal(0.);
218  }
219  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
220  } else if (var_name == "qc") {
221  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
222  if (S.nComp() > RhoQ2_comp) {
223  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
224  {
225  const Box& bx = mfi.tilebox();
226  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoQ2_comp);
227  }
228  } else {
229  Temp.setVal(0.);
230  }
231  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
232  } else if (var_name == "velocity") {
233  MultiFab Vel(S.boxArray(), S.DistributionMap(), 3, m_out_rad);
234  average_face_to_cellcenter(Vel,0,Array<const MultiFab*,3>{&xvel,&yvel,&zvel});
235  bndry.copyFrom(Vel, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
236  } else {
237  //Print() << "Trying to write planar output for " << var_name << std::endl;
238  Error("Don't know how to output this variable");
239  }
240 
241  for (OrientationIter oit; oit != nullptr; ++oit) {
242  auto ori = oit();
243  if (ori.coordDir() < 2) {
244  std::string facename = Concatenate(filename + '_', ori, 1);
245  br_shift(oit, bndry, bndry_shifted);
246  bndry_shifted[ori].write(facename);
247  }
248  }
249 
250  } // loop over num_vars
251 
252  // Writing time.dat
253  if (ParallelDescriptor::IOProcessor()) {
254  std::ofstream oftime(m_time_file, std::ios::out | std::ios::app);
255  oftime << std::setprecision(17) << t_step << ' ' << time << '\n';
256  oftime.close();
257  }
258 }
struct @21 out
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#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
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
void br_shift(OrientationIter oit, const BndryRegister &b1, BndryRegister &b2)
Definition: ERF_WriteBndryPlanes.cpp:18
const int m_out_rad
Definition: ERF_WriteBndryPlanes.H:51
const int m_extent_rad
Definition: ERF_WriteBndryPlanes.H:52
@ xvel
Definition: ERF_IndexDefines.H:157
@ cons
Definition: ERF_IndexDefines.H:156
@ zvel
Definition: ERF_IndexDefines.H:159
@ yvel
Definition: ERF_IndexDefines.H:158
void erf_derrhodivide(const Box &bx, FArrayBox &derfab, const FArrayBox &datfab, const int scalar_index)
Definition: ERF_Derive.cpp:18
void erf_dertemp(const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const Geometry &, Real, const int *, const int)
Definition: ERF_Derive.cpp:91
void erf_dermoisttemp(const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const Geometry &, Real, const int *, const int)
Definition: ERF_Derive.cpp:113
Here is the call graph for this function:

Member Data Documentation

◆ bndry_lev

int WriteBndryPlanes::bndry_lev = 0
staticprivate

Referenced by write_planes(), and WriteBndryPlanes().

◆ m_extent_rad

const int WriteBndryPlanes::m_extent_rad = 0
private

Referenced by write_planes().

◆ m_filename

std::string WriteBndryPlanes::m_filename {""}
private

File name for IO.

Referenced by write_planes(), and WriteBndryPlanes().

◆ m_geom

amrex::Vector<amrex::Geometry>& WriteBndryPlanes::m_geom
private

Geometry objects for all levels.

Referenced by write_planes(), and WriteBndryPlanes().

◆ m_in_rad

int WriteBndryPlanes::m_in_rad = 1
private

controls extents on native bndry output

Referenced by write_planes(), and WriteBndryPlanes().

◆ m_in_times

amrex::Vector<amrex::Real> WriteBndryPlanes::m_in_times
private

Timestep and times to be stored in time.dat.

◆ m_in_timesteps

amrex::Vector<int> WriteBndryPlanes::m_in_timesteps
private

◆ m_out_rad

const int WriteBndryPlanes::m_out_rad = 1
private

Referenced by write_planes().

◆ m_time_file

std::string WriteBndryPlanes::m_time_file {""}
private

File name for Native time file.

Referenced by write_planes(), and WriteBndryPlanes().

◆ m_var_names

amrex::Vector<std::string> WriteBndryPlanes::m_var_names
private

Variables for IO.

Referenced by write_planes(), and WriteBndryPlanes().

◆ target_box

amrex::Box WriteBndryPlanes::target_box
private

IO output box region.

Referenced by write_planes(), and WriteBndryPlanes().


The documentation for this class was generated from the following files: