5 #include "AMReX_MultiFab.H"
6 #include "AMReX_BCRec.H"
15 void ChopGrids2D (amrex::BoxArray& ba,
const amrex::Box& domain,
int target_size);
27 void make_J (
const amrex::Geometry& geom,
28 amrex::MultiFab& z_phys_nd,
29 amrex::MultiFab& detJ_cc);
32 amrex::MultiFab& z_phys_nd,
41 amrex::MultiFab& z_phys_nd,
42 amrex::MultiFab& z_phys_cc);
48 amrex::MultiFab& yvel_out,
49 amrex::MultiFab& zvel_out,
50 const amrex::MultiFab& cons_in,
51 const amrex::MultiFab& xmom_in,
52 const amrex::MultiFab& ymom_in,
53 const amrex::MultiFab& zmom_in,
54 const amrex::Box& domain,
55 const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
61 const amrex::IntVect & xvel_ngrow,
62 const amrex::MultiFab& yvel_in,
63 const amrex::IntVect & yvel_ngrow,
64 const amrex::MultiFab& zvel_in,
65 const amrex::IntVect & zvel_ngrow,
66 const amrex::MultiFab& cons_in,
67 amrex::MultiFab& xmom_out,
68 amrex::MultiFab& ymom_out,
69 amrex::MultiFab& zmom_out,
70 const amrex::Box& domain,
71 const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
78 amrex::MultiFab&
xmom, amrex::MultiFab&
ymom, amrex::MultiFab&
zmom,
79 const amrex::Box& domain,
const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
86 const amrex::Box& domain,
92 const int& set_width=0,
93 const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0),
94 const bool get_int_ng=
false);
101 const amrex::Box& domain,
102 const int& set_width,
107 const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0));
118 const amrex::Geometry& geom,
119 amrex::Vector<amrex::MultiFab>& S_rhs,
120 amrex::Vector<amrex::MultiFab>& S_old_data,
121 amrex::Vector<amrex::MultiFab>& S_cur_data,
122 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xlo,
123 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xhi,
124 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_ylo,
125 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_yhi);
134 const int& set_width,
135 const amrex::Geometry& geom,
140 amrex::Vector<amrex::BCRec>& domain_bcs_type,
141 amrex::Vector<amrex::MultiFab>& S_rhs_f,
142 amrex::Vector<amrex::MultiFab>& S_data_f);
150 amrex::MultiFab* vel_t_avg,
151 amrex::MultiFab&
xvel,
152 amrex::MultiFab&
yvel,
153 amrex::MultiFab&
zvel);
173 const int& set_width_x,
174 const int& set_width_y,
175 const amrex::Dim3&
dom_lo,
176 const amrex::Dim3&
dom_hi,
177 const amrex::Box& bx_xlo,
178 const amrex::Box& bx_xhi,
179 const amrex::Box& bx_ylo,
180 const amrex::Box& bx_yhi,
181 const amrex::Array4<const amrex::Real>& arr_xlo,
182 const amrex::Array4<const amrex::Real>& arr_xhi,
183 const amrex::Array4<const amrex::Real>& arr_ylo,
184 const amrex::Array4<const amrex::Real>& arr_yhi,
185 const amrex::Array4<const amrex::Real>& data_arr,
186 const amrex::Array4<amrex::Real>& rhs_arr)
188 int Spec_z_x = set_width_x;
189 int Spec_z_y = set_width_y;
190 amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
193 int j_lo = std::min(j-
dom_lo.y,width-1);
194 int j_hi = std::min(
dom_hi.y-j,width-1);
195 int jj = std::min(j_lo,j_hi);
196 int n_ind_y = jj + 1;
197 int n_ind_x = i -
dom_lo.x + 1;
198 if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
199 amrex::Real tend = ( arr_xlo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
200 rhs_arr(i,j,k,n+icomp) = tend;
203 bx_xhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
206 int j_lo = std::min(j-
dom_lo.y,width-1);
207 int j_hi = std::min(
dom_hi.y-j,width-1);
208 int jj = std::min(j_lo,j_hi);
209 int n_ind_y = jj + 1;
210 int n_ind_x =
dom_hi.x - i + 1;
211 if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
212 amrex::Real tend = ( arr_xhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
213 rhs_arr(i,j,k,n+icomp) = tend;
217 amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
220 int n_ind = j -
dom_lo.y + 1;
221 if (n_ind <= Spec_z_y) {
222 amrex::Real tend = ( arr_ylo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
223 rhs_arr(i,j,k,n+icomp) = tend;
226 bx_yhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
229 int n_ind =
dom_hi.y - j + 1;
230 if (n_ind <= Spec_z_y) {
231 amrex::Real tend = ( arr_yhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
232 rhs_arr(i,j,k,n+icomp) = tend;
266 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& dx,
267 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& ProbLo,
268 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& ProbHi,
271 const amrex::Box& bx_xlo,
272 const amrex::Box& bx_xhi,
273 const amrex::Box& bx_ylo,
274 const amrex::Box& bx_yhi,
275 const amrex::Array4<const amrex::Real>& arr_xlo,
276 const amrex::Array4<const amrex::Real>& arr_xhi,
277 const amrex::Array4<const amrex::Real>& arr_ylo,
278 const amrex::Array4<const amrex::Real>& arr_yhi,
279 const amrex::Array4<const amrex::Real>& data_arr,
280 const amrex::Array4<amrex::Real>& rhs_arr)
283 amrex::IntVect iv = bx_xlo.type();
286 amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
291 amrex::Real x_end = ProbLo[0] + (width + ioff) * dx[0];
292 amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
293 amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
295 amrex::Real eta_lo = (
y < y_end ) ? (y_end -
y) / (y_end - ProbLo[1]) : 0.0;
296 amrex::Real eta_hi = (
y > y_strt) ? (
y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
306 amrex::Real delta_xp = arr_xlo(i+1,j ,k,n) - d_ip1;
307 amrex::Real delta_xm = arr_xlo(i-1,j ,k,n) - d_im1;
308 amrex::Real delta_yp = arr_xlo(i ,j+1,k,n) - d_jp1;
309 amrex::Real delta_ym = arr_xlo(i ,j-1,k,n) - d_jm1;
310 amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
311 amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
312 rhs_arr(i,j,k,n+icomp) += Temp;
314 bx_xhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
319 amrex::Real x_strt = ProbHi[0] - (width + ioff) * dx[0];
320 amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
321 amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
323 amrex::Real eta_lo = (
y < y_end ) ? (y_end -
y) / (y_end - ProbLo[1]) : 0.0;
324 amrex::Real eta_hi = (
y > y_strt) ? (
y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
334 amrex::Real delta_xp = arr_xhi(i+1,j ,k,n) - d_ip1;
335 amrex::Real delta_xm = arr_xhi(i-1,j ,k,n) - d_im1;
336 amrex::Real delta_yp = arr_xhi(i ,j+1,k,n) - d_jp1;
337 amrex::Real delta_ym = arr_xhi(i ,j-1,k,n) - d_jm1;
338 amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
339 amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
340 rhs_arr(i,j,k,n+icomp) += Temp;
343 amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
347 amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
357 amrex::Real delta_xp = arr_ylo(i+1,j ,k,n) - d_ip1;
358 amrex::Real delta_xm = arr_ylo(i-1,j ,k,n) - d_im1;
359 amrex::Real delta_yp = arr_ylo(i ,j+1,k,n) - d_jp1;
360 amrex::Real delta_ym = arr_ylo(i ,j-1,k,n) - d_jm1;
361 amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
362 amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
363 rhs_arr(i,j,k,n+icomp) += Temp;
365 bx_yhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
369 amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
379 amrex::Real delta_xp = arr_yhi(i+1,j ,k,n) - d_ip1;
380 amrex::Real delta_xm = arr_yhi(i-1,j ,k,n) - d_im1;
381 amrex::Real delta_yp = arr_yhi(i ,j+1,k,n) - d_jp1;
382 amrex::Real delta_ym = arr_yhi(i ,j-1,k,n) - d_jm1;
383 amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
384 amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
385 rhs_arr(i,j,k,n+icomp) += Temp;
396 const amrex::iMultiFab& imask,
397 const int nghost = 0)
399 for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
401 const amrex::Box& bx = mfi.growntilebox(nghost);
404 auto dstFab = dst.array(mfi);
405 const auto maskFab = imask.const_array(mfi);
406 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
408 dstFab(i,j,k) *= maskFab(i,j,k);
418 const amrex::iMultiFab& imask,
419 const int nghost = 0)
421 for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
423 const amrex::Box& bx = mfi.growntilebox(nghost);
426 auto dstFab = dst.array(mfi);
427 const auto maskFab = imask.const_array(mfi);
428 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
430 dstFab(i,j,k) *= (1-maskFab(i,j,k));
438 amrex::Vector<amrex::IntVect>& xfaces,
439 amrex::Vector<amrex::IntVect>& yfaces,
440 amrex::Vector<amrex::IntVect>& zfaces,
441 const amrex::Geometry& geomdata,
442 std::unique_ptr<amrex::MultiFab>& z_phys_cc);
const auto & dom_hi
Definition: ERF_DiffSetup.H:10
const auto & dom_lo
Definition: ERF_DiffSetup.H:9
amrex::Real Real
Definition: ERF_ShocInterface.H:16
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_set_rhs_in_spec_region(const amrex::Real &dt, const int &icomp, const int &num_var, const int &width, const int &set_width_x, const int &set_width_y, const amrex::Dim3 &dom_lo, const amrex::Dim3 &dom_hi, 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:169
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)
void realbdy_compute_interior_ghost_rhs(const amrex::Real &bdy_time_interval, const amrex::Real &time, const amrex::Real &delta_t, const amrex::Real &stop_time, int width, int set_width, const amrex::Geometry &geom, amrex::Vector< amrex::MultiFab > &S_rhs, amrex::Vector< amrex::MultiFab > &S_old_data, 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)
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)
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 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)
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 int &set_width=0, const amrex::IntVect &ng_vect=amrex::IntVect(0, 0, 0), const bool get_int_ng=false)
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))
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_compute_laplacian_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::Real &F2, 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:263
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:417
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:395
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 make_areas(const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &ax, amrex::MultiFab &ay, amrex::MultiFab &az)
Definition: ERF_FillPatcher.H:9
@ ymom
Definition: ERF_IndexDefines.H:160
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159
@ xvel
Definition: ERF_IndexDefines.H:141
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142