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)
 

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::Real > m_in_times
 Timestep and times to be stored in time.dat. More...
 
amrex::Vector< int > m_in_timesteps
 
const 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  // User-specified region is given in physical coordinates, not index space
60  std::vector<Real> box_lo(3), box_hi(3);
61  pp.getarr("bndry_output_box_lo",box_lo,0,2);
62  pp.getarr("bndry_output_box_hi",box_hi,0,2);
63 
64  // If the target area is contained at a finer level, use the finest data possible
65  for (int ilev = 0; ilev < grids.size(); ilev++) {
66 
67  auto const dxi = geom[ilev].InvCellSizeArray();
68  const Box& domain = m_geom[ilev].Domain();
69 
70  // We create the smallest box that contains all of the cell centers
71  // in the physical region specified
72  int ilo = static_cast<int>(Math::floor(box_lo[0] * dxi[0])+.5);
73  int jlo = static_cast<int>(Math::floor(box_lo[1] * dxi[1])+.5);
74  int ihi = static_cast<int>(Math::floor(box_hi[0] * dxi[0])+.5)-1;
75  int jhi = static_cast<int>(Math::floor(box_hi[1] * dxi[1])+.5)-1;
76 
77  // Map this to index space -- for now we do no interpolation
78  target_box.setSmall(IntVect(ilo,jlo,0));
79  target_box.setBig(IntVect(ihi,jhi,domain.bigEnd(2)));
80 
81  // Test if the target box at this level fits in the grids at this level
82  Box gbx = target_box; gbx.grow(IntVect(1,1,0));
83 
84  // Ensure that the box is no larger than can fit in the (periodically grown) domain
85  // at level 0
86  if (ilev == 0) {
87  Box per_grown_domain(domain);
88  int growx = (geom[0].isPeriodic(0)) ? 1 : 0;
89  int growy = (geom[0].isPeriodic(1)) ? 1 : 0;
90  per_grown_domain.grow(IntVect(growx,growy,0));
91  if (!per_grown_domain.contains(gbx))
92  Error("WriteBndryPlanes: Requested box is too large to fill");
93  }
94 
95  if (grids[ilev].contains(gbx)) bndry_lev = ilev;
96  }
97 
98  // The folder "m_filename" will contain the time series of data and the time.dat file
99  pp.get("bndry_output_planes_file", m_filename);
100 
101  m_time_file = m_filename + "/time.dat";
102 
103  if (pp.contains("bndry_output_var_names"))
104  {
105  int num_vars = pp.countval("bndry_output_var_names");
106  m_var_names.resize(num_vars);
107  pp.queryarr("bndry_output_var_names",m_var_names,0,num_vars);
108  }
109 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: Microphysics_Utils.H:183
amrex::Vector< amrex::Geometry > & m_geom
Geometry objects for all levels.
Definition: ERF_WriteBndryPlanes.H:29
static int bndry_lev
Definition: ERF_WriteBndryPlanes.H:46
amrex::Box target_box
IO output box region.
Definition: ERF_WriteBndryPlanes.H:26
amrex::Vector< std::string > m_var_names
Variables for IO.
Definition: ERF_WriteBndryPlanes.H:38
std::string m_time_file
File name for Native time file.
Definition: ERF_WriteBndryPlanes.H:35
std::string m_filename
File name for IO.
Definition: ERF_WriteBndryPlanes.H:32
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 
)

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

const int WriteBndryPlanes::m_in_rad = 1
private

controls extents on native bndry output

Referenced by write_planes().

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