5 #include "AMReX_MultiFab.H"
6 #include "AMReX_BCRec.H"
16 void ChopGrids2D (amrex::BoxArray& ba,
const amrex::Box& domain,
int target_size);
26 cons_to_prim(
const amrex::MultiFab& cons_state, amrex::MultiFab& S_prim,
int ng);
29 make_qt(
const amrex::MultiFab& cons_state, amrex::MultiFab&
qt,
int n_qstate_into_total);
34 void make_J (
const amrex::Geometry& geom,
35 amrex::MultiFab& z_phys_nd,
36 amrex::MultiFab& detJ_cc);
39 amrex::MultiFab& z_phys_nd,
48 amrex::MultiFab& z_phys_nd,
49 amrex::MultiFab& z_phys_cc);
55 amrex::MultiFab& yvel_out,
56 amrex::MultiFab& zvel_out,
57 const amrex::MultiFab& cons_in,
58 const amrex::MultiFab& xmom_in,
59 const amrex::MultiFab& ymom_in,
60 const amrex::MultiFab& zmom_in,
61 const amrex::Box& domain,
62 const amrex::Vector<amrex::BCRec>& domain_bcs_type_h,
63 const amrex::MultiFab* c_vfrac =
nullptr
70 const amrex::IntVect & xvel_ngrow,
71 const amrex::MultiFab& yvel_in,
72 const amrex::IntVect & yvel_ngrow,
73 const amrex::MultiFab& zvel_in,
74 const amrex::IntVect & zvel_ngrow,
75 const amrex::MultiFab& cons_in,
76 amrex::MultiFab& xmom_out,
77 amrex::MultiFab& ymom_out,
78 amrex::MultiFab& zmom_out,
79 const amrex::Box& domain,
80 const amrex::Vector<amrex::BCRec>& domain_bcs_type_h,
81 const amrex::MultiFab* c_vfrac =
nullptr
89 amrex::MultiFab&
xmom, amrex::MultiFab&
ymom, amrex::MultiFab&
zmom,
90 const amrex::Box& domain,
const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
98 amrex::FArrayBox& bdy_data_ylo, amrex::FArrayBox& bdy_data_yhi,
99 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& area_vec,
100 const amrex::Geometry& geom,
109 amrex::FArrayBox& bdy_data_xlo,
110 amrex::FArrayBox& bdy_data_xhi,
111 amrex::FArrayBox& bdy_data_ylo,
112 amrex::FArrayBox& bdy_data_yhi,
113 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& area_vec,
114 const amrex::Geometry& geom,
115 const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
122 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& vels_vec,
123 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& area_vec,
124 const amrex::Geometry& geom);
131 const amrex::Box& domain,
137 const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0),
138 const bool get_int_ng=
false);
145 const amrex::Box& domain,
146 const int& set_width,
151 const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0));
163 const amrex::Geometry& geom,
164 amrex::Vector<amrex::MultiFab>& S_rhs,
165 amrex::Vector<amrex::MultiFab>& S_cur_data,
166 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xlo,
167 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xhi,
168 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_ylo,
169 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_yhi,
170 std::unique_ptr<ReadBndryPlanes>& m_r2d);
179 const int& set_width,
180 const amrex::Geometry& geom,
185 amrex::Vector<amrex::BCRec>& domain_bcs_type,
186 amrex::Vector<amrex::MultiFab>& S_rhs_f,
187 amrex::Vector<amrex::MultiFab>& S_data_f);
195 amrex::MultiFab* vel_t_avg,
196 amrex::MultiFab&
xvel,
197 amrex::MultiFab&
yvel,
198 amrex::MultiFab&
zvel);
228 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>&
dx,
229 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& ProbLo,
230 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& ProbHi,
232 const amrex::Box& bx_xlo,
233 const amrex::Box& bx_xhi,
234 const amrex::Box& bx_ylo,
235 const amrex::Box& bx_yhi,
236 const amrex::Array4<const amrex::Real>& arr_xlo,
237 const amrex::Array4<const amrex::Real>& arr_xhi,
238 const amrex::Array4<const amrex::Real>& arr_ylo,
239 const amrex::Array4<const amrex::Real>& arr_yhi,
240 const amrex::Array4<const amrex::Real>& data_arr,
241 const amrex::Array4<amrex::Real>& rhs_arr)
243 amrex::IntVect iv = bx_xlo.type();
246 amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
260 amrex::Real delta = arr_xlo(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
262 rhs_arr(i,j,k,n+icomp) += Temp;
264 bx_xhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
278 amrex::Real delta = arr_xhi(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
280 rhs_arr(i,j,k,n+icomp) += Temp;
283 amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
291 amrex::Real delta = arr_ylo(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
293 rhs_arr(i,j,k,n+icomp) += Temp;
295 bx_yhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
303 amrex::Real delta = arr_yhi(i ,j ,k,n) - data_arr(i ,j ,k ,n+icomp);
305 rhs_arr(i,j,k,n+icomp) += Temp;
316 const amrex::iMultiFab& imask,
317 const int nghost = 0)
319 for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
321 const amrex::Box& bx = mfi.growntilebox(nghost);
324 auto dstFab = dst.array(mfi);
325 const auto maskFab = imask.const_array(mfi);
326 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
328 dstFab(i,j,k) *= maskFab(i,j,k);
338 const amrex::iMultiFab& imask,
339 const int nghost = 0)
341 for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
343 const amrex::Box& bx = mfi.growntilebox(nghost);
346 auto dstFab = dst.array(mfi);
347 const auto maskFab = imask.const_array(mfi);
348 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
350 dstFab(i,j,k) *= (1-maskFab(i,j,k));
358 amrex::Vector<amrex::IntVect>& xfaces,
359 amrex::Vector<amrex::IntVect>& yfaces,
360 amrex::Vector<amrex::IntVect>& zfaces,
361 const amrex::Geometry& geomdata,
362 std::unique_ptr<amrex::MultiFab>& z_phys_cc);
367 AMREX_GPU_HOST_DEVICE
375 return p0 * ( (
z - z1) * (
z - z2) ) / ( (
z0 - z1) * (
z0 - z2) )
376 + p1 * ( (
z -
z0) * (
z - z2) ) / ( (z1 -
z0) * (z1 - z2) )
377 + p2 * ( (
z -
z0) * (
z - z1) ) / ( (z2 -
z0) * (z2 - z1) );
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
Real z0
Definition: ERF_InitCustomPertVels_ScalarAdvDiff.H:8
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
amrex::Real Real
Definition: ERF_ShocInterface.H:19
void MomentumToVelocity(amrex::MultiFab &xvel_out, amrex::MultiFab &yvel_out, amrex::MultiFab &zvel_out, const amrex::MultiFab &cons_in, const amrex::MultiFab &xmom_in, const amrex::MultiFab &ymom_in, const amrex::MultiFab &zmom_in, const amrex::Box &domain, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h, const amrex::MultiFab *c_vfrac=nullptr)
void compute_influx_outflux_bdy(amrex::FArrayBox &bdy_data_xlo, amrex::FArrayBox &bdy_data_xhi, amrex::FArrayBox &bdy_data_ylo, amrex::FArrayBox &bdy_data_yhi, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &area_vec, const amrex::Geometry &geom, amrex::Real &influx, amrex::Real &outflux, const int n)
void VelocityToMomentum(const amrex::MultiFab &xvel_in, const amrex::IntVect &xvel_ngrow, const amrex::MultiFab &yvel_in, const amrex::IntVect &yvel_ngrow, const amrex::MultiFab &zvel_in, const amrex::IntVect &zvel_ngrow, const amrex::MultiFab &cons_in, amrex::MultiFab &xmom_out, amrex::MultiFab &ymom_out, amrex::MultiFab &zmom_out, const amrex::Box &domain, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h, const amrex::MultiFab *c_vfrac=nullptr)
void realbdy_compute_interior_ghost_rhs(const amrex::Real &total_time, const amrex::Real &delta_t, const amrex::Real &start_bdy_time, const amrex::Real &final_bdy_time, const amrex::Real &bdy_time_interval, const amrex::Real &nudge_factor, int width, const amrex::Geometry &geom, amrex::Vector< amrex::MultiFab > &S_rhs, amrex::Vector< amrex::MultiFab > &S_cur_data, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_xlo, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_xhi, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_ylo, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_yhi, std::unique_ptr< ReadBndryPlanes > &m_r2d)
void cons_to_prim(const amrex::MultiFab &cons_state, amrex::MultiFab &S_prim, int ng)
void WeatherDataInterpolation(const amrex::Real time)
void fine_compute_interior_ghost_rhs(const amrex::Real &time, const amrex::Real &delta_t, const int &width, const int &set_width, const amrex::Geometry &geom, ERFFillPatcher *FPr_c, ERFFillPatcher *FPr_u, ERFFillPatcher *FPr_v, ERFFillPatcher *FPr_w, amrex::Vector< amrex::BCRec > &domain_bcs_type, amrex::Vector< amrex::MultiFab > &S_rhs_f, amrex::Vector< amrex::MultiFab > &S_data_f)
void enforceInOutSolvability_bdy(const amrex::MultiFab &rho0, amrex::FArrayBox &bdy_data_xlo, amrex::FArrayBox &bdy_data_xhi, amrex::FArrayBox &bdy_data_ylo, amrex::FArrayBox &bdy_data_yhi, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &area_vec, const amrex::Geometry &geom, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real quad_interp_1d(amrex::Real z, amrex::Real z0, amrex::Real p0, amrex::Real z1, amrex::Real p1, amrex::Real z2, amrex::Real p2)
Definition: ERF_Utils.H:370
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_compute_relaxation(const int &icomp, const int &num_var, const int &width, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dx, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &ProbLo, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &ProbHi, const amrex::Real &F1, const amrex::Box &bx_xlo, const amrex::Box &bx_xhi, const amrex::Box &bx_ylo, const amrex::Box &bx_yhi, const amrex::Array4< const amrex::Real > &arr_xlo, const amrex::Array4< const amrex::Real > &arr_xhi, const amrex::Array4< const amrex::Real > &arr_ylo, const amrex::Array4< const amrex::Real > &arr_yhi, const amrex::Array4< const amrex::Real > &data_arr, const amrex::Array4< amrex::Real > &rhs_arr)
Definition: ERF_Utils.H:225
void thinbody_wall_dist(std::unique_ptr< amrex::MultiFab > &wdist, amrex::Vector< amrex::IntVect > &xfaces, amrex::Vector< amrex::IntVect > &yfaces, amrex::Vector< amrex::IntVect > &zfaces, const amrex::Geometry &geomdata, std::unique_ptr< amrex::MultiFab > &z_phys_cc)
amrex::BoxArray ERFPostProcessBaseGrids(const amrex::Box &domain, bool decompose_in_z)
void make_zcc(const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc)
void realbdy_bc_bxs_xy(const amrex::Box &bx, const amrex::Box &domain, const int &set_width, amrex::Box &bx_xlo, amrex::Box &bx_xhi, amrex::Box &bx_ylo, amrex::Box &bx_yhi, const amrex::IntVect &ng_vect=amrex::IntVect(0, 0, 0))
void realbdy_interior_bxs_xy(const amrex::Box &bx, const amrex::Box &domain, const int &width, amrex::Box &bx_xlo, amrex::Box &bx_xhi, amrex::Box &bx_ylo, amrex::Box &bx_yhi, const amrex::IntVect &ng_vect=amrex::IntVect(0, 0, 0), const bool get_int_ng=false)
void Time_Avg_Vel_atCC(const amrex::Real &dt, amrex::Real &t_avg_cnt, amrex::MultiFab *vel_t_avg, amrex::MultiFab &xvel, amrex::MultiFab &yvel, amrex::MultiFab &zvel)
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask(amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
Definition: ERF_Utils.H:337
void make_J(const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &detJ_cc)
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyMask(amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
Definition: ERF_Utils.H:315
void ConvertForProjection(const amrex::MultiFab &den_div, const amrex::MultiFab &den_mlt, amrex::MultiFab &xmom, amrex::MultiFab &ymom, amrex::MultiFab &zmom, const amrex::Box &domain, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h)
void ChopGrids2D(amrex::BoxArray &ba, const amrex::Box &domain, int target_size)
void enforceInOutSolvability(int lev, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &vels_vec, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &area_vec, const amrex::Geometry &geom)
void make_areas(const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &ax, amrex::MultiFab &ay, amrex::MultiFab &az)
void make_qt(const amrex::MultiFab &cons_state, amrex::MultiFab &qt, int n_qstate_into_total)
Definition: ERF_FillPatcher.H:9
@ ymom
Definition: ERF_IndexDefines.H:194
@ zmom
Definition: ERF_IndexDefines.H:195
@ xmom
Definition: ERF_IndexDefines.H:193
@ qt
Definition: ERF_Kessler.H:28
@ ng
Definition: ERF_Morrison.H:48
@ xvel
Definition: ERF_IndexDefines.H:175
@ zvel
Definition: ERF_IndexDefines.H:177
@ yvel
Definition: ERF_IndexDefines.H:176
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40