1 #ifndef ERF_TURB_PERT_STRUCT_H_
2 #define ERF_TURB_PERT_STRUCT_H_
5 #include <AMReX_MultiFabUtil.H>
19 Source, Direct, CPM, None
30 if (pert_type == PerturbationType::Source) {
32 }
else if (pert_type == PerturbationType::Direct) {
34 }
else if (pert_type == PerturbationType::CPM) {
37 AMREX_ALWAYS_ASSERT(
pt_type >= 0);
46 const amrex::IntVect& valid_box_lo,
47 const amrex::IntVect& valid_box_hi,
48 const amrex::GpuArray<amrex::Real,3> dx,
49 const amrex::BoxArray& ba,
50 const amrex::DistributionMapping& dm,
51 const int ngrow_state,
52 std::string pp_prefix,
53 const amrex::Vector<amrex::IntVect> refRatio,
63 amrex::ParmParse
pp(pp_prefix);
73 pp.query(
"perturbation_T_intensity",
tpi_Ti);
77 if (
tpi_layers < 0) { amrex::Abort(
"Please provide a valid perturbation layer value (ie. 3-5)"); }
78 if (
tpi_offset < 0) { amrex::Abort(
"Please provide a valid inflow cell offset value for perturbation region (ie. 0-5)"); }
80 if (
tpi_boxDim[i] == 0) { amrex::Abort(
"Please provide valid dimensions for perturbation boxes."); }
82 if (
tpi_nonDim < 0.) { amrex::Abort(
"Please provide a valid nondimensional number (ie. Ri = 0.042)"); }
83 if (
input_Ug < 0.) { amrex::Abort(
"Please provide a valid geostrophic wind speed (ie. Ug = 10.0 m/s)"); }
84 if (
tpi_Tinf < 0.) { amrex::Abort(
"Please provide a valid ambient temperature value (ie. T_0 = T_infty)"); }
85 if (
tpi_Ti < 0.) { amrex::Abort(
"Please provide a valid temperature intensity value (ie. 0-1.0)"); }
88 amrex::Box lo_x_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
89 amrex::Box hi_x_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
90 amrex::Box lo_y_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
91 amrex::Box hi_y_bx(amrex::IntVect(0), amrex::IntVect(1), amrex::IntVect(0));
94 amrex::BoxList tmp_bl;
131 amrex::Box lo_x_lo_y_u = lo_x_bx & lo_y_bx;
132 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)));
136 amrex::Box hi_x_lo_y_u = hi_x_bx & lo_y_bx;
137 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)));
141 amrex::Box lo_x_hi_y_u = lo_x_bx & hi_y_bx;
142 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)));
146 amrex::Box hi_x_hi_y_u = hi_x_bx & hi_y_bx;
147 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)));
156 amrex::BoxArray tmp_ba(tmp_bl);
157 tmp_ba.maxSize(boxSize);
158 pb_ba.push_back(tmp_ba);
161 pb_mag.resize(max_level+1);
162 pb_dir.resize(max_level+1);
166 pb_amp.resize(max_level+1);
185 pb_cell[lev].define(ba, dm, 1, ngrow_state);
224 const amrex::Real dt,
225 amrex::MultiFab& mf_xvel,
226 amrex::MultiFab& mf_yvel,
227 amrex::MultiFab& mf_cons)
233 srand( (
unsigned) time(NULL) );
235 auto m_ixtype = mf_cons.boxArray().ixType();
238 int fix_random_seed = 0;
239 amrex::ParmParse
pp(
"erf");
240 pp.query(
"fix_random_seed", fix_random_seed);
241 if (fix_random_seed) {
243 amrex::InitRandom(1024UL);
249 for (
int boxIdx = 0; boxIdx <
pb_ba[lev].size(); boxIdx++) {
251 bool update_box =
true;
266 if (
pb_mag[lev][boxIdx] !=0.) {
267 amrex::Real interval = 0.0;
270 amrex::Real wind_direction =
pb_dir[lev][boxIdx];
271 if (wind_direction >
PI / 4) { wind_direction =
PI / 2 - wind_direction; }
273 interval = 1.0 / std::cos(wind_direction) * std::min(
tpi_Lpb[lev],
tpi_Wpb[lev]) /
pb_mag[lev][boxIdx];
315 for (
int boxIdx = 0; boxIdx <
pb_ba[lev].size(); boxIdx++) {
325 const amrex::Box& vbx,
327 const amrex::IndexType& m_ixtype,
328 const amrex::Array4<amrex::Real>& src_arr,
329 const amrex::Array4<amrex::Real const>& pert_cell)
331 for (
int boxIdx = 0; boxIdx <
pb_ba[lev].size(); boxIdx++) {
332 amrex::Box pbx = amrex::convert(
pb_ba[lev][boxIdx], m_ixtype);
333 amrex::Box ubx = pbx & vbx;
335 ParallelFor(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
336 src_arr(i,j,k,comp) += pert_cell(i,j,k);
340 src_arr(i,j,k,comp) = (amrex::Real) (boxIdx + 5.);
350 const amrex::Real& interval)
354 amrex::Real
cp = 1004;
355 amrex::Real Ec = 0.2;
358 amrex::Real Um =
pb_mag[lev][boxIdx];
364 int total_ref_ratio = 1;
365 for (
int level = lev; level >= 1; level--) {
366 total_ref_ratio *=
ref_ratio[level-1][2];
378 pb_amp[lev][boxIdx] /= interval;
388 const amrex::IndexType& m_ixtype)
390 for (amrex::MFIter mfi(
pb_cell[lev],
TileNoZ()); mfi.isValid(); ++mfi) {
391 amrex::Box vbx = mfi.validbox();
392 amrex::Box pbx = amrex::convert(
pb_ba[lev][boxIdx], m_ixtype);
393 amrex::Box ubx = pbx & vbx;
395 amrex::Real amp_copy =
pb_amp[lev][boxIdx];
396 const amrex::Array4<amrex::Real>& pert_cell =
pb_cell[lev].array(mfi);
399 amrex::Real rand_number_const =
RandomReal(-1.0, 1.0);
400 ParallelFor(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
401 pert_cell(i,j,k) = rand_number_const * amp_copy;
404 ParallelForRNG(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k,
const amrex::RandomEngine& engine) noexcept {
405 amrex::Real rand_double = amrex::Random(engine);
406 pert_cell(i,j,k) = (rand_double*2.0 - 1.0) * amp_copy;
417 const amrex::IndexType& m_ixtype)
420 for (amrex::MFIter mfi(
pb_cell[lev],
TileNoZ()); mfi.isValid(); ++mfi) {
421 amrex::Box vbx = mfi.validbox();
422 amrex::Box pbx = amrex::convert(
pb_ba[lev][boxIdx], m_ixtype);
423 amrex::Box ubx = pbx & vbx;
425 const amrex::Array4<amrex::Real>& pert_cell =
pb_cell[lev].array(mfi);
426 ParallelFor(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
427 pert_cell(i,j,k) = 0.0;
438 const amrex::BoxArray m_pb_ba =
pb_ba[lev];
439 amrex::Real* m_pb_netZero =
pb_netZero[lev].data();
442 amrex::Vector<amrex::Real>
avg_h(1,0.);
443 amrex::Gpu::DeviceVector<amrex::Real> avg_d(1,0.);
444 amrex::Real* avg = avg_d.data();
447 for (amrex::MFIter mfi(
pb_cell[lev],
TileNoZ()) ; mfi.isValid(); ++mfi) {
448 const amrex::Box& vbx = mfi.validbox();
449 amrex::Box pbx = amrex::convert(m_pb_ba[boxIdx], vbx.ixType());
450 amrex::Box ubx = pbx & vbx;
452 const amrex::Array4<const amrex::Real>& pert_cell =
pb_cell[lev].const_array(mfi);
453 int npts = ubx.numPts();
454 amrex::Real norm = 1.0 / (amrex::Real) npts;
455 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubx, [=]
456 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
457 amrex::Gpu::deviceReduceSum(&avg[0], pert_cell(i,j,k)*norm, handler);
459 amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(),
avg_h.begin());
462 m_pb_netZero[boxIdx] =
avg_h[0];
473 const amrex::BoxArray m_pb_ba =
pb_ba[lev];
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;
480 const amrex::Array4<amrex::Real>& pert_cell =
pb_cell[lev].array(mfi);
481 ParallelFor(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
482 pert_cell(i,j,k) -= adjust;
489 #define USE_VOLUME_AVERAGE
494 amrex::MultiFab& mf_cons,
495 amrex::MultiFab& mf_xvel,
496 amrex::MultiFab& mf_yvel)
500 const amrex::BoxArray m_pb_ba =
pb_ba[lev];
501 amrex::Real* m_pb_mag =
pb_mag[lev].data();
502 amrex::Real* m_pb_dir =
pb_dir[lev].data();
503 m_pb_mag[boxIdx] = 0.;
504 m_pb_dir[boxIdx] = 0.;
510 amrex::Vector<amrex::Real>
avg_h(n_avg,0.);
511 amrex::Gpu::DeviceVector<amrex::Real> avg_d(n_avg,0.);
512 amrex::Real* avg = avg_d.data();
515 for (amrex::MFIter mfi(mf_cons,
TileNoZ()); mfi.isValid(); ++mfi) {
518 const amrex::Box& vbx = mfi.validbox();
521 auto ixtype_u = mf_xvel.boxArray().ixType();
522 amrex::Box vbx_u = amrex::convert(vbx,ixtype_u);
523 amrex::Box pbx_u = amrex::convert(m_pb_ba[boxIdx], ixtype_u);
524 amrex::Box ubx_u = pbx_u & vbx_u;
527 auto ixtype_v = mf_yvel.boxArray().ixType();
528 amrex::Box vbx_v = amrex::convert(vbx,ixtype_v);
529 amrex::Box pbx_v = amrex::convert(m_pb_ba[boxIdx], ixtype_v);
530 amrex::Box ubx_v = pbx_v & vbx_v;
534 const amrex::Array4<const amrex::Real>& xvel_arry = mf_xvel.const_array(mfi);
536 #ifdef USE_VOLUME_AVERAGE
537 int npts = ubx_u.numPts();
538 amrex::Real norm = 1.0 / (amrex::Real) npts;
539 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubx_u, [=]
540 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
541 amrex::Gpu::deviceReduceSum(&avg[0], xvel_arry(i,j,k)*norm, handler);
545 #ifdef USE_SLAB_AVERAGE
546 amrex::Box ubxSlab_lo = makeSlab(ubx_u,2,ubx_u.smallEnd(2));
547 amrex::Box ubxSlab_hi = makeSlab(ubx_u,2,ubx_u.bigEnd(2));
548 int npts_lo = ubxSlab_lo.numPts();
549 int npts_hi = ubxSlab_hi.numPts();
550 amrex::Real norm_lo = 1.0 / (amrex::Real) npts_lo;
551 amrex::Real norm_hi = 1.0 / (amrex::Real) npts_hi;
554 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubxSlab_lo, [=]
555 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
556 amrex::Gpu::deviceReduceSum(&avg[0], xvel_arry(i,j,k)*norm_lo, handler);
560 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubxSlab_hi, [=]
561 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
562 amrex::Gpu::deviceReduceSum(&avg[2], xvel_arry(i,j,k)*norm_hi, handler);
569 const amrex::Array4<const amrex::Real>& yvel_arry = mf_yvel.const_array(mfi);
571 #ifdef USE_VOLUME_AVERAGE
572 int npts = ubx_v.numPts();
573 amrex::Real norm = 1.0 / (amrex::Real) npts;
574 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubx_v, [=]
575 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
576 amrex::Gpu::deviceReduceSum(&avg[1], yvel_arry(i,j,k)*norm, handler);
580 #ifdef USE_SLAB_AVERAGE
581 amrex::Box ubxSlab_lo = makeSlab(ubx_v,2,ubx_v.smallEnd(2));
582 amrex::Box ubxSlab_hi = makeSlab(ubx_v,2,ubx_v.bigEnd(2));
583 int npts_lo = ubxSlab_lo.numPts();
584 int npts_hi = ubxSlab_hi.numPts();
585 amrex::Real norm_lo = 1.0 / (amrex::Real) npts_lo;
586 amrex::Real norm_hi = 1.0 / (amrex::Real) npts_hi;
589 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubxSlab_lo, [=]
590 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
591 amrex::Gpu::deviceReduceSum(&avg[1], yvel_arry(i,j,k)*norm_lo, handler);
595 ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), ubxSlab_hi, [=]
596 AMREX_GPU_DEVICE(
int i,
int j,
int k, amrex::Gpu::Handler
const& handler) noexcept {
597 amrex::Gpu::deviceReduceSum(&avg[3], yvel_arry(i,j,k)*norm_hi, handler);
604 amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(),
avg_h.begin());
607 #ifdef USE_VOLUME_AVERAGE
612 #ifdef USE_SLAB_AVERAGE
641 amrex::Vector<amrex::BoxArray>
pb_ba;
642 amrex::Vector<amrex::Vector<amrex::Real>>
pb_mag;
643 amrex::Vector<amrex::Vector<amrex::Real>>
pb_dir;
678 amrex::Vector<amrex::Vector<amrex::Real>>
pb_amp;
682 amrex::Real
RandomReal (
const amrex::Real min,
const amrex::Real max)
684 amrex::Real r = (amrex::Real) rand() / (amrex::Real) RAND_MAX;
685 return min + r * (max - min);
constexpr amrex::Real PI
Definition: ERF_Constants.H:6
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_ENUM(PerturbationType, Source, Direct, CPM, None)
real(c_double), parameter g
Definition: ERF_module_model_constants.F90:19
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
real(c_double), parameter cp
Definition: ERF_module_model_constants.F90:22
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:22
amrex::Vector< amrex::MultiFab > pb_cell
Definition: ERF_TurbPertStruct.H:647
int tpi_layers
Definition: ERF_TurbPertStruct.H:652
void pseudoRandomPert(const int &boxIdx, const int &lev, const amrex::IndexType &m_ixtype)
Definition: ERF_TurbPertStruct.H:386
amrex::Real tpi_Tinf
Definition: ERF_TurbPertStruct.H:661
amrex::Vector< amrex::Vector< amrex::Real > > pb_local_etime
Definition: ERF_TurbPertStruct.H:677
amrex::Vector< amrex::IntVect > ref_ratio
Definition: ERF_TurbPertStruct.H:672
void calc_tpi_amp(const int &lev, const int &boxIdx, const amrex::Real &interval)
Definition: ERF_TurbPertStruct.H:348
amrex::Real input_Ug
Definition: ERF_TurbPertStruct.H:673
amrex::Vector< amrex::BoxArray > pb_ba
Definition: ERF_TurbPertStruct.H:641
void netZeroBuoyantAdjust(const int &boxIdx, const int &lev)
Definition: ERF_TurbPertStruct.H:469
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:223
amrex::Vector< amrex::Vector< amrex::Real > > pb_mag
Definition: ERF_TurbPertStruct.H:642
amrex::Vector< amrex::Real > tpi_lref
Definition: ERF_TurbPertStruct.H:667
int tpi_offset
Definition: ERF_TurbPertStruct.H:653
amrex::Vector< amrex::Real > tpi_Hpb
Definition: ERF_TurbPertStruct.H:664
void init_tpi_type(const PerturbationType &pert_type)
Definition: ERF_TurbPertStruct.H:28
amrex::Vector< amrex::Vector< amrex::Real > > pb_netZero
Definition: ERF_TurbPertStruct.H:679
amrex::Real tpi_Ti
Definition: ERF_TurbPertStruct.H:660
amrex::Real tpi_pert_adjust
Definition: ERF_TurbPertStruct.H:670
void debug(amrex::Real)
Definition: ERF_TurbPertStruct.H:620
amrex::Vector< amrex::Vector< amrex::Real > > pb_interval
Definition: ERF_TurbPertStruct.H:676
amrex::Vector< amrex::Vector< amrex::Real > > pb_dir
Definition: ERF_TurbPertStruct.H:643
void zero_amp(const int &boxIdx, const int &lev, const amrex::IndexType &m_ixtype)
Definition: ERF_TurbPertStruct.H:415
amrex::Vector< int > tpi_boxDim
Definition: ERF_TurbPertStruct.H:655
amrex::Vector< amrex::Real > tpi_Wpb
Definition: ERF_TurbPertStruct.H:666
amrex::Real tpi_nonDim
Definition: ERF_TurbPertStruct.H:659
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:492
void init_tpi(const int lev, const amrex::IntVect &valid_box_lo, const amrex::IntVect &valid_box_hi, 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)
Definition: ERF_TurbPertStruct.H:45
void netZeroBuoyantAdd(const int &boxIdx, const int &lev)
Definition: ERF_TurbPertStruct.H:434
int pt_type
Definition: ERF_TurbPertStruct.H:638
~TurbulentPerturbation()
Definition: ERF_TurbPertStruct.H:26
amrex::Real RandomReal(const amrex::Real min, const amrex::Real max)
Definition: ERF_TurbPertStruct.H:682
amrex::Vector< amrex::Real > tpi_Lpb
Definition: ERF_TurbPertStruct.H:665
amrex::Real tpi_net_buoyant
Definition: ERF_TurbPertStruct.H:669
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:324
amrex::Vector< int > tpi_direction
Definition: ERF_TurbPertStruct.H:656
amrex::Vector< amrex::Vector< amrex::Real > > pb_amp
Definition: ERF_TurbPertStruct.H:678