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
 
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  const Real* xLo = m_geom[ilev].ProbLo();
68  auto const dxi = geom[ilev].InvCellSizeArray();
69  const Box& domain = m_geom[ilev].Domain();
70 
71  // We create the smallest box that contains all of the cell centers
72  // in the physical region specified
73  int ilo = static_cast<int>(Math::floor((box_lo[0] - xLo[0]) * dxi[0])+.5);
74  int jlo = static_cast<int>(Math::floor((box_lo[1] - xLo[1]) * dxi[1])+.5);
75  int ihi = static_cast<int>(Math::floor((box_hi[0] - xLo[0]) * dxi[0])+.5)-1;
76  int jhi = static_cast<int>(Math::floor((box_hi[1] - xLo[1]) * dxi[1])+.5)-1;
77 
78  // Map this to index space -- for now we do no interpolation
79  target_box.setSmall(IntVect(ilo,jlo,0));
80  target_box.setBig(IntVect(ihi,jhi,domain.bigEnd(2)));
81 
82  // Test if the target box at this level fits in the grids at this level
83  Box gbx = target_box; gbx.grow(IntVect(1,1,0));
84 
85  // Ensure that the box is no larger than can fit in the (periodically grown) domain
86  // at level 0
87  if (ilev == 0) {
88  Box per_grown_domain(domain);
89  int growx = (geom[0].isPeriodic(0)) ? 1 : 0;
90  int growy = (geom[0].isPeriodic(1)) ? 1 : 0;
91  per_grown_domain.grow(IntVect(growx,growy,0));
92  if (!per_grown_domain.contains(gbx))
93  Error("WriteBndryPlanes: Requested box is too large to fill");
94  }
95 
96  if (grids[ilev].contains(gbx)) bndry_lev = ilev;
97  }
98 
99  // The folder "m_filename" will contain the time series of data and the time.dat file
100  pp.get("bndry_output_planes_file", m_filename);
101 
102  m_time_file = m_filename + "/time.dat";
103 
104  if (pp.contains("bndry_output_var_names"))
105  {
106  int num_vars = pp.countval("bndry_output_var_names");
107  m_var_names.resize(num_vars);
108  pp.queryarr("bndry_output_var_names",m_var_names,0,num_vars);
109  }
110 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:233
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Vector< amrex::Geometry > & m_geom
Geometry objects for all levels.
Definition: ERF_WriteBndryPlanes.H:30
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
122 {
123  BL_PROFILE("ERF::WriteBndryPlanes::write_planes");
124 
125  MultiFab& S = vars_new[bndry_lev][Vars::cons];
126  MultiFab& xvel = vars_new[bndry_lev][Vars::xvel];
127  MultiFab& yvel = vars_new[bndry_lev][Vars::yvel];
128  MultiFab& zvel = vars_new[bndry_lev][Vars::zvel];
129 
130  const std::string chkname =
131  m_filename + Concatenate("/bndry_output", t_step);
132 
133  //Print() << "Writing boundary planes at time " << time << std::endl;
134 
135  const std::string level_prefix = "Level_";
136  PreBuildDirectorHierarchy(chkname, level_prefix, 1, true);
137 
138  // note: by using the entire domain box we end up using 1 processor
139  // to hold all boundaries
140  BoxArray ba(target_box);
141  DistributionMapping dm{ba};
142 
143  IntVect new_hi = target_box.bigEnd() - target_box.smallEnd();
144  Box target_box_shifted(IntVect(0,0,0),new_hi);
145  BoxArray ba_shifted(target_box_shifted);
146 
147  for (int i = 0; i < m_var_names.size(); i++)
148  {
149  std::string var_name = m_var_names[i];
150  std::string filename = MultiFabFileFullPrefix(bndry_lev, chkname, level_prefix, var_name);
151 
152  int ncomp;
153  if (var_name == "velocity") {
154  ncomp = AMREX_SPACEDIM;
155  } else {
156  ncomp = 1;
157  }
158 
159  BndryRegister bndry (ba , dm, m_in_rad, m_out_rad, m_extent_rad, ncomp);
160  BndryRegister bndry_shifted(ba_shifted, dm, m_in_rad, m_out_rad, m_extent_rad, ncomp);
161 
162  int nghost = 0;
163  if (var_name == "density")
164  {
165  bndry.copyFrom(S, nghost, Rho_comp, 0, ncomp, m_geom[bndry_lev].periodicity());
166 
167  } else if (var_name == "temperature") {
168  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
169  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
170  {
171  const Box& bx = mfi.tilebox();
172  if (is_moist) {
173  derived::erf_dermoisttemp(bx, Temp[mfi], 0, 1, S[mfi], m_geom[bndry_lev], time, nullptr, bndry_lev);
174  } else {
175  derived::erf_dertemp(bx, Temp[mfi], 0, 1, S[mfi], m_geom[bndry_lev], time, nullptr, bndry_lev);
176  }
177  }
178  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
179  } else if (var_name == "theta") {
180  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
181  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
182  {
183  const Box& bx = mfi.tilebox();
184  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoTheta_comp);
185  }
186  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
187  } else if (var_name == "scalar") {
188  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
189  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
190  {
191  const Box& bx = mfi.tilebox();
192  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoScalar_comp);
193  }
194  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
195  } else if (var_name == "ke") {
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], RhoKE_comp);
201  }
202  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
203  } else if (var_name == "qv") {
204  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
205  if (S.nComp() > RhoQ2_comp) {
206  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
207  {
208  const Box& bx = mfi.tilebox();
209  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoQ1_comp);
210  }
211  } else {
212  Temp.setVal(0.);
213  }
214  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
215  } else if (var_name == "qc") {
216  MultiFab Temp(S.boxArray(),S.DistributionMap(),ncomp,0);
217  if (S.nComp() > RhoQ2_comp) {
218  for (MFIter mfi(Temp, TilingIfNotGPU()); mfi.isValid(); ++mfi)
219  {
220  const Box& bx = mfi.tilebox();
221  derived::erf_derrhodivide(bx, Temp[mfi], S[mfi], RhoQ2_comp);
222  }
223  } else {
224  Temp.setVal(0.);
225  }
226  bndry.copyFrom(Temp, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
227  } else if (var_name == "velocity") {
228  MultiFab Vel(S.boxArray(), S.DistributionMap(), 3, m_out_rad);
229  average_face_to_cellcenter(Vel,0,Array<const MultiFab*,3>{&xvel,&yvel,&zvel});
230  bndry.copyFrom(Vel, nghost, 0, 0, ncomp, m_geom[bndry_lev].periodicity());
231  } else {
232  //Print() << "Trying to write planar output for " << var_name << std::endl;
233  Error("Don't know how to output this variable");
234  }
235 
236  for (OrientationIter oit; oit != nullptr; ++oit) {
237  auto ori = oit();
238  if (ori.coordDir() < 2) {
239  std::string facename = Concatenate(filename + '_', ori, 1);
240  br_shift(oit, bndry, bndry_shifted);
241  bndry_shifted[ori].write(facename);
242  }
243  }
244 
245  } // loop over num_vars
246 
247  // Writing time.dat
248  if (ParallelDescriptor::IOProcessor()) {
249  std::ofstream oftime(m_time_file, std::ios::out | std::ios::app);
250  oftime << t_step << ' ' << time << '\n';
251  oftime.close();
252  }
253 }
struct @19 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
const int m_in_rad
controls extents on native bndry output
Definition: ERF_WriteBndryPlanes.H:50
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
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

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: