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

#include <ERF_TurbPertStruct.H>

Collaboration diagram for TurbulentPerturbation:

Public Member Functions

 ~TurbulentPerturbation ()
 
void init_tpi_type (const int lev, const PerturbationType &pert_type, const int max_level)
 
void init_tpi (const int lev, const amrex::Vector< amrex::BoxArray > &subdomains_lev, const amrex::GpuArray< amrex::Real, 3 > dx, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const int ngrow_state, std::string pp_prefix, const amrex::Vector< amrex::IntVect > refRatio, const int max_level)
 
void calc_tpi_update (const int lev, const amrex::Real dt, amrex::MultiFab &mf_xvel, amrex::MultiFab &mf_yvel, amrex::MultiFab &mf_cons)
 
void apply_tpi (const int &lev, const amrex::Box &vbx, const int &comp, const amrex::IndexType &m_ixtype, const amrex::Array4< amrex::Real > &src_arr, const amrex::Array4< amrex::Real const > &pert_cell)
 
void calc_tpi_amp (const int &lev, const int &boxIdx, const amrex::Real &interval)
 
void pseudoRandomPert (const int &boxIdx, const int &lev, const amrex::IndexType &m_ixtype)
 
void zero_amp (const int &boxIdx, const int &lev, const amrex::IndexType &m_ixtype)
 
void netZeroBuoyantAdd (const int &boxIdx, const int &lev)
 
void netZeroBuoyantAdjust (const int &boxIdx, const int &lev)
 
void calc_tpi_meanMag_perBox (const int &boxIdx, const int &lev, amrex::MultiFab &mf_cons, amrex::MultiFab &mf_xvel, amrex::MultiFab &mf_yvel)
 
void debug (amrex::Real)
 

Public Attributes

amrex::Vector< int > pt_type
 
amrex::Vector< amrex::BoxArray > pb_ba
 
amrex::Vector< amrex::Vector< amrex::Real > > pb_mag
 
amrex::Vector< amrex::Vector< amrex::Real > > pb_dir
 
amrex::Vector< amrex::MultiFab > pb_cell
 

Private Member Functions

amrex::Real RandomReal (const amrex::Real min, const amrex::Real max)
 

Private Attributes

int tpi_layers
 
int tpi_offset
 
amrex::Vector< int > tpi_boxDim
 
amrex::Vector< int > tpi_direction
 
amrex::Real tpi_nonDim
 
amrex::Real tpi_Ti
 
amrex::Real tpi_Tinf
 
amrex::Vector< amrex::Realtpi_Hpb
 
amrex::Vector< amrex::Realtpi_Lpb
 
amrex::Vector< amrex::Realtpi_Wpb
 
amrex::Vector< amrex::Realtpi_lref
 
amrex::Real tpi_net_buoyant
 
amrex::Real tpi_pert_adjust
 
amrex::Vector< amrex::IntVect > ref_ratio
 
amrex::Real input_Ug
 
amrex::Vector< amrex::Vector< amrex::Real > > pb_interval
 
amrex::Vector< amrex::Vector< amrex::Real > > pb_local_etime
 
amrex::Vector< amrex::Vector< amrex::Real > > pb_amp
 
amrex::Vector< amrex::Vector< amrex::Real > > pb_netZero
 

Constructor & Destructor Documentation

◆ ~TurbulentPerturbation()

TurbulentPerturbation::~TurbulentPerturbation ( )
inline
26 {}

Member Function Documentation

◆ apply_tpi()

void TurbulentPerturbation::apply_tpi ( const int &  lev,
const amrex::Box &  vbx,
const int &  comp,
const amrex::IndexType &  m_ixtype,
const amrex::Array4< amrex::Real > &  src_arr,
const amrex::Array4< amrex::Real const > &  pert_cell 
)
inline
357  {
358  for (int boxIdx = 0; boxIdx < pb_ba[lev].size(); boxIdx++) {
359  amrex::Box pbx = amrex::convert(pb_ba[lev][boxIdx], m_ixtype);
360  amrex::Box ubx = pbx & vbx;
361  if (ubx.ok()) {
362  ParallelFor(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
363  src_arr(i,j,k,comp) += pert_cell(i,j,k);
364 
365  // For box region debug only
366  #ifdef INDEX_PERTURB
367  src_arr(i,j,k,comp) = (amrex::Real) (boxIdx + amrex::Real(5.));
368  #endif
369  });
370  }
371  }
372  }
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+myhalf) *dx[0];const Real y=(j+myhalf) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - one)/(two *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(one+deltaT, inv_gm1);const Real T=(one+deltaT) *T_inf;const Real p=std::pow(rho_norm, Gamma)/Gamma *rho_0 *a_inf *a_inf;const Real rho_theta=rho_0 *rho_norm *(T *std::pow(p_0/p, rdOcp));state_pert(i, j, k, RhoTheta_comp)=rho_theta - getRhoThetagivenP(p_hse(i, j, k));const Real r2d_xy=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc));state_pert(i, j, k, RhoScalar_comp)=fourth *(one+std::cos(PI *std::min(r2d_xy, R)/R));})
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Vector< amrex::BoxArray > pb_ba
Definition: ERF_TurbPertStruct.H:661

Referenced by make_sources().

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

◆ calc_tpi_amp()

void TurbulentPerturbation::calc_tpi_amp ( const int &  lev,
const int &  boxIdx,
const amrex::Real interval 
)
inline
378  {
379  pb_amp[lev][boxIdx] = zero; // Safety step
380  if (pt_type[lev] == 2) { // CPM
381  amrex::Real cp = 1004; // specific heat of air [J/(kg K)]
382  amrex::Real Ec = amrex::Real(0.2); // Eckert number
383  pb_amp[lev][boxIdx] = (input_Ug * input_Ug) / (Ec * cp);
384  } else { // box perturbation
385  amrex::Real Um = pb_mag[lev][boxIdx];
386  amrex::Real beta = one/tpi_Tinf; // Thermal expansion coefficient
387 
388  // Pseudo Random temperature the ignores scale when mechanically tripping turbulence
390  // get total refinement ratio on each level
391  int total_ref_ratio = 1;
392  for (int level = lev; level >= 1; level--) {
393  total_ref_ratio *= ref_ratio[level-1][2];
394  }
395  // calculation needs to be scale-aware since the formulation relies on the physical size of the box
396  if (tpi_Ti > zero) g = (tpi_nonDim * Um * Um) / (tpi_Ti * tpi_Hpb[lev]) * 1 / total_ref_ratio;
397 
398  // Ma and Senocak (2023) Eq. 8, solving for delta phi
399  pb_amp[lev][boxIdx] = (tpi_nonDim * Um * Um) / (g * beta * tpi_Hpb[lev]) * 1 / total_ref_ratio;
400 
401  if (pt_type[lev] == 0) {
402  // Performing this step converts the perturbation proportionality into
403  // the forcing term
404  // Ma & Senocak (2023) Eq. 7
405  pb_amp[lev][boxIdx] /= interval;
406  }
407  }
408  }
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:32
amrex::Real beta
Definition: ERF_InitCustomPert_DataAssimilation_ISV.H:10
real(c_double), parameter g
Definition: ERF_module_model_constants.F90:19
real(c_double), parameter cp
Definition: ERF_module_model_constants.F90:22
amrex::Real tpi_Tinf
Definition: ERF_TurbPertStruct.H:681
amrex::Vector< amrex::IntVect > ref_ratio
Definition: ERF_TurbPertStruct.H:692
amrex::Real input_Ug
Definition: ERF_TurbPertStruct.H:693
amrex::Vector< amrex::Vector< amrex::Real > > pb_mag
Definition: ERF_TurbPertStruct.H:662
amrex::Vector< int > pt_type
Definition: ERF_TurbPertStruct.H:658
amrex::Vector< amrex::Real > tpi_Hpb
Definition: ERF_TurbPertStruct.H:684
amrex::Real tpi_Ti
Definition: ERF_TurbPertStruct.H:680
amrex::Real tpi_nonDim
Definition: ERF_TurbPertStruct.H:679
amrex::Vector< amrex::Vector< amrex::Real > > pb_amp
Definition: ERF_TurbPertStruct.H:698

Referenced by calc_tpi_update().

Here is the caller graph for this function:

◆ calc_tpi_meanMag_perBox()

void TurbulentPerturbation::calc_tpi_meanMag_perBox ( const int &  boxIdx,
const int &  lev,
amrex::MultiFab &  mf_cons,
amrex::MultiFab &  mf_xvel,
amrex::MultiFab &  mf_yvel 
)
inline
524  {
525  // Creating local copy of PB box array and magnitude
526  const amrex::BoxArray m_pb_ba = pb_ba[lev];
527  amrex::Real* m_pb_mag = pb_mag[lev].data();
528  amrex::Real* m_pb_dir = pb_dir[lev].data();
529  m_pb_mag[boxIdx] = zero; // Safety step
530  m_pb_dir[boxIdx] = zero; // Safety step
531 
532  // Storage of averages per PB
533  // Index: 0=u (vol/slab_lo), 1=v (vol/slab_lo)
534  // 2=u (slab_hi), 3=v (slab_hi)
535  int n_avg = 4;
536  amrex::Vector<amrex::Real> avg_h(n_avg,zero);
537  amrex::Gpu::DeviceVector<amrex::Real> avg_d(n_avg,zero);
538  amrex::Real* avg = avg_d.data();
539 
540  // Averaging u & v components in single MFIter
541  for (amrex::MFIter mfi(mf_cons, TileNoZ()); mfi.isValid(); ++mfi) {
542 
543  // CC valid box (inherited from mf_cons)
544  const amrex::Box& vbx = mfi.validbox();
545 
546  // Box logic for u velocity
547  auto ixtype_u = mf_xvel.boxArray().ixType();
548  amrex::Box vbx_u = amrex::convert(vbx,ixtype_u);
549  amrex::Box pbx_u = amrex::convert(m_pb_ba[boxIdx], ixtype_u);
550  amrex::Box ubx_u = pbx_u & vbx_u;
551 
552  // Box logic for v velocity
553  auto ixtype_v = mf_yvel.boxArray().ixType();
554  amrex::Box vbx_v = amrex::convert(vbx,ixtype_v);
555  amrex::Box pbx_v = amrex::convert(m_pb_ba[boxIdx], ixtype_v);
556  amrex::Box ubx_v = pbx_v & vbx_v;
557 
558  // Operation over box union (U)
559  if (ubx_u.ok()) {
560  const amrex::Array4<const amrex::Real>& xvel_arry = mf_xvel.const_array(mfi);
561 
562  #ifdef USE_VOLUME_AVERAGE
563  amrex::Real norm = one / static_cast<amrex::Real>(ubx_u.numPts());
564  ParallelFor(amrex::Gpu::KernelInfo().setReduction(true), ubx_u, [=]
565  AMREX_GPU_DEVICE(int i, int j, int k, amrex::Gpu::Handler const& handler) noexcept {
566  amrex::Gpu::deviceReduceSum(&avg[0], xvel_arry(i,j,k)*norm, handler);
567  });
568  #endif // USE_VOLUME_AVERAGE
569 
570  #ifdef USE_SLAB_AVERAGE
571  amrex::Box ubxSlab_lo = makeSlab(ubx_u,2,ubx_u.smallEnd(2));
572  amrex::Box ubxSlab_hi = makeSlab(ubx_u,2,ubx_u.bigEnd(2));
573  amrex::Real norm_lo = one / static_cast<amrex::Real>(ubxSlab_lo.numPts());
574  amrex::Real norm_hi = one / static_cast<amrex::Real>(ubxSlab_hi.numPts());
575 
576  // Average u in the low slab
577  ParallelFor(amrex::Gpu::KernelInfo().setReduction(true), ubxSlab_lo, [=]
578  AMREX_GPU_DEVICE(int i, int j, int k, amrex::Gpu::Handler const& handler) noexcept {
579  amrex::Gpu::deviceReduceSum(&avg[0], xvel_arry(i,j,k)*norm_lo, handler);
580  });
581 
582  // Average u in the high slab
583  ParallelFor(amrex::Gpu::KernelInfo().setReduction(true), ubxSlab_hi, [=]
584  AMREX_GPU_DEVICE(int i, int j, int k, amrex::Gpu::Handler const& handler) noexcept {
585  amrex::Gpu::deviceReduceSum(&avg[2], xvel_arry(i,j,k)*norm_hi, handler);
586  });
587  #endif // USE_SLAB_AVERAGE
588  } // if
589 
590  // Operation over box union (V)
591  if (ubx_v.ok()) {
592  const amrex::Array4<const amrex::Real>& yvel_arry = mf_yvel.const_array(mfi);
593 
594  #ifdef USE_VOLUME_AVERAGE
595  amrex::Real norm = one / static_cast<amrex::Real>(ubx_v.numPts());
596  ParallelFor(amrex::Gpu::KernelInfo().setReduction(true), ubx_v, [=]
597  AMREX_GPU_DEVICE(int i, int j, int k, amrex::Gpu::Handler const& handler) noexcept {
598  amrex::Gpu::deviceReduceSum(&avg[1], yvel_arry(i,j,k)*norm, handler);
599  });
600  #endif // USE_VOLUME_AVERAGE
601 
602  #ifdef USE_SLAB_AVERAGE
603  amrex::Box ubxSlab_lo = makeSlab(ubx_v,2,ubx_v.smallEnd(2));
604  amrex::Box ubxSlab_hi = makeSlab(ubx_v,2,ubx_v.bigEnd(2));
605  amrex::Real norm_lo = one / static_cast<amrex::Real>(ubxSlab_lo.numPts());
606  amrex::Real norm_hi = one / static_cast<amrex::Real>(ubxSlab_hi.numPts());
607 
608  // Average v in the low slab
609  ParallelFor(amrex::Gpu::KernelInfo().setReduction(true), ubxSlab_lo, [=]
610  AMREX_GPU_DEVICE(int i, int j, int k, amrex::Gpu::Handler const& handler) noexcept {
611  amrex::Gpu::deviceReduceSum(&avg[1], yvel_arry(i,j,k)*norm_lo, handler);
612  });
613 
614  // Average v in the high slab
615  ParallelFor(amrex::Gpu::KernelInfo().setReduction(true), ubxSlab_hi, [=]
616  AMREX_GPU_DEVICE(int i, int j, int k, amrex::Gpu::Handler const& handler) noexcept {
617  amrex::Gpu::deviceReduceSum(&avg[3], yvel_arry(i,j,k)*norm_hi, handler);
618  });
619  #endif // USE_SLAB_AVERAGE
620  } // if
621  } // MFIter
622 
623  // Copy from device back to host
624  amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(), avg_h.begin());
625 
626  // Computing the average magnitude within PB
627  #ifdef USE_VOLUME_AVERAGE
628  m_pb_mag[boxIdx] = std::sqrt(avg_h[0]*avg_h[0] + avg_h[1]*avg_h[1]);
629  m_pb_dir[boxIdx] = std::atan(std::abs(avg_h[0]) / std::abs(avg_h[1]+std::numeric_limits<amrex::Real>::epsilon()));
630  #endif
631 
632  #ifdef USE_SLAB_AVERAGE
633  m_pb_mag[boxIdx] = myhalf*( std::sqrt(avg_h[0]*avg_h[0] + avg_h[1]*avg_h[1])
634  + std::sqrt(avg_h[2]*avg_h[2] + avg_h[3]*avg_h[3]));
635  m_pb_dir[boxIdx] = std::atan(std::abs(myhalf*(avg_h[0]+avg_h[2])) / std::abs(myhalf*(avg_h[1]+avg_h[3])+std::numeric_limits<amrex::Real>::epsilon()));
636  #endif
637  }
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
amrex::Vector< amrex::Vector< amrex::Real > > pb_dir
Definition: ERF_TurbPertStruct.H:663

Referenced by calc_tpi_update().

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

◆ calc_tpi_update()

void TurbulentPerturbation::calc_tpi_update ( const int  lev,
const amrex::Real  dt,
amrex::MultiFab &  mf_xvel,
amrex::MultiFab &  mf_yvel,
amrex::MultiFab &  mf_cons 
)
inline
255  {
256  // Resetting the net buoyant force value
258 
259  // Setting random number generator for update interval
260  srand( (unsigned) time(NULL) );
261 
262  auto m_ixtype = mf_cons.boxArray().ixType(); // safety step
263 
264  // Seed the random generator at 1024UL for regression testing
265  int fix_random_seed = 0;
266  amrex::ParmParse pp("erf");
267  pp.query("fix_random_seed", fix_random_seed);
268  if (fix_random_seed) {
269  // We need this one for the ParalleForRNG used in calc_tpi
270  amrex::InitRandom(1024UL, amrex::ParallelDescriptor::NProcs(), 1024UL);
271 
272  // We need this one for the RandomReal below
273  srand(1024UL);
274  }
275 
276  for (int boxIdx = 0; boxIdx < pb_ba[lev].size(); boxIdx++) {
277 
278  bool update_box = true; // initialize flag
279  // Check if the local elapsed time is greater than the update interval and don't go into boxes the rank doesn't own
280  if (pt_type[lev] == 2) {
281  // for CPM, perturbation boxes/cells refresh after boxes have advected box width/length * num_layers (advective time scale)
282  update_box = ((pb_local_etime[lev][boxIdx] >= pb_interval[lev][boxIdx]*tpi_layers) && pb_local_etime[lev][boxIdx] != -one);
283  } else {
284  update_box = ( pb_local_etime[lev][boxIdx] >= pb_interval[lev][boxIdx] && pb_local_etime[lev][boxIdx] != -one );
285  }
286  if ( update_box ) {
287 
288  // Compute mean velocity of each perturbation box
289  calc_tpi_meanMag_perBox(boxIdx, lev, mf_cons, mf_xvel, mf_yvel);
290 
291  // Only the rank owning the box will be able to access the storage location
292  // Done for parallelism to avoid Inf being stored in array
293  if (pb_mag[lev][boxIdx] !=zero) {
294  amrex::Real interval = zero;
295  if (pt_type[lev] == 2) {
296  // Wind direction correction for angled wind
297  amrex::Real wind_direction = pb_dir[lev][boxIdx];
298  if (wind_direction > PI / 4) { wind_direction = PI / 2 - wind_direction; }
299  // CPM only cares about the side length of a box, min call maintains flexibility.
300  interval = one / std::cos(wind_direction) * std::min(tpi_Lpb[lev], tpi_Wpb[lev]) / pb_mag[lev][boxIdx];
301  } else {
302  interval = tpi_lref[lev] / pb_mag[lev][boxIdx];
303  }
304  pb_interval[lev][boxIdx] = RandomReal(amrex::Real(0.9)*interval,amrex::Real(1.1)*interval); // 10% variation
305 
306  // Reset local elapsed time
307  pb_local_etime[lev][boxIdx] = zero;
308  } else {
309  // this box is not on this rank, we shouldn't enter it again.
310  // Technically, all boxes are looped through on the very first step of a simulation and this is when this is set
311  pb_local_etime[lev][boxIdx] = -one;
312  }
313 
314  // Trigger amplitude calculation per perturbation box
315  calc_tpi_amp(lev, boxIdx, pb_interval[lev][boxIdx]);
316 
317  // Trigger random amplitude storage per cell within perturbation box
318  pseudoRandomPert(boxIdx, lev, m_ixtype);
319 
320  } else {
321  // set perturbation amplitudes to 0 for CPM. A little inefficient but leverages as much as existing code as possible.
322  if (pt_type[lev] == 2) { zero_amp(boxIdx, lev, m_ixtype); }
323 
324  // Increase by timestep of level 0 (but only if the box is owned by this rank)
325  if (pb_local_etime[lev][boxIdx] != -one) { pb_local_etime[lev][boxIdx] += dt; }
326  } // if
327 
328  if (pt_type[lev] < 2) { // box perturbation method only
329  // Per iteration operation of net-zero buoyant force check
330  if (pb_mag[lev][boxIdx] !=zero) netZeroBuoyantAdd(boxIdx, lev);
331  tpi_net_buoyant += pb_netZero[lev][boxIdx];
332  }
333  } // for
334 
335  if (pt_type[lev] < 2) { // box perturbation method only
336  // Normalizing the adjustment based on how many boxes there are
337  // the values within the array is already normalized by the number
338  // of cells within each box
340 
341  // Per iteration operation of net-zero buoyant force adjustment
342  for (int boxIdx = 0; boxIdx < pb_ba[lev].size(); boxIdx++) {
343  if (pb_mag[lev][boxIdx] !=zero) netZeroBuoyantAdjust(boxIdx, lev);
344  }
345  }
346  }
constexpr amrex::Real PI
Definition: ERF_Constants.H:16
ParmParse pp("prob")
int tpi_layers
Definition: ERF_TurbPertStruct.H:672
void pseudoRandomPert(const int &boxIdx, const int &lev, const amrex::IndexType &m_ixtype)
Definition: ERF_TurbPertStruct.H:413
amrex::Vector< amrex::Vector< amrex::Real > > pb_local_etime
Definition: ERF_TurbPertStruct.H:697
void calc_tpi_amp(const int &lev, const int &boxIdx, const amrex::Real &interval)
Definition: ERF_TurbPertStruct.H:375
void netZeroBuoyantAdjust(const int &boxIdx, const int &lev)
Definition: ERF_TurbPertStruct.H:495
amrex::Vector< amrex::Real > tpi_lref
Definition: ERF_TurbPertStruct.H:687
amrex::Vector< amrex::Vector< amrex::Real > > pb_netZero
Definition: ERF_TurbPertStruct.H:699
amrex::Real tpi_pert_adjust
Definition: ERF_TurbPertStruct.H:690
amrex::Vector< amrex::Vector< amrex::Real > > pb_interval
Definition: ERF_TurbPertStruct.H:696
void zero_amp(const int &boxIdx, const int &lev, const amrex::IndexType &m_ixtype)
Definition: ERF_TurbPertStruct.H:442
amrex::Vector< amrex::Real > tpi_Wpb
Definition: ERF_TurbPertStruct.H:686
void calc_tpi_meanMag_perBox(const int &boxIdx, const int &lev, amrex::MultiFab &mf_cons, amrex::MultiFab &mf_xvel, amrex::MultiFab &mf_yvel)
Definition: ERF_TurbPertStruct.H:518
void netZeroBuoyantAdd(const int &boxIdx, const int &lev)
Definition: ERF_TurbPertStruct.H:461
amrex::Real RandomReal(const amrex::Real min, const amrex::Real max)
Definition: ERF_TurbPertStruct.H:702
amrex::Vector< amrex::Real > tpi_Lpb
Definition: ERF_TurbPertStruct.H:685
amrex::Real tpi_net_buoyant
Definition: ERF_TurbPertStruct.H:689
Here is the call graph for this function:

◆ debug()

void TurbulentPerturbation::debug ( amrex::Real  )
inline
641  {
642  /*
643  amrex::PrintToFile("BoxPerturbationOutput") << "#################### PB output at time = "
644  << time << " ####################\n";
645  amrex::PrintToFile("BoxPerturbationOutput") << " Using type: " << pt_type << "\n";
646  amrex::PrintToFile("BoxPerturbationOutput") << " Net: " << tpi_net_buoyant << " Adjust : " << tpi_pert_adjust << "\n";
647  for (int i = 0; i < pb_mag.size(); i++) {
648  amrex::PrintToFile("BoxPerturbationOutput") << "[" << i
649  << "] pb_Umag=" << pb_mag[i]
650  << " | pb_interval=" << pb_interval[i]
651  << " (" << pb_local_etime[i]
652  << ") | pb_amp=" << pb_amp[i] << "\n";
653  }
654  amrex::PrintToFile("BoxPerturbationOutput") << "\n";
655  */
656  }

◆ init_tpi()

void TurbulentPerturbation::init_tpi ( const int  lev,
const amrex::Vector< amrex::BoxArray > &  subdomains_lev,
const amrex::GpuArray< amrex::Real, 3 >  dx,
const amrex::BoxArray &  ba,
const amrex::DistributionMapping &  dm,
const int  ngrow_state,
std::string  pp_prefix,
const amrex::Vector< amrex::IntVect >  refRatio,
const int  max_level 
)
inline
62  {
63  ref_ratio = refRatio;
64 
65  amrex::ParmParse pp(pp_prefix);
66 
67  // Reading inputs, and placing assertion for the perturbation inflow to work
68  pp.getarr("perturbation_box_dims",tpi_boxDim);
69  pp.getarr("perturbation_direction",tpi_direction);
70  pp.get("perturbation_layers",tpi_layers);
71  pp.get("perturbation_offset",tpi_offset);
72 
73  pp.get("perturbation_nondimensional",tpi_nonDim);
74 
75  tpi_Tinf = amrex::Real(300.);
76  pp.query("perturbation_T_infinity",tpi_Tinf);
77 
78  tpi_Ti = zero;
79  pp.query("perturbation_T_intensity",tpi_Ti);
80 
81  input_Ug = zero;
82  pp.query("perturbation_Ug",input_Ug);
83 
84  // Check variables message
85  if (tpi_offset < 0) { amrex::Abort("Please provide a valid inflow cell offset value for perturbation region (ie. 0-5)"); }
86  if (tpi_layers < 0) { amrex::Abort("Please provide a valid perturbation layer value (ie. 3-5)"); }
87  if (tpi_nonDim < zero) { amrex::Abort("Please provide a valid nondimensional number (ie. Ri = amrex::Real(0.042))"); }
88  for (int i = 0; i < tpi_boxDim.size(); i++) {
89  if (tpi_boxDim[i] < 3) { amrex::Abort("Please provide valid dimensions for perturbation boxes."); }
90  }
91  if (input_Ug < zero) { amrex::Abort("Please provide a valid geostrophic wind speed (ie. Ug = amrex::Real(10.0) m/s)"); }
92  if (tpi_Tinf < zero) { amrex::Abort("Please provide a valid ambient temperature value (ie. T_0 = T_infty)"); }
93  if (tpi_Ti < zero) { amrex::Abort("Please provide a valid temperature intensity value (ie. 0-one)"); }
94 
95  // Create a temporary box list to accumulate all the perturbation regions after box modification
96  amrex::BoxList tmp_bl;
97 
98  // boxSize for individual boxes
99  amrex::IntVect boxSize(tpi_boxDim[0],tpi_boxDim[1],tpi_boxDim[2]);
100 
101  if (tpi_direction[2] || tpi_direction[5]) { amrex::Abort("Currently not supporting z-direction flow perturbation"); }
102 
103  for (int isub = 0; isub < subdomains_lev.size(); ++isub) {
104  const amrex::BoxArray& subdomain = subdomains_lev[isub];
105  amrex::Box subdomain_box(subdomain.minimalBox());
106 
107  if (subdomain_box.numPts() != subdomain.numPts()) {
108  amrex::Abort("Turbulent perturbations require rectangular subdomains. "
109  "Level " + std::to_string(lev) +
110  ", subdomain " + std::to_string(isub) +
111  " is not a rectangular region fully covered by grids.");
112  }
113 
114  const amrex::IntVect& valid_box_lo = subdomain_box.smallEnd();
115  const amrex::IntVect& valid_box_hi = subdomain_box.bigEnd();
116 
117  // Creating perturbation regions and initializing with generic size.
118  amrex::Box lo_x_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
119  amrex::Box hi_x_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
120  amrex::Box lo_y_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
121  amrex::Box hi_y_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
122 
123  // Starting logic to set the size of the perturbation region(s)
124  //amrex::PrintToFile("BoxPerturbationOutput") << "Setting perturbation region in:";
125  // ***** X-direction perturbation *****
126  if (tpi_direction[0]) { // West
127  lo_x_bx.setSmall(amrex::IntVect(valid_box_lo[0]+tpi_offset, valid_box_lo[1]+tpi_direction[1]*tpi_offset, valid_box_lo[2]));
128  lo_x_bx.setBig (amrex::IntVect(valid_box_lo[0]+(tpi_layers*tpi_boxDim[0]-1)+tpi_offset, valid_box_hi[1]-(tpi_direction[4]*tpi_offset), valid_box_hi[2]));
129  amrex::PrintToFile("BoxPerturbationOutput") << " West face";
130  }
131 
132  if (tpi_direction[3]) { // East
133  hi_x_bx.setSmall(amrex::IntVect(valid_box_hi[0]-((tpi_layers*tpi_boxDim[0]-1)+tpi_offset), valid_box_lo[1]+tpi_direction[1]*tpi_offset, valid_box_lo[2]));
134  hi_x_bx.setBig(amrex::IntVect(valid_box_hi[0]-tpi_offset, valid_box_hi[1]-(tpi_direction[4]*tpi_offset), valid_box_hi[2]));
135  amrex::PrintToFile("BoxPerturbationOutput") << " East face";
136  }
137 
138  // ***** Y-direction Perturbation *****
139  if (tpi_direction[1]) { // North
140  lo_y_bx.setSmall(amrex::IntVect(valid_box_lo[0]+tpi_direction[0]*tpi_offset, valid_box_lo[1]+tpi_offset, valid_box_lo[2]));
141  lo_y_bx.setBig(amrex::IntVect(valid_box_hi[0]-tpi_direction[3]*tpi_offset, valid_box_lo[1]+((tpi_layers*tpi_boxDim[1])-1)+tpi_offset, valid_box_hi[2]));
142  amrex::PrintToFile("BoxPerturbationOutput") << " North face";
143  }
144 
145  if (tpi_direction[4]) { // South
146  hi_y_bx.setSmall(amrex::IntVect(valid_box_lo[0]+tpi_direction[0]*tpi_offset, valid_box_hi[1]-((tpi_layers*tpi_boxDim[1]-1)+tpi_offset), valid_box_lo[2]));
147  hi_y_bx.setBig(amrex::IntVect(valid_box_hi[0]-tpi_direction[3]*tpi_offset, valid_box_hi[1]-tpi_offset, valid_box_hi[2]));
148  amrex::PrintToFile("BoxPerturbationOutput") << " South face";
149  }
150 
151  // Performing box union for intersecting perturbation regions to avoid overlapping sections (double counting at corners)
152  if (tpi_direction[0] && tpi_direction[1]) { // Reshaping South smallEnd
153  amrex::Box lo_x_lo_y_u = lo_x_bx & lo_y_bx;
154  lo_y_bx.setSmall(amrex::IntVect(lo_x_lo_y_u.bigEnd(0)+1, lo_x_lo_y_u.smallEnd(1), lo_x_lo_y_u.smallEnd(2)));
155  }
156 
157  if (tpi_direction[3] && tpi_direction[1]) { // Reshaping South bigEnd
158  amrex::Box hi_x_lo_y_u = hi_x_bx & lo_y_bx;
159  lo_y_bx.setBig(amrex::IntVect(hi_x_lo_y_u.smallEnd(0)-1, hi_x_lo_y_u.bigEnd(1), hi_x_lo_y_u.bigEnd(2)));
160  }
161 
162  if (tpi_direction[0] && tpi_direction[4]) { // Reshaping North smallEnd
163  amrex::Box lo_x_hi_y_u = lo_x_bx & hi_y_bx;
164  hi_y_bx.setSmall(amrex::IntVect(lo_x_hi_y_u.bigEnd(0)+1, lo_x_hi_y_u.smallEnd(1), lo_x_hi_y_u.smallEnd(2)));
165  }
166 
167  if (tpi_direction[3] && tpi_direction[4]) { // Reshaping North bigEnd
168  amrex::Box hi_x_hi_y_u = hi_x_bx & hi_y_bx;
169  hi_y_bx.setBig(amrex::IntVect(hi_x_hi_y_u.smallEnd(0)-1, hi_x_hi_y_u.bigEnd(1), hi_x_hi_y_u.bigEnd(2)));
170  }
171 
172  // Creating structure box array for conserved quantity
173  if (tpi_direction[0]) { tmp_bl.push_back(lo_x_bx); }
174  if (tpi_direction[1]) { tmp_bl.push_back(lo_y_bx); }
175  if (tpi_direction[3]) { tmp_bl.push_back(hi_x_bx); }
176  if (tpi_direction[4]) { tmp_bl.push_back(hi_y_bx); }
177  }
178 
179  //amrex::PrintToFile("BoxPerturbationOutput") << "\nBoxList: " << tmp_bl << "\n";
180  amrex::BoxArray tmp_ba(tmp_bl);
181  tmp_ba.maxSize(boxSize);
182 
183  const int num_levels = max_level + 1;
184  if (pb_ba.size() < num_levels) {
185  pb_ba.resize(num_levels);
186  pb_mag.resize(num_levels);
187  pb_dir.resize(num_levels);
188  pb_netZero.resize(num_levels);
189  pb_interval.resize(num_levels);
190  pb_local_etime.resize(num_levels);
191  pb_amp.resize(num_levels);
192  pb_cell.resize(num_levels);
193  tpi_Lpb.resize(num_levels);
194  tpi_Wpb.resize(num_levels);
195  tpi_Hpb.resize(num_levels);
196  tpi_lref.resize(num_levels);
197  }
198 
199  pb_ba[lev] = tmp_ba;
200 
201  // Initializing mean magnitude and direction vectors
202  pb_mag[lev].resize(pb_ba[lev].size(), zero);
203  pb_dir[lev].resize(pb_ba[lev].size(), zero);
204  pb_netZero[lev].resize(pb_ba[lev].size(), zero);
205 
206  // Set size of vector and initialize
207  pb_interval[lev].resize(pb_ba[lev].size(), -one);
208  pb_local_etime[lev].resize(pb_ba[lev].size(), zero);
209  pb_amp[lev].resize(pb_ba[lev].size(), zero);
210 
211  // Creating data array for perturbation amplitude storage
212  pb_cell[lev].define(ba, dm, 1, ngrow_state); // this is the only place ba is used. Maybe we can print here to determine what's valid...
213  pb_cell[lev].setVal(0.);
214 
215  // Computing perturbation reference length
216  tpi_Lpb[lev] = tpi_boxDim[0]*dx[0];
217  tpi_Wpb[lev] = tpi_boxDim[1]*dx[1];
218  tpi_Hpb[lev] = tpi_boxDim[2]*dx[2];
219  tpi_lref[lev] = std::sqrt(tpi_Lpb[lev]*tpi_Lpb[lev] + tpi_Wpb[lev]*tpi_Wpb[lev]);
220 
223 
224  /*
225  // Function check point message
226  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_box_dims: "
227  << tpi_boxDim[0] << " "
228  << tpi_boxDim[1] << " "
229  << tpi_boxDim[2] << "\n";
230  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_direction: "
231  << tpi_direction[0] << " "
232  << tpi_direction[1] << " "
233  << tpi_direction[2] << "\n\n";
234  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_layers: " << tpi_layers << "\n";
235  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_offset: " << tpi_offset << "\n\n";
236  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_nondimensional: " << tpi_nonDim << "\n";
237  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_T_infinity: " << tpi_Tinf << "\n";
238  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_T_intensity: " << tpi_Ti << "\n";
239  amrex::PrintToFile("BoxPerturbationOutput") << "Reference length per box = " << tpi_lref[lev] << "\n\n";
240  amrex::PrintToFile("BoxPerturbationOutput") << "Turbulent perturbation BoxArray:\n" << pb_ba[lev] << "\n";
241  */
242  }
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
integer, private isub
Definition: ERF_module_mp_morr_two_moment.F90:164
amrex::Vector< amrex::MultiFab > pb_cell
Definition: ERF_TurbPertStruct.H:667
int tpi_offset
Definition: ERF_TurbPertStruct.H:673
amrex::Vector< int > tpi_boxDim
Definition: ERF_TurbPertStruct.H:675
amrex::Vector< int > tpi_direction
Definition: ERF_TurbPertStruct.H:676
Here is the call graph for this function:

◆ init_tpi_type()

void TurbulentPerturbation::init_tpi_type ( const int  lev,
const PerturbationType &  pert_type,
const int  max_level 
)
inline
31  {
32  if (pt_type.size() < max_level + 1) {
33  pt_type.resize(max_level + 1, -1);
34  }
35 
36  if (pert_type == PerturbationType::Source) {
37  pt_type[lev] = 0;
38  } else if (pert_type == PerturbationType::Direct) {
39  pt_type[lev] = 1;
40  } else if (pert_type == PerturbationType::CPM) {
41  pt_type[lev] = 2;
42  } else {
43  pt_type[lev] = -1;
44  }
45  }

◆ netZeroBuoyantAdd()

void TurbulentPerturbation::netZeroBuoyantAdd ( const int &  boxIdx,
const int &  lev 
)
inline
463  {
464  // Creating local copy of PB box array and magnitude
465  const amrex::BoxArray m_pb_ba = pb_ba[lev];
466  amrex::Real* m_pb_netZero = pb_netZero[lev].data();
467 
468  // Create device array for summation
469  amrex::Vector<amrex::Real> avg_h(1,zero);
470  amrex::Gpu::DeviceVector<amrex::Real> avg_d(1,zero);
471  amrex::Real* avg = avg_d.data();
472 
473  // Iterates through the cells of each box and sum the white noise perturbation
474  for (amrex::MFIter mfi(pb_cell[lev], TileNoZ()) ; mfi.isValid(); ++mfi) {
475  const amrex::Box& vbx = mfi.validbox();
476  amrex::Box pbx = amrex::convert(m_pb_ba[boxIdx], vbx.ixType());
477  amrex::Box ubx = pbx & vbx;
478  if (ubx.ok()) {
479  const amrex::Array4<const amrex::Real>& pert_cell = pb_cell[lev].const_array(mfi);
480  amrex::Real norm = one / static_cast<amrex::Real>(ubx.numPts());
481  ParallelFor(amrex::Gpu::KernelInfo().setReduction(true), ubx, [=]
482  AMREX_GPU_DEVICE(int i, int j, int k, amrex::Gpu::Handler const& handler) noexcept {
483  amrex::Gpu::deviceReduceSum(&avg[0], pert_cell(i,j,k)*norm, handler);
484  });
485  amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(), avg_h.begin());
486 
487  // Assigning onto storage array
488  m_pb_netZero[boxIdx] = avg_h[0];
489  }
490  }
491  }

Referenced by calc_tpi_update().

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

◆ netZeroBuoyantAdjust()

void TurbulentPerturbation::netZeroBuoyantAdjust ( const int &  boxIdx,
const int &  lev 
)
inline
497  {
498  // Creating local copy of PB box array and magnitude
499  const amrex::BoxArray m_pb_ba = pb_ba[lev];
500  for (amrex::MFIter mfi(pb_cell[lev], TileNoZ()) ; mfi.isValid(); ++mfi) {
501  const amrex::Box& vbx = mfi.validbox();
502  amrex::Box pbx = amrex::convert(m_pb_ba[boxIdx], vbx.ixType());
503  amrex::Box ubx = pbx & vbx;
504  if (ubx.ok()) {
505  const amrex::Real adjust = tpi_pert_adjust;
506  const amrex::Array4<amrex::Real>& pert_cell = pb_cell[lev].array(mfi);
507  ParallelFor(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
508  pert_cell(i,j,k) -= adjust;
509  });
510  }
511  }
512  }

Referenced by calc_tpi_update().

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

◆ pseudoRandomPert()

void TurbulentPerturbation::pseudoRandomPert ( const int &  boxIdx,
const int &  lev,
const amrex::IndexType &  m_ixtype 
)
inline
416  {
417  for (amrex::MFIter mfi(pb_cell[lev],TileNoZ()); mfi.isValid(); ++mfi) {
418  amrex::Box vbx = mfi.validbox();
419  amrex::Box pbx = amrex::convert(pb_ba[lev][boxIdx], m_ixtype);
420  amrex::Box ubx = pbx & vbx;
421  if (ubx.ok()) {
422  amrex::Real amp_copy = pb_amp[lev][boxIdx];
423  const amrex::Array4<amrex::Real>& pert_cell = pb_cell[lev].array(mfi);
424 
425  if (pt_type[lev] == 2) { // CPM
426  amrex::Real rand_number_const = RandomReal(-one, one);
427  ParallelFor(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
428  pert_cell(i,j,k) = rand_number_const * amp_copy;
429  });
430  } else {
431  ParallelForRNG(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k, const amrex::RandomEngine& engine) noexcept {
432  amrex::Real rand_double = amrex::Random(engine);
433  pert_cell(i,j,k) = (rand_double*two - one) * amp_copy;
434  });
435  }
436 
437  }
438  }
439  }
constexpr amrex::Real two
Definition: ERF_Constants.H:8
ParallelForRNG(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k, const amrex::RandomEngine &engine) noexcept { const Real x=prob_lo_x+(i+myhalf) *dx;const Real y=prob_lo_y+(j+myhalf) *dy;const Real z=z_cc(i, j, k);const Real r=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc)+(z-zc) *(z-zc));if((z<=pert_ref_height) &&(T_0_Pert_Mag !=amrex::Real(0))) { Real rand_double=amrex::Random(engine);state_pert(i, j, k, RhoTheta_comp)=(rand_double *amrex::Real(2) - amrex::Real(1)) *T_0_Pert_Mag;if(!pert_rhotheta) { state_pert(i, j, k, RhoTheta_comp) *=r_hse(i, j, k);} } state_pert(i, j, k, RhoScalar_comp)=A_0 *std::exp(-amrex::Real(10.) *r *r);if(state_pert.nComp() > RhoKE_comp) { if(rhoKE_0 > 0) { state_pert(i, j, k, RhoKE_comp)=rhoKE_0;} else { state_pert(i, j, k, RhoKE_comp)=r_hse(i, j, k) *KE_0;} if(KE_decay_height > 0) { state_pert(i, j, k, RhoKE_comp) *=amrex::max(std::pow(1 - amrex::min(z/KE_decay_height, amrex::Real(1)), KE_decay_order), Real(1e-12));} } })

Referenced by calc_tpi_update().

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

◆ RandomReal()

amrex::Real TurbulentPerturbation::RandomReal ( const amrex::Real  min,
const amrex::Real  max 
)
inlineprivate
703  {
704  amrex::Real r = (amrex::Real) rand() / (amrex::Real) RAND_MAX;
705  return min + r * (max - min);
706  }

Referenced by calc_tpi_update(), and pseudoRandomPert().

Here is the caller graph for this function:

◆ zero_amp()

void TurbulentPerturbation::zero_amp ( const int &  boxIdx,
const int &  lev,
const amrex::IndexType &  m_ixtype 
)
inline
445  {
446 
447  for (amrex::MFIter mfi(pb_cell[lev],TileNoZ()); mfi.isValid(); ++mfi) {
448  amrex::Box vbx = mfi.validbox();
449  amrex::Box pbx = amrex::convert(pb_ba[lev][boxIdx], m_ixtype);
450  amrex::Box ubx = pbx & vbx;
451  if (ubx.ok()) {
452  const amrex::Array4<amrex::Real>& pert_cell = pb_cell[lev].array(mfi);
453  ParallelFor(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
454  pert_cell(i,j,k) = zero;
455  });
456  }
457  }
458  }

Referenced by calc_tpi_update().

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

Member Data Documentation

◆ input_Ug

amrex::Real TurbulentPerturbation::input_Ug
private

Referenced by calc_tpi_amp(), and init_tpi().

◆ pb_amp

amrex::Vector<amrex::Vector<amrex::Real> > TurbulentPerturbation::pb_amp
private

◆ pb_ba

amrex::Vector<amrex::BoxArray> TurbulentPerturbation::pb_ba

◆ pb_cell

amrex::Vector<amrex::MultiFab> TurbulentPerturbation::pb_cell

◆ pb_dir

amrex::Vector<amrex::Vector<amrex::Real> > TurbulentPerturbation::pb_dir

◆ pb_interval

amrex::Vector<amrex::Vector<amrex::Real> > TurbulentPerturbation::pb_interval
private

Referenced by calc_tpi_update(), and init_tpi().

◆ pb_local_etime

amrex::Vector<amrex::Vector<amrex::Real> > TurbulentPerturbation::pb_local_etime
private

Referenced by calc_tpi_update(), and init_tpi().

◆ pb_mag

amrex::Vector<amrex::Vector<amrex::Real> > TurbulentPerturbation::pb_mag

◆ pb_netZero

amrex::Vector<amrex::Vector<amrex::Real> > TurbulentPerturbation::pb_netZero
private

◆ pt_type

amrex::Vector<int> TurbulentPerturbation::pt_type

◆ ref_ratio

amrex::Vector<amrex::IntVect> TurbulentPerturbation::ref_ratio
private

Referenced by calc_tpi_amp(), and init_tpi().

◆ tpi_boxDim

amrex::Vector<int> TurbulentPerturbation::tpi_boxDim
private

Referenced by init_tpi().

◆ tpi_direction

amrex::Vector<int> TurbulentPerturbation::tpi_direction
private

Referenced by init_tpi().

◆ tpi_Hpb

amrex::Vector<amrex::Real> TurbulentPerturbation::tpi_Hpb
private

Referenced by calc_tpi_amp(), and init_tpi().

◆ tpi_layers

int TurbulentPerturbation::tpi_layers
private

Referenced by calc_tpi_update(), and init_tpi().

◆ tpi_Lpb

amrex::Vector<amrex::Real> TurbulentPerturbation::tpi_Lpb
private

Referenced by calc_tpi_update(), and init_tpi().

◆ tpi_lref

amrex::Vector<amrex::Real> TurbulentPerturbation::tpi_lref
private

Referenced by calc_tpi_update(), and init_tpi().

◆ tpi_net_buoyant

amrex::Real TurbulentPerturbation::tpi_net_buoyant
private

Referenced by calc_tpi_update(), and init_tpi().

◆ tpi_nonDim

amrex::Real TurbulentPerturbation::tpi_nonDim
private

Referenced by calc_tpi_amp(), and init_tpi().

◆ tpi_offset

int TurbulentPerturbation::tpi_offset
private

Referenced by init_tpi().

◆ tpi_pert_adjust

amrex::Real TurbulentPerturbation::tpi_pert_adjust
private

◆ tpi_Ti

amrex::Real TurbulentPerturbation::tpi_Ti
private

Referenced by calc_tpi_amp(), and init_tpi().

◆ tpi_Tinf

amrex::Real TurbulentPerturbation::tpi_Tinf
private

Referenced by calc_tpi_amp(), and init_tpi().

◆ tpi_Wpb

amrex::Vector<amrex::Real> TurbulentPerturbation::tpi_Wpb
private

Referenced by calc_tpi_update(), and init_tpi().


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