Skip to content
Snippets Groups Projects
Commit cb28e955 authored by Dennis Gläser's avatar Dennis Gläser
Browse files

neralization] simplify kozeny-karman law

This commit also removes porosity and permeability laws from the
1pncmin spatial parameters as they were not used anyway.
parent 54d13d38
No related branches found
No related tags found
1 merge request!889Feature/simplify perm poro laws
......@@ -35,71 +35,27 @@ namespace Dumux {
* \brief The Kozeny-Carman relationship for the calculation of a porosity-dependent permeability.
* When the porosity is implemented as solution-independent, using this relationship for the
* permeability leads to unnecessary overhead.
*
* \tparam PermeabilityType The type used for the intrinsic permeability
*/
template<class TypeTag>
template<class PermeabilityType>
class PermeabilityKozenyCarman
{
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
using SpatialParams = typename GET_PROP_TYPE(TypeTag, SpatialParams);
using VolumeVariables = typename GET_PROP_TYPE(TypeTag, VolumeVariables);
using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::LocalView;
using SubControlVolume = typename FVElementGeometry::SubControlVolume;
static constexpr int dim = GridView::dimension;
static constexpr int dimWorld = GridView::dimensionworld;
using Element = typename GridView::template Codim<0>::Entity;
using Tensor = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
using PermType = typename SpatialParams::PermeabilityType;
public:
/*!
* \brief The initial parameter distribution.
* \param spatialParams the spatial parameters
*/
void init(const SpatialParams& spatialParams)
{
spatialParamsPtr_ = &spatialParams;
}
/*!
* \brief calculates the permeability for a given sub-control volume
* \param element element
* \param elemSol the element solution
* \param scv sub control volume
* \param refPerm Reference permeability before porosity changes
* \param refPoro The poro corresponding to the reference permeability
* \param poro The porosity for which permeability is to be evaluated
*/
template<class ElementSolution>
PermType evaluatePermeability(const Element& element,
const SubControlVolume& scv,
const ElementSolution& elemSol) const
template<class Scalar>
PermeabilityType evaluatePermeability(PermeabilityType refPerm, Scalar refPoro, Scalar poro) const
{
auto refPoro = spatialParams().referencePorosity(element, scv);
auto poro = spatialParams().porosity(element, scv, elemSol);
using std::pow;
auto factor = pow((1.0 - refPoro)/(1.0 - poro), 2) * pow(poro/refPoro, 3);
return applyFactorToPermeability_(spatialParams().referencePermeability(element, scv), factor);
}
private:
const SpatialParams& spatialParams() const
{ return *spatialParamsPtr_; }
Scalar applyFactorToPermeability_(Scalar k, Scalar factor) const
{ return k*factor; }
Tensor applyFactorToPermeability_(Tensor K, Scalar factor) const
{
Tensor result(K);
for (int i = 0; i < dim; ++i)
result[i][i] *= factor;
return result;
refPerm *= factor;
return refPerm;
}
const SpatialParams* spatialParamsPtr_;
};
} // namespace Dumux
......
......@@ -28,8 +28,6 @@
#include <dumux/material/spatialparams/fv1p.hh>
#include <dumux/material/fluidmatrixinteractions/porosityreactivebed.hh>
#include <dumux/material/fluidmatrixinteractions/permeabilitykozenycarman.hh>
#include <dumux/material/fluidmatrixinteractions/mineralization/effectivesoliddensity.hh>
#include <dumux/material/fluidmatrixinteractions/mineralization/effectivesolidheatcapacity.hh>
......@@ -79,8 +77,6 @@ class ThermoChemSpatialParams : public FVSpatialParamsOneP<TypeTag>
using SubControlVolume = typename FVElementGeometry::SubControlVolume;
using Element = typename GridView::template Codim<0>::Entity;
using PorosityLaw = PorosityReactiveBed<TypeTag>;
using PermeabilityLaw = PermeabilityKozenyCarman<TypeTag>;
using EffectiveSolidRho = EffectiveSolidDensity<TypeTag>;
using EffectiveSolidCp = EffectiveSolidHeatCapacity<TypeTag>;
......@@ -102,32 +98,10 @@ public:
cp_[cPhaseIdx-numPhases] = 934; //[J/kgK] heat capacity of CaO (see Shao et al. 2014)
cp_[hPhaseIdx-numPhases] = 1530; //[J/kgK] heat capacity of Ca(OH)_2 (see Shao et al. 2014)
eps_ = 1e-6;
poroLaw_.init(*this);
permLaw_.init(*this);
effSolRho_.init(*this);
effSolCp_.init(*this);
}
/*!
* \brief Define the reference permeability \f$[m^2]\f$ distribution
*
* \param element The finite element
* \param scv The sub-control volume
*/
Scalar referencePermeability(const Element& element, const SubControlVolume &scv) const
{ return 8.53e-12; }
/*!
* \brief Define the reference porosity \f$[-]\f$ distribution
*
* \param element The finite element
* \param scv The sub-control volume
*/
Scalar referencePorosity(const Element& element, const SubControlVolume &scv) const
{
return 0.8;
}
/*! Intrinsic permeability tensor K \f$[m^2]\f$ depending
* on the position in the domain
*
......@@ -141,19 +115,7 @@ public:
Scalar permeability(const Element& element,
const SubControlVolume& scv,
const ElementSolution& elemSol) const
{ return permLaw_.evaluatePermeability(element, scv, elemSol); }
/*!
* \brief Define the minimum porosity \f$[-]\f$ distribution
*
* \param element The finite element
* \param scv The sub-control volume
*/
Scalar minPorosity(const Element& element, const SubControlVolume &scv) const
{
// return 0.604; //intrinsic porosity of CaO2H2 (see Nagel et al. 2014)
return 0.8;
}
{ return 8.53e-12; }
/*!
* \brief Define the minimum porosity \f$[-]\f$ after clogging caused by mineralization
......@@ -166,10 +128,7 @@ public:
Scalar porosity(const Element& element,
const SubControlVolume& scv,
const ElementSolution& elemSol) const
{
// return poroLaw_.evaluatePorosity(element, scv, elemSol);
return 0.8;
}
{ return 0.8; }
/*!
* \brief Returns the average heat capacity \f$[J / (kg K)]\f$ of solid phases.
......@@ -248,8 +207,6 @@ private:
Scalar lambdaSolid_;
std::array<Scalar, numSPhases> rho_;
std::array<Scalar, numSPhases> cp_;
PorosityLaw poroLaw_;
PermeabilityLaw permLaw_;
EffectiveSolidRho effSolRho_;
EffectiveSolidCp effSolCp_;
};
......
......@@ -83,10 +83,6 @@ class DissolutionSpatialparams : public FVSpatialParams<TypeTag>
using SubControlVolume = typename FVElementGeometry::SubControlVolume;
using Element = typename GridView::template Codim<0>::Entity;
using ModelTraits = typename GET_PROP_TYPE(TypeTag, ModelTraits);
using PorosityLaw = PorosityPrecipitation<Scalar, ModelTraits::numComponents(), ModelTraits::numSPhases()>;
using PermeabilityLaw = PermeabilityKozenyCarman<TypeTag>;
public:
// type used for the permeability (i.e. tensor or scalar)
using PermeabilityType = Scalar;
......@@ -109,9 +105,6 @@ public:
// parameters of Brooks & Corey Law
materialParams_.setPe(pEntry1_);
materialParams_.setLambda(bcLambda1_);
//! Initialize the parameter laws
permLaw_.init(*this);
}
/*!
......@@ -147,16 +140,6 @@ public:
const ElementSolution& elemSol) const
{ return poroLaw_.evaluatePorosity(element, scv, elemSol, referencePorosity_, 1e-5); }
/*!
* \brief Define the reference permeability \f$[m^2]\f$ distribution
*
* \param element The finite element
* \param scv The sub-control volume
*/
PermeabilityType referencePermeability(const Element& element, const SubControlVolume &scv) const
{ return referencePermeability_; }
/*! Intrinsic permeability tensor K \f$[m^2]\f$ depending
* on the position in the domain
*
......@@ -169,7 +152,10 @@ public:
PermeabilityType permeability(const Element& element,
const SubControlVolume& scv,
const ElementSolution& elemSol) const
{ return permLaw_.evaluatePermeability(element, scv, elemSol); }
{
const auto poro = porosity(element, scv, elemSol);
return permLaw_.evaluatePermeability(referencePermeability_, referencePorosity_, poro);
}
Scalar solidity(const SubControlVolume &scv) const
{ return 1.0 - porosityAtPos(scv.center()); }
......@@ -188,8 +174,11 @@ private:
MaterialLawParams materialParams_;
using ModelTraits = typename GET_PROP_TYPE(TypeTag, ModelTraits);
using PorosityLaw = PorosityPrecipitation<Scalar, ModelTraits::numComponents(), ModelTraits::numSPhases()>;
PorosityLaw poroLaw_;
PermeabilityLaw permLaw_;
PermeabilityKozenyCarman<PermeabilityType> permLaw_;
Scalar solubilityLimit_;
Scalar referencePorosity_;
PermeabilityType referencePermeability_ = 0.0;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment