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
355  {
356  for (int boxIdx = 0; boxIdx < pb_ba[lev].size(); boxIdx++) {
357  amrex::Box pbx = amrex::convert(pb_ba[lev][boxIdx], m_ixtype);
358  amrex::Box ubx = pbx & vbx;
359  if (ubx.ok()) {
360  ParallelFor(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
361  src_arr(i,j,k,comp) += pert_cell(i,j,k);
362 
363  // For box region debug only
364  #ifdef INDEX_PERTURB
365  src_arr(i,j,k,comp) = (amrex::Real) (boxIdx + amrex::Real(5.));
366  #endif
367  });
368  }
369  }
370  }
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+0.5) *dx[0];const Real y=(j+0.5) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - 1.0)/(2.0 *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(1.0+deltaT, inv_gm1);const Real T=(1.0+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)=0.25 *(1.0+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:659

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
376  {
377  pb_amp[lev][boxIdx] = zero; // Safety step
378  if (pt_type[lev] == 2) { // CPM
379  amrex::Real cp = 1004; // specific heat of air [J/(kg K)]
380  amrex::Real Ec = amrex::Real(0.2); // Eckert number
381  pb_amp[lev][boxIdx] = (input_Ug * input_Ug) / (Ec * cp);
382  } else { // box perturbation
383  amrex::Real Um = pb_mag[lev][boxIdx];
384  amrex::Real beta = one/tpi_Tinf; // Thermal expansion coefficient
385 
386  // Pseudo Random temperature the ignores scale when mechanically tripping turbulence
388  // get total refinement ratio on each level
389  int total_ref_ratio = 1;
390  for (int level = lev; level >= 1; level--) {
391  total_ref_ratio *= ref_ratio[level-1][2];
392  }
393  // calculation needs to be scale-aware since the formulation relies on the physical size of the box
394  if (tpi_Ti > zero) g = (tpi_nonDim * Um * Um) / (tpi_Ti * tpi_Hpb[lev]) * 1 / total_ref_ratio;
395 
396  // Ma and Senocak (2023) Eq. 8, solving for delta phi
397  pb_amp[lev][boxIdx] = (tpi_nonDim * Um * Um) / (g * beta * tpi_Hpb[lev]) * 1 / total_ref_ratio;
398 
399  if (pt_type[lev] == 0) {
400  // Performing this step converts the perturbation proportionality into
401  // the forcing term
402  // Ma & Senocak (2023) Eq. 7
403  pb_amp[lev][boxIdx] /= interval;
404  }
405  }
406  }
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:679
amrex::Vector< amrex::IntVect > ref_ratio
Definition: ERF_TurbPertStruct.H:690
amrex::Real input_Ug
Definition: ERF_TurbPertStruct.H:691
amrex::Vector< amrex::Vector< amrex::Real > > pb_mag
Definition: ERF_TurbPertStruct.H:660
amrex::Vector< int > pt_type
Definition: ERF_TurbPertStruct.H:656
amrex::Vector< amrex::Real > tpi_Hpb
Definition: ERF_TurbPertStruct.H:682
amrex::Real tpi_Ti
Definition: ERF_TurbPertStruct.H:678
amrex::Real tpi_nonDim
Definition: ERF_TurbPertStruct.H:677
amrex::Vector< amrex::Vector< amrex::Real > > pb_amp
Definition: ERF_TurbPertStruct.H:696

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

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

◆ debug()

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

◆ 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  // Initialization for some 0 dependent terms
64  tpi_Ti = zero;
65  tpi_Tinf = amrex::Real(300.);
66 
67  ref_ratio = refRatio;
68 
69  amrex::ParmParse pp(pp_prefix);
70 
71  // Reading inputs, and placing assertion for the perturbation inflow to work
72  pp.queryarr("perturbation_box_dims",tpi_boxDim);
73  pp.queryarr("perturbation_direction",tpi_direction);
74  pp.query("perturbation_layers",tpi_layers);
75  pp.query("perturbation_offset",tpi_offset);
76 
77  pp.query("perturbation_nondimensional",tpi_nonDim);
78  pp.query("perturbation_T_infinity",tpi_Tinf);
79  pp.query("perturbation_T_intensity",tpi_Ti);
80  pp.query("perturbation_Ug",input_Ug);
81 
82  // Check variables message
83  if (tpi_layers < 0) { amrex::Abort("Please provide a valid perturbation layer value (ie. 3-5)"); }
84  if (tpi_offset < 0) { amrex::Abort("Please provide a valid inflow cell offset value for perturbation region (ie. 0-5)"); }
85  for (int i = 0; i < tpi_boxDim.size(); i++) {
86  if (tpi_boxDim[i] < 3) { amrex::Abort("Please provide valid dimensions for perturbation boxes."); }
87  }
88  if (tpi_nonDim < zero) { amrex::Abort("Please provide a valid nondimensional number (ie. Ri = amrex::Real(0.042))"); }
89  if (input_Ug < zero) { amrex::Abort("Please provide a valid geostrophic wind speed (ie. Ug = amrex::Real(10.0) m/s)"); }
90  if (tpi_Tinf < zero) { amrex::Abort("Please provide a valid ambient temperature value (ie. T_0 = T_infty)"); }
91  if (tpi_Ti < zero) { amrex::Abort("Please provide a valid temperature intensity value (ie. 0-one)"); }
92 
93  // Create a temporary box list to accumulate all the perturbation regions after box modification
94  amrex::BoxList tmp_bl;
95 
96  // boxSize for individual boxes
97  amrex::IntVect boxSize(tpi_boxDim[0],tpi_boxDim[1],tpi_boxDim[2]);
98 
99  if (tpi_direction[2] || tpi_direction[5]) { amrex::Abort("Currently not supporting z-direction flow perturbation"); }
100 
101  for (int isub = 0; isub < subdomains_lev.size(); ++isub) {
102  const amrex::BoxArray& subdomain = subdomains_lev[isub];
103  amrex::Box subdomain_box(subdomain.minimalBox());
104 
105  if (subdomain_box.numPts() != subdomain.numPts()) {
106  amrex::Abort("Turbulent perturbations require rectangular subdomains. "
107  "Level " + std::to_string(lev) +
108  ", subdomain " + std::to_string(isub) +
109  " is not a rectangular region fully covered by grids.");
110  }
111 
112  const amrex::IntVect& valid_box_lo = subdomain_box.smallEnd();
113  const amrex::IntVect& valid_box_hi = subdomain_box.bigEnd();
114 
115  // Creating perturbation regions and initializing with generic size.
116  amrex::Box lo_x_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
117  amrex::Box hi_x_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
118  amrex::Box lo_y_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
119  amrex::Box hi_y_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
120 
121  // Starting logic to set the size of the perturbation region(s)
122  //amrex::PrintToFile("BoxPerturbationOutput") << "Setting perturbation region in:";
123  // ***** X-direction perturbation *****
124  if (tpi_direction[0]) { // West
125  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]));
126  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]));
127  amrex::PrintToFile("BoxPerturbationOutput") << " West face";
128  }
129 
130  if (tpi_direction[3]) { // East
131  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]));
132  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]));
133  amrex::PrintToFile("BoxPerturbationOutput") << " East face";
134  }
135 
136  // ***** Y-direction Perturbation *****
137  if (tpi_direction[1]) { // North
138  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]));
139  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]));
140  amrex::PrintToFile("BoxPerturbationOutput") << " North face";
141  }
142 
143  if (tpi_direction[4]) { // South
144  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]));
145  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]));
146  amrex::PrintToFile("BoxPerturbationOutput") << " South face";
147  }
148 
149  // Performing box union for intersecting perturbation regions to avoid overlapping sections (double counting at corners)
150  if (tpi_direction[0] && tpi_direction[1]) { // Reshaping South smallEnd
151  amrex::Box lo_x_lo_y_u = lo_x_bx & lo_y_bx;
152  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)));
153  }
154 
155  if (tpi_direction[3] && tpi_direction[1]) { // Reshaping South bigEnd
156  amrex::Box hi_x_lo_y_u = hi_x_bx & lo_y_bx;
157  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)));
158  }
159 
160  if (tpi_direction[0] && tpi_direction[4]) { // Reshaping North smallEnd
161  amrex::Box lo_x_hi_y_u = lo_x_bx & hi_y_bx;
162  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)));
163  }
164 
165  if (tpi_direction[3] && tpi_direction[4]) { // Reshaping North bigEnd
166  amrex::Box hi_x_hi_y_u = hi_x_bx & hi_y_bx;
167  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)));
168  }
169 
170  // Creating structure box array for conserved quantity
171  if (tpi_direction[0]) { tmp_bl.push_back(lo_x_bx); }
172  if (tpi_direction[1]) { tmp_bl.push_back(lo_y_bx); }
173  if (tpi_direction[3]) { tmp_bl.push_back(hi_x_bx); }
174  if (tpi_direction[4]) { tmp_bl.push_back(hi_y_bx); }
175  }
176 
177  //amrex::PrintToFile("BoxPerturbationOutput") << "\nBoxList: " << tmp_bl << "\n";
178  amrex::BoxArray tmp_ba(tmp_bl);
179  tmp_ba.maxSize(boxSize);
180 
181  const int num_levels = max_level + 1;
182  if (pb_ba.size() < num_levels) {
183  pb_ba.resize(num_levels);
184  pb_mag.resize(num_levels);
185  pb_dir.resize(num_levels);
186  pb_netZero.resize(num_levels);
187  pb_interval.resize(num_levels);
188  pb_local_etime.resize(num_levels);
189  pb_amp.resize(num_levels);
190  pb_cell.resize(num_levels);
191  tpi_Lpb.resize(num_levels);
192  tpi_Wpb.resize(num_levels);
193  tpi_Hpb.resize(num_levels);
194  tpi_lref.resize(num_levels);
195  }
196 
197  pb_ba[lev] = tmp_ba;
198 
199  // Initializing mean magnitude and direction vectors
200  pb_mag[lev].resize(pb_ba[lev].size(), zero);
201  pb_dir[lev].resize(pb_ba[lev].size(), zero);
202  pb_netZero[lev].resize(pb_ba[lev].size(), zero);
203 
204  // Set size of vector and initialize
205  pb_interval[lev].resize(pb_ba[lev].size(), -one);
206  pb_local_etime[lev].resize(pb_ba[lev].size(), zero);
207  pb_amp[lev].resize(pb_ba[lev].size(), zero);
208 
209  // Creating data array for perturbation amplitude storage
210  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...
211  pb_cell[lev].setVal(0.);
212 
213  // Computing perturbation reference length
214  tpi_Lpb[lev] = tpi_boxDim[0]*dx[0];
215  tpi_Wpb[lev] = tpi_boxDim[1]*dx[1];
216  tpi_Hpb[lev] = tpi_boxDim[2]*dx[2];
217  tpi_lref[lev] = sqrt(tpi_Lpb[lev]*tpi_Lpb[lev] + tpi_Wpb[lev]*tpi_Wpb[lev]);
218 
221 
222  /*
223  // Function check point message
224  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_box_dims: "
225  << tpi_boxDim[0] << " "
226  << tpi_boxDim[1] << " "
227  << tpi_boxDim[2] << "\n";
228  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_direction: "
229  << tpi_direction[0] << " "
230  << tpi_direction[1] << " "
231  << tpi_direction[2] << "\n\n";
232  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_layers: " << tpi_layers << "\n";
233  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_offset: " << tpi_offset << "\n\n";
234  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_nondimensional: " << tpi_nonDim << "\n";
235  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_T_infinity: " << tpi_Tinf << "\n";
236  amrex::PrintToFile("BoxPerturbationOutput") << "perturbation_T_intensity: " << tpi_Ti << "\n";
237  amrex::PrintToFile("BoxPerturbationOutput") << "Reference length per box = " << tpi_lref[lev] << "\n\n";
238  amrex::PrintToFile("BoxPerturbationOutput") << "Turbulent perturbation BoxArray:\n" << pb_ba[lev] << "\n";
239  */
240  }
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:665
int tpi_offset
Definition: ERF_TurbPertStruct.H:671
amrex::Vector< int > tpi_boxDim
Definition: ERF_TurbPertStruct.H:673
amrex::Vector< int > tpi_direction
Definition: ERF_TurbPertStruct.H:674
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
461  {
462  // Creating local copy of PB box array and magnitude
463  const amrex::BoxArray m_pb_ba = pb_ba[lev];
464  amrex::Real* m_pb_netZero = pb_netZero[lev].data();
465 
466  // Create device array for summation
467  amrex::Vector<amrex::Real> avg_h(1,zero);
468  amrex::Gpu::DeviceVector<amrex::Real> avg_d(1,zero);
469  amrex::Real* avg = avg_d.data();
470 
471  // Iterates through the cells of each box and sum the white noise perturbation
472  for (amrex::MFIter mfi(pb_cell[lev], TileNoZ()) ; mfi.isValid(); ++mfi) {
473  const amrex::Box& vbx = mfi.validbox();
474  amrex::Box pbx = amrex::convert(m_pb_ba[boxIdx], vbx.ixType());
475  amrex::Box ubx = pbx & vbx;
476  if (ubx.ok()) {
477  const amrex::Array4<const amrex::Real>& pert_cell = pb_cell[lev].const_array(mfi);
478  amrex::Real norm = one / static_cast<amrex::Real>(ubx.numPts());
479  ParallelFor(amrex::Gpu::KernelInfo().setReduction(true), ubx, [=]
480  AMREX_GPU_DEVICE(int i, int j, int k, amrex::Gpu::Handler const& handler) noexcept {
481  amrex::Gpu::deviceReduceSum(&avg[0], pert_cell(i,j,k)*norm, handler);
482  });
483  amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(), avg_h.begin());
484 
485  // Assigning onto storage array
486  m_pb_netZero[boxIdx] = avg_h[0];
487  }
488  }
489  }

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

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
414  {
415  for (amrex::MFIter mfi(pb_cell[lev],TileNoZ()); mfi.isValid(); ++mfi) {
416  amrex::Box vbx = mfi.validbox();
417  amrex::Box pbx = amrex::convert(pb_ba[lev][boxIdx], m_ixtype);
418  amrex::Box ubx = pbx & vbx;
419  if (ubx.ok()) {
420  amrex::Real amp_copy = pb_amp[lev][boxIdx];
421  const amrex::Array4<amrex::Real>& pert_cell = pb_cell[lev].array(mfi);
422 
423  if (pt_type[lev] == 2) { // CPM
424  amrex::Real rand_number_const = RandomReal(-one, one);
425  ParallelFor(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
426  pert_cell(i,j,k) = rand_number_const * amp_copy;
427  });
428  } else {
429  ParallelForRNG(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k, const amrex::RandomEngine& engine) noexcept {
430  amrex::Real rand_double = amrex::Random(engine);
431  pert_cell(i,j,k) = (rand_double*two - one) * amp_copy;
432  });
433  }
434 
435  }
436  }
437  }
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 *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
701  {
702  amrex::Real r = (amrex::Real) rand() / (amrex::Real) RAND_MAX;
703  return min + r * (max - min);
704  }

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
443  {
444 
445  for (amrex::MFIter mfi(pb_cell[lev],TileNoZ()); mfi.isValid(); ++mfi) {
446  amrex::Box vbx = mfi.validbox();
447  amrex::Box pbx = amrex::convert(pb_ba[lev][boxIdx], m_ixtype);
448  amrex::Box ubx = pbx & vbx;
449  if (ubx.ok()) {
450  const amrex::Array4<amrex::Real>& pert_cell = pb_cell[lev].array(mfi);
451  ParallelFor(ubx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
452  pert_cell(i,j,k) = zero;
453  });
454  }
455  }
456  }

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: