1 #ifndef ERF_TURB_PERT_STRUCT_H_
2 #define ERF_TURB_PERT_STRUCT_H_
5 #include <AMReX_MultiFabUtil.H>
19 Source, Direct, CPM, None
29 const PerturbationType& pert_type,
32 if (
pt_type.size() < max_level + 1) {
33 pt_type.resize(max_level + 1, -1);
36 if (pert_type == PerturbationType::Source) {
38 }
else if (pert_type == PerturbationType::Direct) {
40 }
else if (pert_type == PerturbationType::CPM) {
53 const amrex::Vector<amrex::BoxArray>& subdomains_lev,
54 const amrex::GpuArray<amrex::Real,3>
dx,
55 const amrex::BoxArray& ba,
56 const amrex::DistributionMapping& dm,
57 const int ngrow_state,
58 std::string pp_prefix,
59 const amrex::Vector<amrex::IntVect> refRatio,
65 amrex::ParmParse
pp(pp_prefix);
79 pp.query(
"perturbation_T_intensity",
tpi_Ti);
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))"); }
89 if (
tpi_boxDim[i] < 3) { amrex::Abort(
"Please provide valid dimensions for perturbation boxes."); }
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)"); }
96 amrex::BoxList tmp_bl;
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());
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.");
114 const amrex::IntVect& valid_box_lo = subdomain_box.smallEnd();
115 const amrex::IntVect& valid_box_hi = subdomain_box.bigEnd();
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));
129 amrex::PrintToFile(
"BoxPerturbationOutput") <<
" West face";
135 amrex::PrintToFile(
"BoxPerturbationOutput") <<
" East face";
142 amrex::PrintToFile(
"BoxPerturbationOutput") <<
" North face";
148 amrex::PrintToFile(
"BoxPerturbationOutput") <<
" South face";
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)));
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)));
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)));
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)));
180 amrex::BoxArray tmp_ba(tmp_bl);
181 tmp_ba.maxSize(boxSize);
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);
191 pb_amp.resize(num_levels);
212 pb_cell[lev].define(ba, dm, 1, ngrow_state);
252 amrex::MultiFab& mf_xvel,
253 amrex::MultiFab& mf_yvel,
254 amrex::MultiFab& mf_cons)
260 srand( (
unsigned) time(NULL) );
262 auto m_ixtype = mf_cons.boxArray().ixType();
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) {
270 amrex::InitRandom(1024UL, amrex::ParallelDescriptor::NProcs(), 1024UL);
276 for (
int boxIdx = 0; boxIdx <
pb_ba[lev].size(); boxIdx++) {
278 bool update_box =
true;
298 if (wind_direction >
PI / 4) { wind_direction =
PI / 2 - wind_direction; }
342 for (
int boxIdx = 0; boxIdx <
pb_ba[lev].size(); boxIdx++) {
352 const amrex::Box& vbx,
354 const amrex::IndexType& m_ixtype,
355 const amrex::Array4<amrex::Real>& src_arr,
356 const amrex::Array4<amrex::Real const>& pert_cell)
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;
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);
391 int total_ref_ratio = 1;
392 for (
int level = lev; level >= 1; level--) {
393 total_ref_ratio *=
ref_ratio[level-1][2];
405 pb_amp[lev][boxIdx] /= interval;
415 const amrex::IndexType& m_ixtype)
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;
423 const amrex::Array4<amrex::Real>& pert_cell =
pb_cell[lev].array(mfi);
427 ParallelFor(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
428 pert_cell(i,j,k) = rand_number_const * amp_copy;
431 ParallelForRNG(ubx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k,
const amrex::RandomEngine& engine) noexcept {
433 pert_cell(i,j,k) = (rand_double*
two -
one) * amp_copy;
444 const amrex::IndexType& m_ixtype)
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;
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;
465 const amrex::BoxArray m_pb_ba =
pb_ba[lev];
469 amrex::Vector<amrex::Real> avg_h(1,
zero);
470 amrex::Gpu::DeviceVector<amrex::Real> avg_d(1,
zero);
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;
479 const amrex::Array4<const amrex::Real>& pert_cell =
pb_cell[lev].const_array(mfi);
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);
485 amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(), avg_h.begin());
488 m_pb_netZero[boxIdx] = avg_h[0];
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;
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;
515 #define USE_VOLUME_AVERAGE
520 amrex::MultiFab& mf_cons,
521 amrex::MultiFab& mf_xvel,
522 amrex::MultiFab& mf_yvel)
526 const amrex::BoxArray m_pb_ba =
pb_ba[lev];
529 m_pb_mag[boxIdx] =
zero;
530 m_pb_dir[boxIdx] =
zero;
536 amrex::Vector<amrex::Real> avg_h(n_avg,
zero);
537 amrex::Gpu::DeviceVector<amrex::Real> avg_d(n_avg,
zero);
541 for (amrex::MFIter mfi(mf_cons,
TileNoZ()); mfi.isValid(); ++mfi) {
544 const amrex::Box& vbx = mfi.validbox();
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;
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;
560 const amrex::Array4<const amrex::Real>& xvel_arry = mf_xvel.const_array(mfi);
562 #ifdef USE_VOLUME_AVERAGE
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);
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));
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);
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);
592 const amrex::Array4<const amrex::Real>& yvel_arry = mf_yvel.const_array(mfi);
594 #ifdef USE_VOLUME_AVERAGE
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);
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));
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);
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);
624 amrex::Gpu::copy(amrex::Gpu::deviceToHost, avg_d.begin(), avg_d.end(), avg_h.begin());
627 #ifdef USE_VOLUME_AVERAGE
628 m_pb_mag[boxIdx] = std::sqrt(avg_h[0]*avg_h[0] + avg_h[1]*avg_h[1]);
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]));
661 amrex::Vector<amrex::BoxArray>
pb_ba;
662 amrex::Vector<amrex::Vector<amrex::Real>>
pb_mag;
663 amrex::Vector<amrex::Vector<amrex::Real>>
pb_dir;
698 amrex::Vector<amrex::Vector<amrex::Real>>
pb_amp;
705 return min + r * (max - min);
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real PI
Definition: ERF_Constants.H:16
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:32
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));} } })
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
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 beta
Definition: ERF_InitCustomPert_DataAssimilation_ISV.H:10
amrex::Real Real
Definition: ERF_ShocInterface.H:19
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
integer, private isub
Definition: ERF_module_mp_morr_two_moment.F90:164
Definition: ERF_TurbPertStruct.H:22
amrex::Vector< amrex::MultiFab > pb_cell
Definition: ERF_TurbPertStruct.H:667
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::Real tpi_Tinf
Definition: ERF_TurbPertStruct.H:681
amrex::Vector< amrex::Vector< amrex::Real > > pb_local_etime
Definition: ERF_TurbPertStruct.H:697
amrex::Vector< amrex::IntVect > ref_ratio
Definition: ERF_TurbPertStruct.H:692
void calc_tpi_amp(const int &lev, const int &boxIdx, const amrex::Real &interval)
Definition: ERF_TurbPertStruct.H:375
amrex::Real input_Ug
Definition: ERF_TurbPertStruct.H:693
amrex::Vector< amrex::BoxArray > pb_ba
Definition: ERF_TurbPertStruct.H:661
void netZeroBuoyantAdjust(const int &boxIdx, const int &lev)
Definition: ERF_TurbPertStruct.H:495
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:250
amrex::Vector< amrex::Vector< amrex::Real > > pb_mag
Definition: ERF_TurbPertStruct.H:662
amrex::Vector< amrex::Real > tpi_lref
Definition: ERF_TurbPertStruct.H:687
amrex::Vector< int > pt_type
Definition: ERF_TurbPertStruct.H:658
int tpi_offset
Definition: ERF_TurbPertStruct.H:673
amrex::Vector< amrex::Real > tpi_Hpb
Definition: ERF_TurbPertStruct.H:684
void init_tpi_type(const int lev, const PerturbationType &pert_type, const int max_level)
Definition: ERF_TurbPertStruct.H:28
amrex::Vector< amrex::Vector< amrex::Real > > pb_netZero
Definition: ERF_TurbPertStruct.H:699
amrex::Real tpi_Ti
Definition: ERF_TurbPertStruct.H:680
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)
Definition: ERF_TurbPertStruct.H:52
amrex::Real tpi_pert_adjust
Definition: ERF_TurbPertStruct.H:690
void debug(amrex::Real)
Definition: ERF_TurbPertStruct.H:640
amrex::Vector< amrex::Vector< amrex::Real > > pb_interval
Definition: ERF_TurbPertStruct.H:696
amrex::Vector< amrex::Vector< amrex::Real > > pb_dir
Definition: ERF_TurbPertStruct.H:663
void zero_amp(const int &boxIdx, const int &lev, const amrex::IndexType &m_ixtype)
Definition: ERF_TurbPertStruct.H:442
amrex::Vector< int > tpi_boxDim
Definition: ERF_TurbPertStruct.H:675
amrex::Vector< amrex::Real > tpi_Wpb
Definition: ERF_TurbPertStruct.H:686
amrex::Real tpi_nonDim
Definition: ERF_TurbPertStruct.H:679
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
~TurbulentPerturbation()
Definition: ERF_TurbPertStruct.H:26
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
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:351
amrex::Vector< int > tpi_direction
Definition: ERF_TurbPertStruct.H:676
amrex::Vector< amrex::Vector< amrex::Real > > pb_amp
Definition: ERF_TurbPertStruct.H:698