mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
blackoil, ebos: implement non-trivial boundary conditions
with this, it is possible do define fluxes of conservation quantities over the domain boundaries by specifying the thermodynamic state on the boundary when using the black-oil model. The main motivation is are thermal fluxes which are required to maintain geothermal temperature gradients over time.
This commit is contained in:
parent
815be1451b
commit
54c96aa1c2
@ -44,6 +44,7 @@ NEW_PROP_TAG(FluidSystem);
|
||||
NEW_PROP_TAG(GridView);
|
||||
NEW_PROP_TAG(Scalar);
|
||||
NEW_PROP_TAG(MaterialLaw);
|
||||
NEW_PROP_TAG(EnableEnergy);
|
||||
}
|
||||
|
||||
/*!
|
||||
@ -65,10 +66,6 @@ class EclEquilInitializer
|
||||
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
|
||||
typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
|
||||
|
||||
typedef Opm::BlackOilFluidState<Scalar,
|
||||
FluidSystem,
|
||||
/*enableTemperature=*/true> ScalarFluidState;
|
||||
|
||||
enum { numPhases = FluidSystem::numPhases };
|
||||
enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
|
||||
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
|
||||
@ -80,6 +77,12 @@ class EclEquilInitializer
|
||||
enum { waterCompIdx = FluidSystem::waterCompIdx };
|
||||
|
||||
enum { dimWorld = GridView::dimensionworld };
|
||||
enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
|
||||
|
||||
typedef Opm::BlackOilFluidState<Scalar,
|
||||
FluidSystem,
|
||||
/*enableTemperature=*/true,
|
||||
/*enableEnthalpy=*/enableEnergy> ScalarFluidState;
|
||||
|
||||
public:
|
||||
template <class EclMaterialLawManager>
|
||||
|
@ -349,12 +349,110 @@ protected:
|
||||
}
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Update the required gradients for boundary faces
|
||||
*/
|
||||
template <class FluidState>
|
||||
void calculateBoundaryGradients_(const ElementContext& elemCtx,
|
||||
unsigned scvfIdx,
|
||||
unsigned timeIdx,
|
||||
const FluidState& exFluidState)
|
||||
{
|
||||
bool enableBoundaryMassFlux = false;
|
||||
if (!enableBoundaryMassFlux)
|
||||
return;
|
||||
|
||||
const auto& problem = elemCtx.problem();
|
||||
const auto& stencil = elemCtx.stencil(timeIdx);
|
||||
const auto& scvf = stencil.boundaryFace(scvfIdx);
|
||||
|
||||
interiorDofIdx_ = scvf.interiorIndex();
|
||||
|
||||
Scalar trans = problem.transmissibilityBoundary(elemCtx, scvfIdx);
|
||||
Scalar faceArea = scvf.area();
|
||||
|
||||
// estimate the gravity correction: for performance reasons we use a simplified
|
||||
// approach for this flux module that assumes that gravity is constant and always
|
||||
// acts into the downwards direction. (i.e., no centrifuge experiments, sorry.)
|
||||
Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
|
||||
|
||||
const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
|
||||
|
||||
// this is quite hacky because the dune grid interface does not provide a
|
||||
// cellCenterDepth() method (so we ask the problem to provide it). The "good"
|
||||
// solution would be to take the Z coordinate of the element centroids, but since
|
||||
// ECL seems to like to be inconsistent on that front, it needs to be done like
|
||||
// here...
|
||||
Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx_, timeIdx);
|
||||
Scalar zEx = scvf.integrationPos()[dimWorld - 1];
|
||||
|
||||
// the distances from the DOF's depths. (i.e., the additional depth of the
|
||||
// exterior DOF)
|
||||
Scalar distZ = zIn - zEx;
|
||||
|
||||
for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
|
||||
if (!FluidSystem::phaseIsActive(phaseIdx))
|
||||
continue;
|
||||
|
||||
// do the gravity correction: compute the hydrostatic pressure for the
|
||||
// integration position
|
||||
const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
|
||||
const auto& rhoEx = exFluidState.density(phaseIdx);
|
||||
Evaluation rhoAvg = (rhoIn + rhoEx)/2;
|
||||
|
||||
const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
|
||||
Evaluation pressureExterior = exFluidState.pressure(phaseIdx);
|
||||
pressureExterior += rhoAvg*(distZ*g);
|
||||
|
||||
pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
|
||||
|
||||
// decide the upstream index for the phase. for this we make sure that the
|
||||
// degree of freedom which is regarded upstream if both pressures are equal
|
||||
// is always the same: if the pressure is equal, the DOF with the lower
|
||||
// global index is regarded to be the upstream one.
|
||||
if (pressureDifference_[phaseIdx] > 0.0) {
|
||||
upIdx_[phaseIdx] = -1;
|
||||
dnIdx_[phaseIdx] = interiorDofIdx_;
|
||||
}
|
||||
else {
|
||||
upIdx_[phaseIdx] = interiorDofIdx_;
|
||||
dnIdx_[phaseIdx] = -1;
|
||||
}
|
||||
|
||||
// this is slightly hacky because in the automatic differentiation case, it
|
||||
// only works for the element centered finite volume method. for ebos this
|
||||
// does not matter, though.
|
||||
unsigned upstreamIdx = upstreamIndex_(phaseIdx);
|
||||
const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
|
||||
if (upstreamIdx == interiorDofIdx_)
|
||||
volumeFlux_[phaseIdx] =
|
||||
pressureDifference_[phaseIdx]*up.mobility(phaseIdx)*(-trans/faceArea);
|
||||
else {
|
||||
// compute the phase mobility using the material law parameters of the
|
||||
// interior element. TODO: this could probably be done more efficiently
|
||||
const auto& matParams =
|
||||
elemCtx.problem().materialLawParams(elemCtx,
|
||||
interiorDofIdx_,
|
||||
/*timeIdx=*/0);
|
||||
typename FluidState::Scalar kr[numPhases];
|
||||
MaterialLaw::relativePermeabilities(kr, matParams, exFluidState);
|
||||
|
||||
const auto& mob = kr[phaseIdx]/exFluidState.viscosity(phaseIdx);
|
||||
volumeFlux_[phaseIdx] =
|
||||
pressureDifference_[phaseIdx]*mob*(-trans/faceArea);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Update the volumetric fluxes for all fluid phases on the interior faces of the context
|
||||
*/
|
||||
void calculateFluxes_(const ElementContext& elemCtx OPM_UNUSED, unsigned scvfIdx OPM_UNUSED, unsigned timeIdx OPM_UNUSED)
|
||||
{ }
|
||||
|
||||
void calculateBoundaryFluxes_(const ElementContext& elemCtx OPM_UNUSED, unsigned scvfIdx OPM_UNUSED, unsigned timeIdx OPM_UNUSED)
|
||||
{}
|
||||
|
||||
private:
|
||||
Implementation& asImp_()
|
||||
{ return *static_cast<Implementation*>(this); }
|
||||
|
@ -150,7 +150,7 @@ public:
|
||||
typedef typename EclMaterialLawManager::MaterialLaw type;
|
||||
};
|
||||
|
||||
// Set the material law for heat heat storage in rock
|
||||
// Set the material law for energy storage in rock
|
||||
SET_PROP(EclBaseProblem, SolidEnergyLaw)
|
||||
{
|
||||
private:
|
||||
@ -163,7 +163,7 @@ public:
|
||||
typedef typename EclThermalLawManager::SolidEnergyLaw type;
|
||||
};
|
||||
|
||||
// Set the material law for heat conduction
|
||||
// Set the material law for thermal conduction
|
||||
SET_PROP(EclBaseProblem, ThermalConductionLaw)
|
||||
{
|
||||
private:
|
||||
@ -324,8 +324,9 @@ class EclProblem : public GET_PROP_TYPE(TypeTag, BaseProblem)
|
||||
typedef BlackOilPolymerModule<TypeTag> PolymerModule;
|
||||
|
||||
typedef Opm::BlackOilFluidState<Scalar,
|
||||
FluidSystem,
|
||||
/*enableTemperature=*/true> InitialFluidState;
|
||||
FluidSystem,
|
||||
/*enableTemperature=*/true,
|
||||
/*enableEnthalpy=*/enableEnergy> InitialFluidState;
|
||||
|
||||
typedef Opm::MathToolbox<Evaluation> Toolbox;
|
||||
typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
|
||||
@ -768,6 +769,17 @@ public:
|
||||
return pffDofData_.get(context.element(), toDofLocalIdx).transmissibility;
|
||||
}
|
||||
|
||||
/*!
|
||||
* \copydoc EclTransmissiblity::transmissibilityBoundary
|
||||
*/
|
||||
template <class Context>
|
||||
Scalar transmissibilityBoundary(const Context& elemCtx,
|
||||
unsigned boundaryFaceIdx) const
|
||||
{
|
||||
unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
|
||||
return transmissibilities_.transmissibilityBoundary(elemIdx, boundaryFaceIdx);
|
||||
}
|
||||
|
||||
/*!
|
||||
* \copydoc EclTransmissiblity::thermalHalfTransmissibility
|
||||
*/
|
||||
@ -784,11 +796,10 @@ public:
|
||||
* \copydoc EclTransmissiblity::thermalHalfTransmissibility
|
||||
*/
|
||||
template <class Context>
|
||||
Scalar thermalHalfTransmissibilityBoundary(const Context& context,
|
||||
unsigned dofIdx,
|
||||
Scalar thermalHalfTransmissibilityBoundary(const Context& elemCtx,
|
||||
unsigned boundaryFaceIdx) const
|
||||
{
|
||||
unsigned elemIdx = context.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
||||
unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
|
||||
return transmissibilities_.thermalHalfTransBoundary(elemIdx, boundaryFaceIdx);
|
||||
}
|
||||
|
||||
@ -888,7 +899,7 @@ public:
|
||||
{ return materialLawManager_->materialLawParams(globalDofIdx); }
|
||||
|
||||
/*!
|
||||
* \brief Return the parameters for the heat storage law of the rock
|
||||
* \brief Return the parameters for the energy storage law of the rock
|
||||
*/
|
||||
template <class Context>
|
||||
const SolidEnergyLawParams&
|
||||
@ -1069,10 +1080,23 @@ public:
|
||||
*/
|
||||
template <class Context>
|
||||
void boundary(BoundaryRateVector& values,
|
||||
const Context& context OPM_UNUSED,
|
||||
unsigned spaceIdx OPM_UNUSED,
|
||||
unsigned timeIdx OPM_UNUSED) const
|
||||
{ values.setNoFlow(); }
|
||||
const Context& context,
|
||||
unsigned spaceIdx,
|
||||
unsigned timeIdx) const
|
||||
{
|
||||
if (!enableEnergy)
|
||||
values.setNoFlow();
|
||||
else {
|
||||
// in the energy case we need to specify a non-trivial boundary condition
|
||||
// because the geothermal gradient needs to be maintained. for this, we
|
||||
// simply assume the initial temperature at the boundary and specify the
|
||||
// thermal flow accordingly. in this context, "thermal flow" means energy
|
||||
// flow due to a temerature gradient while assuming no-flow for mass
|
||||
unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
|
||||
unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
|
||||
values.setThermalFlow(context, spaceIdx, timeIdx, initialFluidStates_[globalDofIdx]);
|
||||
}
|
||||
}
|
||||
|
||||
/*!
|
||||
* \copydoc FvBaseProblem::initial
|
||||
|
@ -162,10 +162,14 @@ public:
|
||||
trans_.clear();
|
||||
trans_.reserve(numElements*3*1.05);
|
||||
|
||||
transBoundary_.clear();
|
||||
|
||||
// if energy is enabled, let's do the same for the "thermal half transmissibilities"
|
||||
if (enableEnergy) {
|
||||
thermalHalfTrans_->clear();
|
||||
thermalHalfTrans_->reserve(numElements*3*1.05);
|
||||
|
||||
thermalHalfTransBoundary_.clear();
|
||||
}
|
||||
|
||||
// compute the transmissibilities for all intersections
|
||||
@ -176,13 +180,54 @@ public:
|
||||
|
||||
auto isIt = gridView.ibegin(elem);
|
||||
const auto& isEndIt = gridView.iend(elem);
|
||||
unsigned boundaryIsIdx = 0;
|
||||
for (; isIt != isEndIt; ++ isIt) {
|
||||
// store intersection, this might be costly
|
||||
const auto& intersection = *isIt;
|
||||
|
||||
// ignore boundary intersections for now (TODO?)
|
||||
if (!intersection.neighbor())
|
||||
// deal with grid boundaries
|
||||
if (!intersection.neighbor()) {
|
||||
// compute the transmissibilty for the boundary intersection
|
||||
const auto& geometry = intersection.geometry();
|
||||
const auto& faceCenterInside = geometry.center();
|
||||
|
||||
auto faceAreaNormal = intersection.centerUnitOuterNormal();
|
||||
faceAreaNormal *= geometry.volume();
|
||||
|
||||
Scalar transBoundaryIs;
|
||||
computeHalfTrans_(transBoundaryIs,
|
||||
faceAreaNormal,
|
||||
intersection.indexInInside(),
|
||||
distanceVector_(faceCenterInside,
|
||||
intersection.indexInInside(),
|
||||
elemIdx,
|
||||
axisCentroids),
|
||||
permeability_[elemIdx]);
|
||||
|
||||
// normally there would be two half-transmissibilities that would be
|
||||
// averaged. on the grid boundary there only is the half
|
||||
// transmissibility of the interior element.
|
||||
transBoundary_[std::make_pair(elemIdx, boundaryIsIdx)] = transBoundaryIs;
|
||||
|
||||
// for boundary intersections we also need to compute the thermal
|
||||
// half transmissibilities
|
||||
if (enableEnergy) {
|
||||
const auto& n = intersection.centerUnitOuterNormal();
|
||||
Scalar A = intersection.geometry().volume();
|
||||
|
||||
const auto& inPos = elem.geometry().center();
|
||||
const auto& outPos = intersection.geometry().center();
|
||||
const auto& d = outPos - inPos;
|
||||
|
||||
Scalar thermalHalfTrans =
|
||||
A * std::abs(n*d)/(d*d);
|
||||
|
||||
thermalHalfTransBoundary_[std::make_pair(elemIdx, boundaryIsIdx)] = thermalHalfTrans;
|
||||
}
|
||||
|
||||
++ boundaryIsIdx;
|
||||
continue;
|
||||
}
|
||||
|
||||
const auto& outsideElem = intersection.outside();
|
||||
unsigned outsideElemIdx = elemMapper.index(outsideElem);
|
||||
@ -217,10 +262,16 @@ public:
|
||||
DimVector faceCenterOutside;
|
||||
DimVector faceAreaNormal;
|
||||
|
||||
typename std::is_same< Grid, Dune::CpGrid> :: type isCpGrid;
|
||||
computeFaceProperties( intersection, elemIdx, insideFaceIdx, outsideElemIdx, outsideFaceIdx,
|
||||
faceCenterInside, faceCenterOutside, faceAreaNormal,
|
||||
isCpGrid );
|
||||
typename std::is_same<Grid, Dune::CpGrid>::type isCpGrid;
|
||||
computeFaceProperties(intersection,
|
||||
elemIdx,
|
||||
insideFaceIdx,
|
||||
outsideElemIdx,
|
||||
outsideFaceIdx,
|
||||
faceCenterInside,
|
||||
faceCenterOutside,
|
||||
faceAreaNormal,
|
||||
isCpGrid);
|
||||
|
||||
Scalar halfTrans1;
|
||||
Scalar halfTrans2;
|
||||
@ -303,6 +354,12 @@ public:
|
||||
Scalar transmissibility(unsigned elemIdx1, unsigned elemIdx2) const
|
||||
{ return trans_.at(isId_(elemIdx1, elemIdx2)); }
|
||||
|
||||
/*!
|
||||
* \brief Return the transmissibility for a given boundary segment.
|
||||
*/
|
||||
Scalar transmissibilityBoundary(unsigned elemIdx, unsigned boundaryFaceIdx) const
|
||||
{ return transBoundary_.at(std::make_pair(elemIdx, boundaryFaceIdx)); }
|
||||
|
||||
/*!
|
||||
* \brief Return the thermal "half transmissibility" for the intersection between two
|
||||
* elements.
|
||||
@ -320,6 +377,9 @@ public:
|
||||
Scalar thermalHalfTrans(unsigned insideElemIdx, unsigned outsideElemIdx) const
|
||||
{ return thermalHalfTrans_->at(isId_(insideElemIdx, outsideElemIdx)); }
|
||||
|
||||
Scalar thermalHalfTransBoundary(unsigned insideElemIdx, unsigned boundaryFaceIdx) const
|
||||
{ return thermalHalfTransBoundary_.at(std::make_pair(insideElemIdx, boundaryFaceIdx)); }
|
||||
|
||||
private:
|
||||
template <class Intersection>
|
||||
void computeFaceProperties( const Intersection& intersection,
|
||||
@ -330,7 +390,7 @@ private:
|
||||
DimVector& faceCenterInside,
|
||||
DimVector& faceCenterOutside,
|
||||
DimVector& faceAreaNormal,
|
||||
std::false_type ) const
|
||||
/*isCpGrid=*/std::false_type) const
|
||||
{
|
||||
// default implementation for DUNE grids
|
||||
const auto& geometry = intersection.geometry();
|
||||
@ -350,7 +410,7 @@ private:
|
||||
DimVector& faceCenterInside,
|
||||
DimVector& faceCenterOutside,
|
||||
DimVector& faceAreaNormal,
|
||||
std::true_type ) const
|
||||
/*isCpGrid=*/std::true_type) const
|
||||
{
|
||||
int faceIdx = intersection.id();
|
||||
faceCenterInside = vanguard_.grid().faceCenterEcl(insideElemIdx,insideFaceIdx);
|
||||
@ -537,6 +597,8 @@ private:
|
||||
const Vanguard& vanguard_;
|
||||
std::vector<DimMatrix> permeability_;
|
||||
std::unordered_map<std::uint64_t, Scalar> trans_;
|
||||
std::map<std::pair<unsigned, unsigned>, Scalar> transBoundary_;
|
||||
std::map<std::pair<unsigned, unsigned>, Scalar> thermalHalfTransBoundary_;
|
||||
Opm::ConditionalStorage<enableEnergy,
|
||||
std::unordered_map<std::uint64_t, Scalar> > thermalHalfTrans_;
|
||||
};
|
||||
|
Loading…
Reference in New Issue
Block a user