Skip to content

Commit

Permalink
Multi-component diffusion (#39)
Browse files Browse the repository at this point in the history
* Add a Multi-Component diffusionOp to class.

* add ncomp to diffusionOp constructor, default to 1

* getMCDiffusionOp function. if ncomp != m_ncomp, reset the operator.

* Reset mcdiffusionOp during regrid.

* Enable DiffusionOp::diffuse_scalar, DiffusionOp::computeDiffFluxes and DiffusionOp::computeDiffLap
with multiple components.

* Switch to getMCDiffusionOp when dealing with species.

* NewStateRedist was removed from AMReX-Hydro
  • Loading branch information
esclapez authored and nickwimer committed Mar 24, 2022
1 parent 99e0410 commit c42dc3e
Show file tree
Hide file tree
Showing 6 changed files with 60 additions and 45 deletions.
3 changes: 2 additions & 1 deletion Source/DiffusionOp.H
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,7 @@ class DiffusionOp

// Methods

DiffusionOp (PeleLM* a_pelelm);
DiffusionOp (PeleLM* a_pelelm, int ncomp = 1);

void diffuse_scalar (amrex::Vector<amrex::MultiFab*> const& a_phi,
int phi_comp,
Expand Down Expand Up @@ -92,6 +92,7 @@ class DiffusionOp
#endif

int m_verbose = 0;
int m_ncomp = 1;

// Options to control MLMG behavior
int m_mg_verbose = 0;
Expand Down
54 changes: 30 additions & 24 deletions Source/DiffusionOp.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -13,8 +13,8 @@ using namespace amrex;
//---------------------------------------------------------------------------------------
// Diffusion Operator

DiffusionOp::DiffusionOp (PeleLM* a_pelelm)
: m_pelelm(a_pelelm)
DiffusionOp::DiffusionOp (PeleLM* a_pelelm, int ncomp)
: m_pelelm(a_pelelm), m_ncomp(ncomp)
{
readParameters();

Expand Down Expand Up @@ -42,12 +42,12 @@ DiffusionOp::DiffusionOp (PeleLM* a_pelelm)
m_scal_apply_op.reset(new MLEBABecLap(m_pelelm->Geom(0,m_pelelm->finestLevel()),
m_pelelm->boxArray(0,m_pelelm->finestLevel()),
m_pelelm->DistributionMap(0,m_pelelm->finestLevel()),
info_apply, ebfactVec));
info_apply, ebfactVec, m_ncomp));
#else
m_scal_apply_op.reset(new MLABecLaplacian(m_pelelm->Geom(0,m_pelelm->finestLevel()),
m_pelelm->boxArray(0,m_pelelm->finestLevel()),
m_pelelm->DistributionMap(0,m_pelelm->finestLevel()),
info_apply));
info_apply, {}, m_ncomp));
#endif
m_scal_apply_op->setMaxOrder(m_mg_maxorder);

Expand All @@ -56,12 +56,12 @@ DiffusionOp::DiffusionOp (PeleLM* a_pelelm)
m_scal_solve_op.reset(new MLEBABecLap(m_pelelm->Geom(0,m_pelelm->finestLevel()),
m_pelelm->boxArray(0,m_pelelm->finestLevel()),
m_pelelm->DistributionMap(0,m_pelelm->finestLevel()),
info_solve, ebfactVec));
info_solve, ebfactVec, m_ncomp));
#else
m_scal_solve_op.reset(new MLABecLaplacian(m_pelelm->Geom(0,m_pelelm->finestLevel()),
m_pelelm->boxArray(0,m_pelelm->finestLevel()),
m_pelelm->DistributionMap(0,m_pelelm->finestLevel()),
info_solve));
info_solve, {}, m_ncomp));
#endif
m_scal_solve_op->setMaxOrder(m_mg_maxorder);

Expand All @@ -70,12 +70,12 @@ DiffusionOp::DiffusionOp (PeleLM* a_pelelm)
m_gradient_op.reset(new MLEBABecLap(m_pelelm->Geom(0,m_pelelm->finestLevel()),
m_pelelm->boxArray(0,m_pelelm->finestLevel()),
m_pelelm->DistributionMap(0,m_pelelm->finestLevel()),
info_apply, ebfactVec));
info_apply, ebfactVec, m_ncomp));
#else
m_gradient_op.reset(new MLABecLaplacian(m_pelelm->Geom(0,m_pelelm->finestLevel()),
m_pelelm->boxArray(0,m_pelelm->finestLevel()),
m_pelelm->DistributionMap(0,m_pelelm->finestLevel()),
info_apply));
info_apply, {}, m_ncomp));
#endif
m_gradient_op->setMaxOrder(m_mg_maxorder);
m_gradient_op->setScalars(0.0,1.0);
Expand Down Expand Up @@ -107,6 +107,7 @@ void DiffusionOp::diffuse_scalar(Vector<MultiFab*> const& a_phi, int phi_comp,

//----------------------------------------------------------------
// Checks
AMREX_ASSERT(m_ncomp == 1 || m_ncomp == ncomp);
AMREX_ASSERT(a_phi[0]->nComp() >= phi_comp+ncomp);
AMREX_ASSERT(a_rhs[0]->nComp() >= rhs_comp+ncomp);
if ( have_fluxes ) {
Expand Down Expand Up @@ -168,8 +169,8 @@ void DiffusionOp::diffuse_scalar(Vector<MultiFab*> const& a_phi, int phi_comp,
}

//----------------------------------------------------------------
// Solve and get fluxes on a per component basis
for (int comp = 0; comp < ncomp; ++comp) {
// Solve and get fluxes on a m_ncomp component basis
for (int comp = 0; comp < ncomp; comp+=m_ncomp) {

// Aliases
Vector<Array<MultiFab*,AMREX_SPACEDIM>> fluxes(finest_level+1);
Expand All @@ -184,12 +185,13 @@ void DiffusionOp::diffuse_scalar(Vector<MultiFab*> const& a_phi, int phi_comp,
for (int lev = 0; lev <= finest_level; ++lev) {
if (have_fluxes) {
for (int idim = 0; idim < AMREX_SPACEDIM; idim++ ) {
fluxes[lev][idim] = new MultiFab(*a_flux[lev][idim],amrex::make_alias,flux_comp+comp,1);
fluxes[lev][idim] = new MultiFab(*a_flux[lev][idim],amrex::make_alias,flux_comp+comp,m_ncomp);
}
}

if (have_bcoeff) {
Array<MultiFab,AMREX_SPACEDIM> bcoeff_ec = m_pelelm->getDiffusivity(lev, bcoeff_comp+comp, 1, {a_bcrec[comp]}, *a_bcoeff[lev]);
Vector<BCRec> subBCRec = {a_bcrec.begin()+comp,a_bcrec.begin()+comp+m_ncomp};
Array<MultiFab,AMREX_SPACEDIM> bcoeff_ec = m_pelelm->getDiffusivity(lev, bcoeff_comp+comp, m_ncomp, subBCRec, *a_bcoeff[lev]);
#ifdef AMREX_USE_EB
m_scal_solve_op->setBCoeffs(lev, GetArrOfConstPtrs(bcoeff_ec), MLMG::Location::FaceCentroid);
#else
Expand All @@ -199,8 +201,8 @@ void DiffusionOp::diffuse_scalar(Vector<MultiFab*> const& a_phi, int phi_comp,
m_scal_solve_op->setBCoeffs(lev, 1.0);
}

component.emplace_back(phi[lev],amrex::make_alias,comp,1);
rhs.emplace_back(*a_rhs[lev],amrex::make_alias,rhs_comp+comp,1);
component.emplace_back(phi[lev],amrex::make_alias,comp,m_ncomp);
rhs.emplace_back(*a_rhs[lev],amrex::make_alias,rhs_comp+comp,m_ncomp);
m_scal_solve_op->setLevelBC(lev, &component[lev]);
}

Expand Down Expand Up @@ -276,6 +278,7 @@ void DiffusionOp::computeDiffLap(Vector<MultiFab*> const& a_laps, int lap_comp,

//----------------------------------------------------------------
// Checks
AMREX_ASSERT(m_ncomp == 1 || m_ncomp == ncomp);
AMREX_ASSERT(a_laps[0]->nComp() >= lap_comp+ncomp);
AMREX_ASSERT(a_phi[0]->nComp() >= phi_comp+ncomp);
AMREX_ASSERT(a_bcoeff[0]->nComp() >= bcoeff_comp+ncomp);
Expand All @@ -300,16 +303,17 @@ void DiffusionOp::computeDiffLap(Vector<MultiFab*> const& a_laps, int lap_comp,
Real beta = -1.0;
m_scal_apply_op->setScalars(alpha, beta);

for (int comp = 0; comp < ncomp; ++comp) {
for (int comp = 0; comp < ncomp; comp+=m_ncomp) {

// Component based vector of data
Vector<MultiFab> laps;
Vector<MultiFab> component;

for (int lev = 0; lev <= finest_level; ++lev) {
laps.emplace_back(*a_laps[lev],amrex::make_alias,lap_comp+comp,1);
component.emplace_back(phi[lev],amrex::make_alias,comp,1);
Array<MultiFab,AMREX_SPACEDIM> bcoeff_ec = m_pelelm->getDiffusivity(lev, bcoeff_comp+comp, 1, {a_bcrec[comp]}, *a_bcoeff[lev]);
laps.emplace_back(*a_laps[lev],amrex::make_alias,lap_comp+comp,m_ncomp);
component.emplace_back(phi[lev],amrex::make_alias,comp,m_ncomp);
Vector<BCRec> subBCRec = {a_bcrec.begin()+comp,a_bcrec.begin()+comp+m_ncomp};
Array<MultiFab,AMREX_SPACEDIM> bcoeff_ec = m_pelelm->getDiffusivity(lev, bcoeff_comp+comp, m_ncomp, subBCRec, *a_bcoeff[lev]);

#ifdef AMREX_USE_EB
m_scal_apply_op->setBCoeffs(lev, GetArrOfConstPtrs(bcoeff_ec), MLMG::Location::FaceCentroid);
Expand Down Expand Up @@ -337,6 +341,7 @@ void DiffusionOp::computeDiffFluxes(Vector<Array<MultiFab*,AMREX_SPACEDIM>> cons

//----------------------------------------------------------------
// Checks
AMREX_ASSERT(m_ncomp == 1 || m_ncomp == ncomp);
AMREX_ASSERT(a_flux[0][0]->nComp() >= flux_comp+ncomp);
AMREX_ASSERT(a_phi[0]->nComp() >= phi_comp+ncomp);
AMREX_ASSERT(a_bcoeff[0]->nComp() >= bcoeff_comp+ncomp);
Expand Down Expand Up @@ -382,8 +387,8 @@ void DiffusionOp::computeDiffFluxes(Vector<Array<MultiFab*,AMREX_SPACEDIM>> cons
Real beta = -1.0;
m_scal_apply_op->setScalars(alpha, beta);

// Get fluxes on a per component basis
for (int comp = 0; comp < ncomp; ++comp) {
// Get fluxes on a m_ncomp component(s) basis
for (int comp = 0; comp < ncomp; comp+=m_ncomp) {

// Component based vector of data
Vector<Array<MultiFab*,AMREX_SPACEDIM>> fluxes(finest_level+1);
Expand All @@ -396,12 +401,13 @@ void DiffusionOp::computeDiffFluxes(Vector<Array<MultiFab*,AMREX_SPACEDIM>> cons

for (int lev = 0; lev <= finest_level; ++lev) {
for (int idim = 0; idim < AMREX_SPACEDIM; idim++ ) {
fluxes[lev][idim] = new MultiFab(*a_flux[lev][idim],amrex::make_alias,flux_comp+comp,1);
fluxes[lev][idim] = new MultiFab(*a_flux[lev][idim],amrex::make_alias,flux_comp+comp,m_ncomp);
}
component.emplace_back(phi[lev],amrex::make_alias,comp,1);
Array<MultiFab,AMREX_SPACEDIM> bcoeff_ec = m_pelelm->getDiffusivity(lev, bcoeff_comp+comp, 1, {a_bcrec[comp]}, *a_bcoeff[lev]);
component.emplace_back(phi[lev],amrex::make_alias,comp,m_ncomp);
Vector<BCRec> subBCRec = {a_bcrec.begin()+comp,a_bcrec.begin()+comp+m_ncomp};
Array<MultiFab,AMREX_SPACEDIM> bcoeff_ec = m_pelelm->getDiffusivity(lev, bcoeff_comp+comp, m_ncomp, subBCRec, *a_bcoeff[lev]);
laps.emplace_back(a_phi[lev]->boxArray(), a_phi[lev]->DistributionMap(),
1, 1, MFInfo(), a_phi[lev]->Factory());
m_ncomp, 1, MFInfo(), a_phi[lev]->Factory());
#ifdef AMREX_USE_EB
m_scal_apply_op->setBCoeffs(lev, GetArrOfConstPtrs(bcoeff_ec), MLMG::Location::FaceCentroid);
#else
Expand Down
4 changes: 3 additions & 1 deletion Source/PeleLM.H
Original file line number Diff line number Diff line change
Expand Up @@ -818,6 +818,7 @@ class PeleLM : public amrex::AmrCore {
amrex::Vector<int> fetchAdvTypeArray(int scomp, int ncomp);
amrex::Vector<int> fetchDiffTypeArray(int scomp, int ncomp);
DiffusionOp* getDiffusionOp ();
DiffusionOp* getMCDiffusionOp (int ncomp = 1);
DiffusionTensorOp* getDiffusionTensorOp ();
//-----------------------------------------------------------------------------

Expand Down Expand Up @@ -885,6 +886,7 @@ class PeleLM : public amrex::AmrCore {

// Linear solvers
std::unique_ptr<DiffusionOp> m_diffusion_op;
std::unique_ptr<DiffusionOp> m_mcdiffusion_op;
std::unique_ptr<DiffusionTensorOp> m_diffusionTensor_op;
std::unique_ptr<Hydro::MacProjector> macproj;
int m_macProjNeedReset = 0;
Expand Down Expand Up @@ -1022,7 +1024,7 @@ class PeleLM : public amrex::AmrCore {
int m_refine_cutcells = 1;
amrex::Vector<amrex::Real> coveredState_h;
amrex::Gpu::DeviceVector<amrex::Real> coveredState_d;
std::string m_adv_redist_type = "NewStateRedist";
std::string m_adv_redist_type = "StateRedist";
std::string m_diff_redist_type = "FluxRedist";
#endif
//-----------------------------------------------------------------------------
Expand Down
31 changes: 19 additions & 12 deletions Source/PeleLMDiffusion.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,13 @@ PeleLM::getDiffusionOp()
return m_diffusion_op.get();
}

DiffusionOp*
PeleLM::getMCDiffusionOp(int ncomp)
{
if (!m_mcdiffusion_op || m_mcdiffusion_op->m_ncomp != ncomp) m_mcdiffusion_op.reset(new DiffusionOp(this,ncomp));
return m_mcdiffusion_op.get();
}

DiffusionTensorOp*
PeleLM::getDiffusionTensorOp ()
{
Expand Down Expand Up @@ -165,11 +172,11 @@ void PeleLM::computeDifferentialDiffusionFluxes(const TimeStamp &a_time,
// Get the species diffusion fluxes from the DiffusionOp
// Don't average down just yet
int do_avgDown = 0;
getDiffusionOp()->computeDiffFluxes(a_fluxes, 0,
GetVecOfConstPtrs(getSpeciesVect(a_time)), 0,
GetVecOfConstPtrs(getDensityVect(a_time)),
GetVecOfConstPtrs(getDiffusivityVect(a_time)), 0, bcRecSpec,
NUM_SPECIES, -1.0, do_avgDown);
getMCDiffusionOp(NUM_SPECIES)->computeDiffFluxes(a_fluxes, 0,
GetVecOfConstPtrs(getSpeciesVect(a_time)), 0,
GetVecOfConstPtrs(getDensityVect(a_time)),
GetVecOfConstPtrs(getDiffusivityVect(a_time)), 0, bcRecSpec,
NUM_SPECIES, -1.0, do_avgDown);

// Add the wbar term
// TODO: might need to do an average_down of the wbar fluxes
Expand Down Expand Up @@ -603,13 +610,13 @@ void PeleLM::differentialDiffusionUpdate(std::unique_ptr<AdvanceAdvData> &advDat
// Solve for \widetilda{rhoY^{np1,kp1}}
// -> return the uncorrected fluxes^{np1,kp1}
// -> and the partially updated species (not including wbar or flux correction)
getDiffusionOp()->diffuse_scalar(getSpeciesVect(AmrNewTime), 0,
GetVecOfConstPtrs(advData->Forcing), 0,
GetVecOfArrOfPtrs(fluxes), 0,
GetVecOfConstPtrs(getDensityVect(AmrNewTime)), // this is the acoeff of LinOp
GetVecOfConstPtrs(getDensityVect(AmrNewTime)), // this triggers proper scaling by density
GetVecOfConstPtrs(getDiffusivityVect(AmrNewTime)), 0, bcRecSpec,
NUM_SPECIES, 0, m_dt);
getMCDiffusionOp(NUM_SPECIES)->diffuse_scalar(getSpeciesVect(AmrNewTime), 0,
GetVecOfConstPtrs(advData->Forcing), 0,
GetVecOfArrOfPtrs(fluxes), 0,
GetVecOfConstPtrs(getDensityVect(AmrNewTime)), // this is the acoeff of LinOp
GetVecOfConstPtrs(getDensityVect(AmrNewTime)), // this triggers proper scaling by density
GetVecOfConstPtrs(getDiffusivityVect(AmrNewTime)), 0, bcRecSpec,
NUM_SPECIES, 0, m_dt);

// Add lagged Wbar term
// Computed in computeDifferentialDiffusionTerms at t^{n} if first SDC iteration, t^{np1,k} otherwise
Expand Down
10 changes: 3 additions & 7 deletions Source/PeleLMEB.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -71,8 +71,7 @@ void PeleLM::redistributeAofS(int a_lev,
// FluxRedistribute
Box gbx = bx;

if (m_adv_redist_type == "StateRedist" ||
m_adv_redist_type == "NewStateRedist")
if (m_adv_redist_type == "StateRedist")
gbx.grow(3);
else if (m_adv_redist_type == "FluxRedist")
gbx.grow(2);
Expand Down Expand Up @@ -144,8 +143,7 @@ void PeleLM::redistributeDiff(int a_lev,
// FluxRedistribute
Box gbx = bx;

if (m_diff_redist_type == "StateRedist" ||
m_diff_redist_type == "NewStateRedist")
if (m_diff_redist_type == "StateRedist")
gbx.grow(3);
else if (m_diff_redist_type == "FluxRedist")
gbx.grow(2);
Expand Down Expand Up @@ -236,9 +234,7 @@ void PeleLM::initialRedistribution()
{
// Redistribute the initial solution if adv/diff scheme uses State or NewState
if (m_adv_redist_type == "StateRedist" ||
m_adv_redist_type == "NewStateRedist" ||
m_diff_redist_type == "StateRedist" ||
m_diff_redist_type == "NewStateRedist") {
m_diff_redist_type == "StateRedist") {

for (int lev = 0; lev <= finest_level; ++lev) {

Expand Down
3 changes: 3 additions & 0 deletions Source/PeleLMRegrid.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -87,6 +87,7 @@ void PeleLM::MakeNewLevelFromCoarse( int lev,

// DiffusionOp will be recreated
m_diffusion_op.reset();
m_mcdiffusion_op.reset();
m_diffusionTensor_op.reset();

// Trigger MacProj reset
Expand Down Expand Up @@ -178,6 +179,7 @@ void PeleLM::RemakeLevel( int lev,

// DiffusionOp will be recreated
m_diffusion_op.reset();
m_mcdiffusion_op.reset();
m_diffusionTensor_op.reset();

// Trigger MacProj reset
Expand All @@ -195,6 +197,7 @@ void PeleLM::ClearLevel(int lev) {
m_dmapChem[lev].reset();
m_factory[lev].reset();
m_diffusion_op.reset();
m_mcdiffusion_op.reset();
m_diffusionTensor_op.reset();
macproj.reset();
#ifdef PELE_USE_EFIELD
Expand Down

0 comments on commit c42dc3e

Please sign in to comment.