1 #ifndef ERF_TURB_PERT_STRUCT_H_
2 #define ERF_TURB_PERT_STRUCT_H_
5 #include <AMReX_MultiFabUtil.H>
30 const amrex::IntVect& nx,
31 const amrex::GpuArray<amrex::Real,3> dx,
32 const amrex::BoxArray& ba,
33 const amrex::DistributionMapping& dm,
34 const int ngrow_state,
35 std::string pp_prefix)
42 amrex::ParmParse
pp(pp_prefix);
52 pp.query(
"perturbation_T_intensity",
tpi_Ti);
55 if (
tpi_layers < 0) { amrex::Abort(
"Please provide a valid perturbation layer value (ie. 3-5)"); }
56 if (
tpi_offset < 0) { amrex::Abort(
"Please provide a valid inflow cell offset value for perturbation region (ie. 0-5)"); }
58 if (
tpi_boxDim[i] == 0) { amrex::Abort(
"Please provide valid dimensions for perturbation boxes."); }
60 if (
tpi_nonDim < 0.) { amrex::Abort(
"Please provide a valid nondimensional number (ie. Ri = 0.042)"); }
61 if (
tpi_Tinf < 0.) { amrex::Abort(
"Please provide a valid ambient temperature value (ie. T_0 = T_infty)"); }
62 if (
tpi_Ti < 0.) { amrex::Abort(
"Please provide a valid temperature intensity value (ie. 0-1.0)"); }
65 amrex::Box lo_x_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
66 amrex::Box hi_x_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
67 amrex::Box lo_y_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
68 amrex::Box hi_y_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
71 amrex::BoxList tmp_bl;
108 amrex::Box lo_x_lo_y_u = lo_x_bx & lo_y_bx;
109 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)));
113 amrex::Box hi_x_lo_y_u = hi_x_bx & lo_y_bx;
114 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)));
118 amrex::Box lo_x_hi_y_u = lo_x_bx & hi_y_bx;
119 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)));
123 amrex::Box hi_x_hi_y_u = hi_x_bx & hi_y_bx;
124 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)));
133 amrex::BoxArray tmp_ba(tmp_bl);
134 tmp_ba.maxSize(boxSize);
135 pb_ba.push_back(tmp_ba);
147 pb_cell.define(ba, dm, 1, ngrow_state);
186 const amrex::Real dt,
187 amrex::MultiFab& mf_xvel,
188 amrex::MultiFab& mf_yvel,
189 amrex::MultiFab& mf_cons)
195 srand( (
unsigned) time(NULL) );
197 auto m_ixtype = mf_cons.boxArray().ixType();
199 for (
int boxIdx = 0; boxIdx <
pb_ba[lev].size(); boxIdx++) {
209 if (
pb_mag[boxIdx] !=0.) {
238 for (
int boxIdx = 0; boxIdx <
pb_ba[lev].size(); boxIdx++) {
247 const amrex::Box& vbx,
249 const amrex::IndexType& m_ixtype,
250 const amrex::Array4<amrex::Real>& src_arr,
251 const amrex::Array4<amrex::Real const>& pert_cell)
253 for (
int boxIdx = 0; boxIdx <
pb_ba[lev].size(); boxIdx++) {
254 amrex::Box pbx = amrex::convert(
pb_ba[lev][boxIdx], m_ixtype);
255 amrex::Box ubx = pbx & vbx;
257 ParallelFor(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
258 src_arr(i,j,k,comp) += pert_cell(i,j,k);
262 src_arr(i,j,k,comp) = (amrex::Real) (boxIdx + 5.);
271 const amrex::Real& interval)
273 amrex::Real Um =
pb_mag[boxIdx];
289 pb_amp[boxIdx] /= interval;
298 const amrex::IndexType& m_ixtype)
301 int fix_random_seed = 0;
302 amrex::ParmParse
pp(
"erf");
pp.query(
"fix_random_seed", fix_random_seed);
303 if (fix_random_seed) {
304 amrex::InitRandom(1024UL);
307 for (amrex::MFIter mfi(
pb_cell,
TileNoZ()); mfi.isValid(); ++mfi) {
308 amrex::Box vbx = mfi.validbox();
309 amrex::Box pbx = amrex::convert(
pb_ba[lev][boxIdx], m_ixtype);
310 amrex::Box ubx = pbx & vbx;
312 amrex::Real amp_copy =
pb_amp[boxIdx];
313 const amrex::Array4<amrex::Real>& pert_cell =
pb_cell.array(mfi);
314 ParallelForRNG(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k,
const amrex::RandomEngine& engine) noexcept {
315 amrex::Real rand_double = amrex::Random(engine);
316 pert_cell(i,j,k) = (rand_double*2.0 - 1.0) * amp_copy;
327 const amrex::BoxArray m_pb_ba =
pb_ba[lev];
331 amrex::Vector<amrex::Real>
avg_h(1,0.);
332 amrex::Gpu::DeviceVector<amrex::Real> avg_d(1,0.);
333 amrex::Real* avg = avg_d.data();
336 for (amrex::MFIter mfi(
pb_cell,
TileNoZ()) ; mfi.isValid(); ++mfi) {
337 const amrex::Box& vbx = mfi.validbox();
338 amrex::Box pbx = amrex::convert(m_pb_ba[boxIdx], vbx.ixType());
339 amrex::Box ubx = pbx & vbx;
341 const amrex::Array4<const amrex::Real>& pert_cell =
pb_cell.const_array(mfi);
342 int npts = ubx.numPts();
343 amrex::Real norm = 1.0 / (amrex::Real) npts;
344 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubx, [=]
345 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
346 amrex::Gpu::deviceReduceSum(&avg[0], pert_cell(i,j,k)*norm, handler);
348 amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(),
avg_h.begin());
351 m_pb_netZero[boxIdx] =
avg_h[0];
362 const amrex::BoxArray m_pb_ba =
pb_ba[lev];
363 for (amrex::MFIter mfi(
pb_cell,
TileNoZ()) ; mfi.isValid(); ++mfi) {
364 const amrex::Box& vbx = mfi.validbox();
365 amrex::Box pbx = amrex::convert(m_pb_ba[boxIdx], vbx.ixType());
366 amrex::Box ubx = pbx & vbx;
369 const amrex::Array4<amrex::Real>& pert_cell =
pb_cell.array(mfi);
370 ParallelFor(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
371 pert_cell(i,j,k) -= adjust;
378 #define USE_VOLUME_AVERAGE
383 amrex::MultiFab& mf_cons,
384 amrex::MultiFab& mf_xvel,
385 amrex::MultiFab& mf_yvel)
389 const amrex::BoxArray m_pb_ba =
pb_ba[lev];
391 m_pb_mag[boxIdx] = 0.;
397 amrex::Vector<amrex::Real>
avg_h(n_avg,0.);
398 amrex::Gpu::DeviceVector<amrex::Real> avg_d(n_avg,0.);
399 amrex::Real* avg = avg_d.data();
402 for (amrex::MFIter mfi(mf_cons,
TileNoZ()); mfi.isValid(); ++mfi) {
405 const amrex::Box& vbx = mfi.validbox();
408 auto ixtype_u = mf_xvel.boxArray().ixType();
409 amrex::Box vbx_u = amrex::convert(vbx,ixtype_u);
410 amrex::Box pbx_u = amrex::convert(m_pb_ba[boxIdx], ixtype_u);
411 amrex::Box ubx_u = pbx_u & vbx_u;
414 auto ixtype_v = mf_yvel.boxArray().ixType();
415 amrex::Box vbx_v = amrex::convert(vbx,ixtype_v);
416 amrex::Box pbx_v = amrex::convert(m_pb_ba[boxIdx], ixtype_v);
417 amrex::Box ubx_v = pbx_v & vbx_v;
421 const amrex::Array4<const amrex::Real>& xvel_arry = mf_xvel.const_array(mfi);
423 #ifdef USE_VOLUME_AVERAGE
424 int npts = ubx_u.numPts();
425 amrex::Real norm = 1.0 / (amrex::Real) npts;
426 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubx_u, [=]
427 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
428 amrex::Gpu::deviceReduceSum(&avg[0], xvel_arry(i,j,k)*norm, handler);
432 #ifdef USE_SLAB_AVERAGE
433 amrex::Box ubxSlab_lo = makeSlab(ubx_u,2,ubx_u.smallEnd(2));
434 amrex::Box ubxSlab_hi = makeSlab(ubx_u,2,ubx_u.bigEnd(2));
435 int npts_lo = ubxSlab_lo.numPts();
436 int npts_hi = ubxSlab_hi.numPts();
437 amrex::Real norm_lo = 1.0 / (amrex::Real) npts_lo;
438 amrex::Real norm_hi = 1.0 / (amrex::Real) npts_hi;
441 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubxSlab_lo, [=]
442 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
443 amrex::Gpu::deviceReduceSum(&avg[0], xvel_arry(i,j,k)*norm_lo, handler);
447 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubxSlab_hi, [=]
448 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
449 amrex::Gpu::deviceReduceSum(&avg[2], xvel_arry(i,j,k)*norm_hi, handler);
456 const amrex::Array4<const amrex::Real>& yvel_arry = mf_yvel.const_array(mfi);
458 #ifdef USE_VOLUME_AVERAGE
459 int npts = ubx_v.numPts();
460 amrex::Real norm = 1.0 / (amrex::Real) npts;
461 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubx_v, [=]
462 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
463 amrex::Gpu::deviceReduceSum(&avg[1], yvel_arry(i,j,k)*norm, handler);
467 #ifdef USE_SLAB_AVERAGE
468 amrex::Box ubxSlab_lo = makeSlab(ubx_v,2,ubx_v.smallEnd(2));
469 amrex::Box ubxSlab_hi = makeSlab(ubx_v,2,ubx_v.bigEnd(2));
470 int npts_lo = ubxSlab_lo.numPts();
471 int npts_hi = ubxSlab_hi.numPts();
472 amrex::Real norm_lo = 1.0 / (amrex::Real) npts_lo;
473 amrex::Real norm_hi = 1.0 / (amrex::Real) npts_hi;
476 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubxSlab_lo, [=]
477 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
478 amrex::Gpu::deviceReduceSum(&avg[1], yvel_arry(i,j,k)*norm_lo, handler);
482 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubxSlab_hi, [=]
483 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
484 amrex::Gpu::deviceReduceSum(&avg[3], yvel_arry(i,j,k)*norm_hi, handler);
491 amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(),
avg_h.begin());
494 #ifdef USE_VOLUME_AVERAGE
498 #ifdef USE_SLAB_AVERAGE
530 amrex::Vector<amrex::BoxArray>
pb_ba;
567 amrex::Real
RandomReal (
const amrex::Real min,
const amrex::Real max)
569 amrex::Real r = (amrex::Real) rand() / (amrex::Real) RAND_MAX;
570 return min + r * (max - min);
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:219
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real avg_h(amrex::Real XXXm, amrex::Real XXXp)
Definition: ERF_EBUtils.H:10
Definition: ERF_TurbPertStruct.H:18
int tpi_layers
Definition: ERF_TurbPertStruct.H:540
void pseudoRandomPert(const int &boxIdx, const int &lev, const amrex::IndexType &m_ixtype)
Definition: ERF_TurbPertStruct.H:296
amrex::Real tpi_Tinf
Definition: ERF_TurbPertStruct.H:549
amrex::Real tpi_Wpb
Definition: ERF_TurbPertStruct.H:554
amrex::Real tpi_Lpb
Definition: ERF_TurbPertStruct.H:553
amrex::Real * get_pb_mag()
Definition: ERF_TurbPertStruct.H:505
amrex::Vector< amrex::Real > pb_local_etime
Definition: ERF_TurbPertStruct.H:562
amrex::Vector< amrex::BoxArray > pb_ba
Definition: ERF_TurbPertStruct.H:530
amrex::Vector< amrex::Real > pb_amp
Definition: ERF_TurbPertStruct.H:563
void netZeroBuoyantAdjust(const int &boxIdx, const int &lev)
Definition: ERF_TurbPertStruct.H:358
void calc_tpi_update(const int lev, const amrex::Real dt, amrex::MultiFab &mf_xvel, amrex::MultiFab &mf_yvel, amrex::MultiFab &mf_cons)
Definition: ERF_TurbPertStruct.H:185
void calc_tpi_amp(const int &boxIdx, const amrex::Real &interval)
Definition: ERF_TurbPertStruct.H:270
amrex::MultiFab pb_cell
Definition: ERF_TurbPertStruct.H:535
amrex::Vector< amrex::Real > pb_netZero
Definition: ERF_TurbPertStruct.H:564
int tpi_offset
Definition: ERF_TurbPertStruct.H:541
amrex::Real tpi_Ti
Definition: ERF_TurbPertStruct.H:548
amrex::Real tpi_pert_adjust
Definition: ERF_TurbPertStruct.H:558
amrex::Vector< amrex::Real > pb_interval
Definition: ERF_TurbPertStruct.H:561
void debug(amrex::Real)
Definition: ERF_TurbPertStruct.H:509
amrex::Real * get_pb_netZero()
Definition: ERF_TurbPertStruct.H:506
amrex::Vector< amrex::Real > pb_mag
Definition: ERF_TurbPertStruct.H:531
amrex::Real tpi_lref
Definition: ERF_TurbPertStruct.H:555
amrex::Vector< int > tpi_boxDim
Definition: ERF_TurbPertStruct.H:543
amrex::Real tpi_nonDim
Definition: ERF_TurbPertStruct.H:547
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:381
void init_tpi(const int lev, const amrex::IntVect &nx, const amrex::GpuArray< amrex::Real, 3 > dx, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const int ngrow_state, std::string pp_prefix)
Definition: ERF_TurbPertStruct.H:29
void netZeroBuoyantAdd(const int &boxIdx, const int &lev)
Definition: ERF_TurbPertStruct.H:323
int pt_type
Definition: ERF_TurbPertStruct.H:527
~TurbulentPerturbation()
Definition: ERF_TurbPertStruct.H:22
amrex::Real RandomReal(const amrex::Real min, const amrex::Real max)
Definition: ERF_TurbPertStruct.H:567
amrex::Real tpi_net_buoyant
Definition: ERF_TurbPertStruct.H:557
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)
Definition: ERF_TurbPertStruct.H:246
amrex::Real tpi_Hpb
Definition: ERF_TurbPertStruct.H:552
amrex::Vector< int > tpi_direction
Definition: ERF_TurbPertStruct.H:544