Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

EB mask operations #21

Merged
merged 6 commits into from
Oct 18, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
93 changes: 93 additions & 0 deletions Exec/RegTests/EB_EnclosedVortex/input.3d-regt
Original file line number Diff line number Diff line change
@@ -0,0 +1,93 @@
#----------------------DOMAIN DEFINITION------------------------
geometry.is_periodic = 0 0 0 # For each dir, 0: non-perio, 1: periodic
geometry.coord_sys = 0 # 0 => cart, 1 => RZ
geometry.prob_lo = -0.02 -0.02 -0.01 # x_lo y_lo (z_lo)
geometry.prob_hi = 0.02 0.02 0.01 # x_hi y_hi (z_hi)

# >>>>>>>>>>>>> BC FLAGS <<<<<<<<<<<<<<<<
# Interior, Inflow, Outflow, Symmetry,
# SlipWallAdiab, NoSlipWallAdiab, SlipWallIsotherm, NoSlipWallIsotherm
peleLM.lo_bc = NoSlipWallAdiab NoSlipWallAdiab Symmetry
peleLM.hi_bc = NoSlipWallAdiab NoSlipWallAdiab Symmetry


#-------------------------AMR CONTROL----------------------------
amr.n_cell = 64 64 32 # Level 0 number of cells in each direction
amr.v = 1 # AMR verbose
amr.max_level = 2 # maximum level number allowed
amr.ref_ratio = 2 2 2 2 2 # refinement ratio
amr.regrid_int = 5 # how often to regrid
amr.n_error_buf = 2 2 2 2 # number of buffer cells in error est
amr.grid_eff = 0.7 # what constitutes an efficient grid
amr.blocking_factor = 16 # block factor in grid generation (min box size)
amr.max_grid_size = 128 # max box size


#--------------------------- Problem -------------------------------
prob.P_mean = 101325.0
prob.rvort = 0.003
prob.forcevort = 0.2

#-------------------------PeleLM CONTROL----------------------------
peleLM.v = 3
peleLM.incompressible = 0
peleLM.rho = 1.17
peleLM.mu = 0.0
peleLM.use_wbar = 1
peleLM.sdc_iterMax = 1
peleLM.floor_species = 0
peleLM.num_divu_iter = 0
peleLM.num_init_iter = 1
peleLM.do_react = 0

peleLM.do_temporals = 1
peleLM.temporal_int = 2
peleLM.mass_balance = 1

#amr.restart = chk00005
amr.check_int = 5
amr.plot_int = 100
amr.max_step = 10
amr.dt_shrink = 0.1
amr.stop_time = 1.0
#amr.stop_time = 1.00
amr.cfl = 0.15
amr.derive_plot_vars = avg_pressure mag_vort mass_fractions

# --------------- INPUTS TO CHEMISTRY REACTOR ---------------
peleLM.chem_integrator = "ReactorNull"
#peleLM.use_typ_vals_chem = 1 # Use species/temp typical values in CVODE
#ode.rtol = 1.0e-6 # Relative tolerance of the chemical solve
#ode.atol = 1.0e-5 # Absolute tolerance factor applied on typical values
cvode.solve_type = denseAJ_direct # CVODE Linear solve type (for Newton direction)
cvode.max_order = 4 # CVODE max BDF order.

#------------------------- EB SETUP -----------------------------
eb2.geom_type = sphere
eb2.sphere_radius = 0.018
eb2.sphere_center = 0.0 0.0 0.0
eb2.sphere_has_fluid_inside = 1
eb2.small_volfrac = 1.e-4

#--------------------REFINEMENT CONTROL------------------------
#amr.refinement_indicators = temp
#amr.temp.max_level = 1
#amr.temp.value_greater = 305
#amr.temp.field_name = temp

#amr.refinement_indicators = magVort
#amr.magVort.max_level = 1
#amr.magVort.value_greater = 500.0
#amr.magVort.field_name = mag_vort

#--------------------LINEAR SOLVER CONTROL------------------------
nodal_proj.verbose = 1
nodal_proj.rtol = 1.0e-10
nodal_proj.atol = 1.0e-10
nodal_proj.mg_max_coarsening_level = 2

#fabarray.mfiter_tile_size = 1024 1024 1024

amrex.fpe_trap_invalid = 1
amrex.fpe_trap_zero = 1
amrex.fpe_trap_overflow = 1
105 changes: 80 additions & 25 deletions Source/PeleLMAdvection.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -330,7 +330,10 @@ void PeleLM::computeScalarAdvTerms(std::unique_ptr<AdvanceAdvData> &advData)

//----------------------------------------------------------------
#ifdef AMREX_USE_EB
// Get EBFact & areafrac
const auto& ebfact = EBFactory(lev);
Array< const MultiCutFab*,AMREX_SPACEDIM> areafrac;
areafrac = ebfact.getAreaFrac();
#endif

// Get the species edge state and advection term
Expand Down Expand Up @@ -434,19 +437,44 @@ void PeleLM::computeScalarAdvTerms(std::unique_ptr<AdvanceAdvData> &advData)

Box const& bx = mfi.tilebox();

#ifdef AMREX_USE_EB
auto const& flagfab = ebfact.getMultiEBCellFlagFab()[mfi];
#endif

// Edge states
for (int idim = 0; idim <AMREX_SPACEDIM; idim++) {
const Box& ebx = amrex::surroundingNodes(bx,idim);
auto const& rho_ed = edgeState[idim].array(mfi,0);
auto const& rhoY_ed = edgeState[idim].array(mfi,1);
amrex::ParallelFor(ebx, [rho_ed, rhoY_ed]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
rho_ed(i,j,k) = 0.0;
for (int n = 0; n < NUM_SPECIES; n++) {
rho_ed(i,j,k) += rhoY_ed(i,j,k,n);
}
});
const Box& ebx = amrex::surroundingNodes(bx,idim);
auto const& rho_ed = edgeState[idim].array(mfi,0);
auto const& rhoY_ed = edgeState[idim].array(mfi,1);
#ifdef AMREX_USE_EB
if (flagfab.getType(ebx) == FabType::covered) { // Covered boxes
amrex::ParallelFor(ebx, [rho_ed]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
rho_ed(i,j,k) = 0.0;
});
} else if (flagfab.getType(ebx) != FabType::regular ) { // EB containing boxes
const auto& afrac = areafrac[idim]->array(mfi);
amrex::ParallelFor(ebx, [rho_ed, rhoY_ed, afrac]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
rho_ed(i,j,k) = 0.0;
if (afrac(i,j,k) > 0.0) { // Uncovered faces
for (int n = 0; n < NUM_SPECIES; n++) {
rho_ed(i,j,k) += rhoY_ed(i,j,k,n);
}
}
});
} else // Regular boxes
#endif
amrex::ParallelFor(ebx, [rho_ed, rhoY_ed]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
rho_ed(i,j,k) = 0.0;
for (int n = 0; n < NUM_SPECIES; n++) {
rho_ed(i,j,k) += rhoY_ed(i,j,k,n);
}
});
}
}

Expand Down Expand Up @@ -488,6 +516,7 @@ void PeleLM::computeScalarAdvTerms(std::unique_ptr<AdvanceAdvData> &advData)
m_advection_type);
}


// Get the edge RhoH states
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
Expand All @@ -496,21 +525,43 @@ void PeleLM::computeScalarAdvTerms(std::unique_ptr<AdvanceAdvData> &advData)

Box const& bx = mfi.tilebox();

#ifdef AMREX_USE_EB
auto const& flagfab = ebfact.getMultiEBCellFlagFab()[mfi];
#endif

for (int idim = 0; idim <AMREX_SPACEDIM; idim++) {
const Box& ebx = amrex::surroundingNodes(bx,idim);
auto const& rho = edgeState[idim].const_array(mfi,0);
auto const& rhoY = edgeState[idim].const_array(mfi,1);
auto const& T = edgeState[idim].const_array(mfi,NUM_SPECIES+2);
auto const& rhoHm = edgeState[idim].array(mfi,NUM_SPECIES+1);
amrex::ParallelFor(ebx, [rho, rhoY, T, rhoHm]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
if (rho(i,j,k) <= 0.0) { // Covered faces
rhoHm(i,j,k) = 0.0;
} else {
getRHmixGivenTY( i, j, k, rho, rhoY, T, rhoHm );
}
});
const Box& ebx = amrex::surroundingNodes(bx,idim);
auto const& rho = edgeState[idim].const_array(mfi,0);
auto const& rhoY = edgeState[idim].const_array(mfi,1);
auto const& T = edgeState[idim].const_array(mfi,NUM_SPECIES+2);
auto const& rhoHm = edgeState[idim].array(mfi,NUM_SPECIES+1);
#ifdef AMREX_USE_EB
if (flagfab.getType(ebx) == FabType::covered) { // Covered boxes
amrex::ParallelFor(ebx, [rhoHm]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
rhoHm(i,j,k) = 0.0;
});
} else if (flagfab.getType(ebx) != FabType::regular ) { // EB containing boxes
const auto& afrac = areafrac[idim]->array(mfi);
amrex::ParallelFor(ebx, [rho, rhoY, T, rhoHm, afrac]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
if (afrac(i,j,k) <= 0.0) { // Covered faces
rhoHm(i,j,k) = 0.0;
} else {
getRHmixGivenTY( i, j, k, rho, rhoY, T, rhoHm );
}
});
} else // Regular boxes
#endif
{
amrex::ParallelFor(ebx, [rho, rhoY, T, rhoHm]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
getRHmixGivenTY( i, j, k, rho, rhoY, T, rhoHm );
});
}
}
}

Expand Down Expand Up @@ -552,6 +603,9 @@ void PeleLM::computeScalarAdvTerms(std::unique_ptr<AdvanceAdvData> &advData)
is_velocity, fluxes_are_area_weighted,
m_advection_type);
}
#ifdef AMREX_USE_EB
EB_set_covered_faces(GetArrOfPtrs(fluxes[lev]),0.);
#endif
}

//----------------------------------------------------------------
Expand Down Expand Up @@ -625,6 +679,7 @@ void PeleLM::computeScalarAdvTerms(std::unique_ptr<AdvanceAdvData> &advData)
1,
bcRecRhoH_d.dataPtr(),
geom[lev]);
EB_set_covered(advData->AofS[lev],0.0);
#else
//----------------------------------------------------------------
// Otherwise go directly into AofS
Expand Down
15 changes: 15 additions & 0 deletions Source/PeleLMDiffusion.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -134,6 +134,21 @@ void PeleLM::computeDifferentialDiffusionTerms(const TimeStamp &a_time,
fluxDivergence(GetVecOfPtrs(diffData->Dwbar), 0, GetVecOfArrOfPtrs(diffData->wbar_fluxes), 0, NUM_SPECIES, 1, -1.0);
#endif
}

#ifdef AMREX_USE_EB
// Set EB-covered diffusion terms here to avoid having to mask operations using this data
// later on
for (int lev = 0; lev <= finest_level; ++lev) {
if (a_time == AmrOldTime) {
EB_set_covered(diffData->Dn[lev],0.0);
} else {
EB_set_covered(diffData->Dnp1[lev],0.0);
}
if (!is_init && m_use_wbar) {
EB_set_covered(diffData->Dwbar[lev],0.0);
}
}
#endif
}

void PeleLM::computeDifferentialDiffusionFluxes(const TimeStamp &a_time,
Expand Down
42 changes: 38 additions & 4 deletions Source/PeleLMEos.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -84,12 +84,24 @@ void PeleLM::calcDivU(int is_init,
}
}

//----------------------------------------------------------------
#ifdef AMREX_USE_EB
// Get EBFact
const auto& ebfact = EBFactory(lev);
#endif

#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (MFIter mfi(ldata_p->divu, TilingIfNotGPU()); mfi.isValid(); ++mfi)
{
const Box& bx = mfi.tilebox();

#ifdef AMREX_USE_EB
auto const& flagfab = ebfact.getMultiEBCellFlagFab()[mfi];
auto const& flag = flagfab.const_array();
#endif

auto const& rhoY = ldata_p->species.const_array(mfi);
auto const& T = ldata_p->temp.const_array(mfi);
auto const& SpecD = ( a_time == AmrOldTime ) ? diffData->Dn[lev].const_array(mfi,0)
Expand All @@ -101,11 +113,33 @@ void PeleLM::calcDivU(int is_init,
auto const& r = (m_do_react) ? RhoYdot.const_array(mfi)
: ldata_p->species.const_array(mfi); // Dummy unused Array4
auto const& divu = ldata_p->divu.array(mfi);
amrex::ParallelFor(bx, [ rhoY, T, SpecD, Fourier, DiffDiff, r, divu, do_react=m_do_react]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept

#ifdef AMREX_USE_EB
if (flagfab.getType(bx) == FabType::covered) { // Covered boxes
amrex::ParallelFor(bx, [divu]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
divu(i,j,k) = 0.0;
});
} else if (flagfab.getType(bx) != FabType::regular ) { // EB containing boxes
amrex::ParallelFor(bx, [ rhoY, T, SpecD, Fourier, DiffDiff, r, divu, do_react=m_do_react, flag]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
if ( flag(i,j,k).isCovered() ) {
divu(i,j,k) = 0.0;
} else {
compute_divu( i, j, k, rhoY, T, SpecD, Fourier, DiffDiff, r, divu, do_react );
}
});
} else
#endif
{
compute_divu( i, j, k, rhoY, T, SpecD, Fourier, DiffDiff, r, divu, do_react );
});
amrex::ParallelFor(bx, [ rhoY, T, SpecD, Fourier, DiffDiff, r, divu, do_react=m_do_react]
AMREX_GPU_DEVICE (int i, int j, int k) noexcept
{
compute_divu( i, j, k, rhoY, T, SpecD, Fourier, DiffDiff, r, divu, do_react );
});
}
}
}

Expand Down
38 changes: 21 additions & 17 deletions Source/PeleLMUtils.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -345,23 +345,27 @@ void PeleLM::advFluxDivergence(int a_lev,
}
#ifdef AMREX_USE_EB
else {
AMREX_D_TERM(auto const& apx_arr = ebfact.getAreaFrac()[0]->const_array(mfi);,
auto const& apy_arr = ebfact.getAreaFrac()[1]->const_array(mfi);,
auto const& apz_arr = ebfact.getAreaFrac()[2]->const_array(mfi););
ParallelFor(bx, ncomp, [=]
AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
{
if (!l_conserv_d[n] && vfrac_arr(i,j,k) > 0.) {
Real qwsum = AMREX_D_TERM( apx_arr(i,j,k) * facex(i,j,k,n) + apx_arr(i+1,j,k) * facex(i+1,j,k,n),
+ apy_arr(i,j,k) * facey(i,j,k,n) + apy_arr(i,j+1,k) * facey(i,j+1,k,n),
+ apz_arr(i,j,k) * facez(i,j,k,n) + apz_arr(i,j,k+1) * facez(i,j,k+1,n));
Real areasum = AMREX_D_TERM( apx_arr(i,j,k) + apx_arr(i+1,j,k),
+ apy_arr(i,j,k) + apy_arr(i,j+1,k),
+ apz_arr(i,j,k) + apz_arr(i,j,k+1));
// Note that because we define adv update as MINUS div(u q), here we add q div (u)
div_arr(i,j,k,n) += qwsum / areasum * divu_arr(i,j,k);
}
});
if (flagfab.getType(bx) == FabType::covered) {
AMREX_PARALLEL_FOR_4D(bx, ncomp, i, j, k, n, {div_arr(i,j,k,n) = 0.0;});
} else {
AMREX_D_TERM(auto const& apx_arr = ebfact.getAreaFrac()[0]->const_array(mfi);,
auto const& apy_arr = ebfact.getAreaFrac()[1]->const_array(mfi);,
auto const& apz_arr = ebfact.getAreaFrac()[2]->const_array(mfi););
ParallelFor(bx, ncomp, [=]
AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
{
if (!l_conserv_d[n] && vfrac_arr(i,j,k) > 0.) {
Real qwsum = AMREX_D_TERM( apx_arr(i,j,k) * facex(i,j,k,n) + apx_arr(i+1,j,k) * facex(i+1,j,k,n),
+ apy_arr(i,j,k) * facey(i,j,k,n) + apy_arr(i,j+1,k) * facey(i,j+1,k,n),
+ apz_arr(i,j,k) * facez(i,j,k,n) + apz_arr(i,j,k+1) * facez(i,j,k+1,n));
Real areasum = AMREX_D_TERM( apx_arr(i,j,k) + apx_arr(i+1,j,k),
+ apy_arr(i,j,k) + apy_arr(i,j+1,k),
+ apz_arr(i,j,k) + apz_arr(i,j,k+1));
// Note that because we define adv update as MINUS div(u q), here we add q div (u)
div_arr(i,j,k,n) += qwsum / areasum * divu_arr(i,j,k);
}
});
}
}
#endif
}
Expand Down