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,
56 const amrex::MultiFab* c_vfrac =
nullptr
63 const amrex::IntVect & xvel_ngrow,
64 const amrex::MultiFab& yvel_in,
65 const amrex::IntVect & yvel_ngrow,
66 const amrex::MultiFab& zvel_in,
67 const amrex::IntVect & zvel_ngrow,
68 const amrex::MultiFab& cons_in,
69 amrex::MultiFab& xmom_out,
70 amrex::MultiFab& ymom_out,
71 amrex::MultiFab& zmom_out,
72 const amrex::Box& domain,
73 const amrex::Vector<amrex::BCRec>& domain_bcs_type_h,
74 const amrex::MultiFab* c_vfrac =
nullptr
82 amrex::MultiFab&
xmom, amrex::MultiFab&
ymom, amrex::MultiFab&
zmom,
83 const amrex::Box& domain,
const amrex::Vector<amrex::BCRec>& domain_bcs_type_h);
90 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& vels_vec,
91 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& area_vec,
92 const amrex::Geometry& geom);
99 const amrex::Box& domain,
105 const int& set_width=0,
106 const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0),
107 const bool get_int_ng=
false);
114 const amrex::Box& domain,
115 const int& set_width,
120 const amrex::IntVect& ng_vect=amrex::IntVect(0,0,0));
131 const amrex::Geometry& geom,
132 amrex::Vector<amrex::MultiFab>& S_rhs,
133 amrex::Vector<amrex::MultiFab>& S_old_data,
134 amrex::Vector<amrex::MultiFab>& S_cur_data,
135 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xlo,
136 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xhi,
137 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_ylo,
138 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_yhi);
147 const int& set_width,
148 const amrex::Geometry& geom,
153 amrex::Vector<amrex::BCRec>& domain_bcs_type,
154 amrex::Vector<amrex::MultiFab>& S_rhs_f,
155 amrex::Vector<amrex::MultiFab>& S_data_f);
163 amrex::MultiFab* vel_t_avg,
164 amrex::MultiFab&
xvel,
165 amrex::MultiFab&
yvel,
166 amrex::MultiFab&
zvel);
186 const int& set_width_x,
187 const int& set_width_y,
188 const amrex::Dim3&
dom_lo,
189 const amrex::Dim3&
dom_hi,
190 const amrex::Box& bx_xlo,
191 const amrex::Box& bx_xhi,
192 const amrex::Box& bx_ylo,
193 const amrex::Box& bx_yhi,
194 const amrex::Array4<const amrex::Real>& arr_xlo,
195 const amrex::Array4<const amrex::Real>& arr_xhi,
196 const amrex::Array4<const amrex::Real>& arr_ylo,
197 const amrex::Array4<const amrex::Real>& arr_yhi,
198 const amrex::Array4<const amrex::Real>& data_arr,
199 const amrex::Array4<amrex::Real>& rhs_arr)
201 int Spec_z_x = set_width_x;
202 int Spec_z_y = set_width_y;
203 amrex::ParallelFor(bx_xlo, 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 = i -
dom_lo.x + 1;
211 if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
212 amrex::Real tend = ( arr_xlo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
213 rhs_arr(i,j,k,n+icomp) = tend;
216 bx_xhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
219 int j_lo = std::min(j-
dom_lo.y,width-1);
220 int j_hi = std::min(
dom_hi.y-j,width-1);
221 int jj = std::min(j_lo,j_hi);
222 int n_ind_y = jj + 1;
223 int n_ind_x =
dom_hi.x - i + 1;
224 if ((n_ind_x <= Spec_z_x) || (n_ind_y <= Spec_z_y)) {
225 amrex::Real tend = ( arr_xhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
226 rhs_arr(i,j,k,n+icomp) = tend;
230 amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
233 int n_ind = j -
dom_lo.y + 1;
234 if (n_ind <= Spec_z_y) {
235 amrex::Real tend = ( arr_ylo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
236 rhs_arr(i,j,k,n+icomp) = tend;
239 bx_yhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
242 int n_ind =
dom_hi.y - j + 1;
243 if (n_ind <= Spec_z_y) {
244 amrex::Real tend = ( arr_yhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
245 rhs_arr(i,j,k,n+icomp) = tend;
279 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& dx,
280 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& ProbLo,
281 const amrex::GpuArray<amrex::Real,AMREX_SPACEDIM>& ProbHi,
284 const amrex::Box& bx_xlo,
285 const amrex::Box& bx_xhi,
286 const amrex::Box& bx_ylo,
287 const amrex::Box& bx_yhi,
288 const amrex::Array4<const amrex::Real>& arr_xlo,
289 const amrex::Array4<const amrex::Real>& arr_xhi,
290 const amrex::Array4<const amrex::Real>& arr_ylo,
291 const amrex::Array4<const amrex::Real>& arr_yhi,
292 const amrex::Array4<const amrex::Real>& data_arr,
293 const amrex::Array4<amrex::Real>& rhs_arr)
296 amrex::IntVect iv = bx_xlo.type();
299 amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
304 amrex::Real x_end = ProbLo[0] + (width + ioff) * dx[0];
305 amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
306 amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
308 amrex::Real eta_lo = (
y < y_end ) ? (y_end -
y) / (y_end - ProbLo[1]) : 0.0;
309 amrex::Real eta_hi = (
y > y_strt) ? (
y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
319 amrex::Real delta_xp = arr_xlo(i+1,j ,k,n) - d_ip1;
320 amrex::Real delta_xm = arr_xlo(i-1,j ,k,n) - d_im1;
321 amrex::Real delta_yp = arr_xlo(i ,j+1,k,n) - d_jp1;
322 amrex::Real delta_ym = arr_xlo(i ,j-1,k,n) - d_jm1;
323 amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
324 amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
325 rhs_arr(i,j,k,n+icomp) += Temp;
327 bx_xhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
332 amrex::Real x_strt = ProbHi[0] - (width + ioff) * dx[0];
333 amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
334 amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
336 amrex::Real eta_lo = (
y < y_end ) ? (y_end -
y) / (y_end - ProbLo[1]) : 0.0;
337 amrex::Real eta_hi = (
y > y_strt) ? (
y - y_strt) / (ProbHi[1] - y_strt) : 0.0;
347 amrex::Real delta_xp = arr_xhi(i+1,j ,k,n) - d_ip1;
348 amrex::Real delta_xm = arr_xhi(i-1,j ,k,n) - d_im1;
349 amrex::Real delta_yp = arr_xhi(i ,j+1,k,n) - d_jp1;
350 amrex::Real delta_ym = arr_xhi(i ,j-1,k,n) - d_jm1;
351 amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
352 amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
353 rhs_arr(i,j,k,n+icomp) += Temp;
356 amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
360 amrex::Real y_end = ProbLo[1] + (width + joff) * dx[1];
370 amrex::Real delta_xp = arr_ylo(i+1,j ,k,n) - d_ip1;
371 amrex::Real delta_xm = arr_ylo(i-1,j ,k,n) - d_im1;
372 amrex::Real delta_yp = arr_ylo(i ,j+1,k,n) - d_jp1;
373 amrex::Real delta_ym = arr_ylo(i ,j-1,k,n) - d_jm1;
374 amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
375 amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
376 rhs_arr(i,j,k,n+icomp) += Temp;
378 bx_yhi, num_var, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
382 amrex::Real y_strt = ProbHi[1] - (width + joff) * dx[1];
392 amrex::Real delta_xp = arr_yhi(i+1,j ,k,n) - d_ip1;
393 amrex::Real delta_xm = arr_yhi(i-1,j ,k,n) - d_im1;
394 amrex::Real delta_yp = arr_yhi(i ,j+1,k,n) - d_jp1;
395 amrex::Real delta_ym = arr_yhi(i ,j-1,k,n) - d_jm1;
396 amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
397 amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
398 rhs_arr(i,j,k,n+icomp) += Temp;
409 const amrex::iMultiFab& imask,
410 const int nghost = 0)
412 for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
414 const amrex::Box& bx = mfi.growntilebox(nghost);
417 auto dstFab = dst.array(mfi);
418 const auto maskFab = imask.const_array(mfi);
419 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
421 dstFab(i,j,k) *= maskFab(i,j,k);
431 const amrex::iMultiFab& imask,
432 const int nghost = 0)
434 for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
436 const amrex::Box& bx = mfi.growntilebox(nghost);
439 auto dstFab = dst.array(mfi);
440 const auto maskFab = imask.const_array(mfi);
441 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
443 dstFab(i,j,k) *= (1-maskFab(i,j,k));
451 amrex::Vector<amrex::IntVect>& xfaces,
452 amrex::Vector<amrex::IntVect>& yfaces,
453 amrex::Vector<amrex::IntVect>& zfaces,
454 const amrex::Geometry& geomdata,
455 std::unique_ptr<amrex::MultiFab>& z_phys_cc);
460 AMREX_GPU_HOST_DEVICE
468 return p0 * ( (
z - z1) * (
z - z2) ) / ( (z0 - z1) * (z0 - z2) )
469 + p1 * ( (
z - z0) * (
z - z2) ) / ( (z1 - z0) * (z1 - z2) )
470 + p2 * ( (
z - z0) * (
z - z1) ) / ( (z2 - z0) * (z2 - z1) );
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
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)
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:182
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 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)
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:463
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 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_elapsed, 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 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:276
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:430
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:408
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)
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
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40