Skip to content
Snippets Groups Projects
Commit 2e7c73b9 authored by Ned Coltman's avatar Ned Coltman
Browse files

[tpfa] Deprecate the effective law methods for the tpfa flux laws

parent 936c3dea
No related branches found
No related tags found
1 merge request!1684Improve effective laws
......@@ -185,23 +185,21 @@ public:
const auto insideScvIdx = scvf.insideScvIdx();
const auto& insideScv = fvGeometry.scv(insideScvIdx);
const auto& insideVolVars = elemVolVars[insideScvIdx];
const auto insideD = getEffectiveDiffusionCoefficient_(insideVolVars, phaseIdx, compIdx); //, Dune::index_constant<hasEffDiff>{});
const auto insideD = insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx);
const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideD, insideVolVars.extrusionFactor());
// for the boundary (dirichlet) or at branching points we only need ti
if (scvf.boundary() || scvf.numOutsideScvs() > 1)
{
tij = scvf.area()*ti;
}
// otherwise we compute a tpfa harmonic mean
else
{
const auto outsideScvIdx = scvf.outsideScvIdx();
const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
const auto& outsideVolVars = elemVolVars[outsideScvIdx];
const auto outsideD = outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx);
const auto outsideD = getEffectiveDiffusionCoefficient_(outsideVolVars, phaseIdx, compIdx); //, Dune::index_constant<hasEffDiff>{});
Scalar tj;
if (dim == dimWorld)
......@@ -268,6 +266,21 @@ private:
}
return rho/(scvf.numOutsideScvs()+1);
}
static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars, const int phaseIdx, const int compIdx)
{
if constexpr (Dumux::Deprecated::hasEffDiffCoeff<VolumeVariables>)
return volVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx);
else
{
// TODO: remove this else clause after release 3.2!
using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
return EffDiffModel::effectiveDiffusivity(volVars.porosity(),
volVars.saturation(phaseIdx),
volVars.diffusionCoefficient(phaseIdx, compIdx));
}
}
};
} // end namespace Dumux
......
......@@ -144,7 +144,7 @@ public:
const auto& insideScv = fvGeometry.scv(insideScvIdx);
const auto& insideVolVars = elemVolVars[insideScvIdx];
const auto insideLambda = insideVolVars.effectiveThermalConductivity();
const auto insideLambda = getEffectiveThermalConductivity_(insideVolVars, problem, element, fvGeometry, insideScv);
const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
// for the boundary (dirichlet) or at branching points we only need ti
......@@ -159,7 +159,7 @@ public:
const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
const auto& outsideVolVars = elemVolVars[outsideScvIdx];
const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
const auto outsideLambda = getEffectiveThermalConductivity_(outsideVolVars, problem, element, fvGeometry, outsideScv);
Scalar tj;
if (dim == dimWorld)
// assume the normal vector from outside is anti parallel so we save flipping a vector
......@@ -205,6 +205,27 @@ private:
}
return sumTempTi/sumTi;
}
template <class VolumeVariables>
static Scalar getEffectiveThermalConductivity_(const VolumeVariables& volVars,
const Problem& problem,
const Element& element,
const FVElementGeometry& fvGeometry,
const SubControlVolume& scv)
{
if constexpr (Dumux::Deprecated::hasEffTherCond<VolumeVariables>)
return volVars.effectiveThermalConductivity();
else
{
// TODO: remove this fallback after release 3.2!
return Deprecated::template effectiveThermalConductivity<ThermalConductivityModel>(volVars,
problem.spatialParams(),
element,
fvGeometry,
scv);
}
}
};
} // end namespace Dumux
......
......@@ -116,53 +116,37 @@ public:
if (phaseIdx != sPhaseIdx)
{
//when number of energyEq for the fluid are smaller than numPhases that means that we need an effecitve law
if (numEnergyEqFluid < ModelTraits::numFluidPhases())
{
insideLambda += insideVolVars.effectiveThermalConductivity();
}
else //numEnergyEqFluid >1
{
if (numEnergyEqFluid < ModelTraits::numFluidPhases())
insideLambda += getEffectiveThermalConductivity_(insideVolVars, problem, element, fvGeometry, insideScv);
else // numEnergyEqFluid >1
insideLambda += insideVolVars.fluidThermalConductivity(phaseIdx)*insideVolVars.saturation(phaseIdx)*insideVolVars.porosity();
}
}
//solid phase
else
{
else // solid phase
insideLambda += insideVolVars.solidThermalConductivity()*(1.0-insideVolVars.porosity());
}
const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
// for the boundary (dirichlet) or at branching points we only need ti
if (scvf.boundary() || scvf.numOutsideScvs() > 1)
{
tij = scvf.area()*ti;
}
// otherwise we compute a tpfa harmonic mean
else
else // otherwise we compute a tpfa harmonic mean
{
const auto outsideScvIdx = scvf.outsideScvIdx();
const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
const auto& outsideVolVars = elemVolVars[outsideScvIdx];
// effective diffusion tensors
if (phaseIdx != sPhaseIdx)
{
//when number of energyEq for the fluid are smaller than numPhases that means that we need an effecitve law
if (numEnergyEqFluid < ModelTraits::numFluidPhases())
// effective diffusion tensors
if (phaseIdx != sPhaseIdx)
{
outsideLambda += outsideVolVars.effectiveThermalConductivity();
//when number of energyEq for the fluid are smaller than numPhases that means that we need an effecitve law
if (numEnergyEqFluid < ModelTraits::numFluidPhases())
outsideLambda += getEffectiveThermalConductivity_(outsideVolVars, problem, element, fvGeometry, outsideScv);
else
outsideLambda += outsideVolVars.fluidThermalConductivity(phaseIdx)*outsideVolVars.saturation(phaseIdx)*outsideVolVars.porosity();
}
else
{
outsideLambda += outsideVolVars.fluidThermalConductivity(phaseIdx)*outsideVolVars.saturation(phaseIdx)*outsideVolVars.porosity();
}
}
//solid phase
else
{
outsideLambda +=outsideVolVars.solidThermalConductivity()*(1.0-outsideVolVars.porosity());
}
else //solid phase
outsideLambda +=outsideVolVars.solidThermalConductivity()*(1.0-outsideVolVars.porosity());
Scalar tj;
if (dim == dimWorld)
// assume the normal vector from outside is anti parallel so we save flipping a vector
......@@ -178,6 +162,27 @@ public:
}
return tij;
}
private:
template <class VolumeVariables>
static Scalar getEffectiveThermalConductivity_(const VolumeVariables& volVars,
const Problem& problem,
const Element& element,
const FVElementGeometry& fvGeometry,
const SubControlVolume& scv)
{
if constexpr (Dumux::Deprecated::hasEffTherCond<VolumeVariables>)
return volVars.effectiveThermalConductivity();
else
{
// TODO: remove this fallback after release 3.2!
return Deprecated::template effectiveThermalConductivity<ThermalConductivityModel>(volVars,
problem.spatialParams(),
element,
fvGeometry,
scv);
}
}
};
} // end namespace Dumux
......
......@@ -229,7 +229,7 @@ private:
const auto xi = volVars.moleFraction(phaseIdx, compIIdx);
const auto Mn = FluidSystem::molarMass(numComponents-1);
Scalar tin = volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1);
Scalar tin = getEffectiveDiffusionCoefficient_(volVars, phaseIdx, compIIdx, numComponents-1, problem, element, scv);
// set the entries of the diffusion matrix of the diagonal
reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
......@@ -244,7 +244,7 @@ private:
const auto xj = volVars.moleFraction(phaseIdx, compJIdx);
const auto Mi = FluidSystem::molarMass(compIIdx);
const auto Mj = FluidSystem::molarMass(compJIdx);
Scalar tij = volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
Scalar tij = getEffectiveDiffusionCoefficient_(volVars, phaseIdx, compIIdx, compJIdx, problem, element, scv);
reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
if (compJIdx < numComponents-1)
......@@ -254,6 +254,61 @@ private:
return reducedDiffusionMatrix;
}
static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars,
const int phaseIdx,
const int compIIdx,
const int compJIdx,
const Problem& problem,
const Element& element,
const SubControlVolume& scv)
{
if constexpr (Dumux::Deprecated::hasEffDiffCoeff<VolumeVariables>)
return volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
else
{
// TODO: remove this else clause after release 3.2!// effective diffusion tensors
using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
const auto tinInside = getDiffusionCoefficient(phaseIdx, compIIdx, compJIdx, problem, element, volVars, scv);
return EffDiffModel::effectiveDiffusivity(volVars.porosity(), volVars.saturation(phaseIdx), tinInside);
}
}
template <class T = TypeTag, typename std::enable_if_t<GetPropType<T, Properties::FluidSystem>::isTracerFluidSystem(), int> =0 >
[[deprecated("Signature deprecated. Will be removed after 3.2!")]]
static Scalar getDiffusionCoefficient(const int phaseIdx,
const int compIIdx,
const int compJIdx,
const Problem& problem,
const Element& element,
const VolumeVariables& volVars,
const SubControlVolume& scv)
{
return FluidSystem::binaryDiffusionCoefficient(compIIdx,
compJIdx,
problem,
element,
scv);
}
template <class T = TypeTag, typename std::enable_if_t<!GetPropType<T, Properties::FluidSystem>::isTracerFluidSystem(), int> =0 >
[[deprecated("Signature deprecated. Will be removed after 3.2!")]]
static Scalar getDiffusionCoefficient(const int phaseIdx,
const int compIIdx,
const int compJIdx,
const Problem& problem,
const Element& element,
const VolumeVariables& volVars,
const SubControlVolume& scv)
{
auto fluidState = volVars.fluidState();
typename FluidSystem::ParameterCache paramCache;
paramCache.updateAll(fluidState);
return FluidSystem::binaryDiffusionCoefficient(fluidState,
paramCache,
phaseIdx,
compIIdx,
compJIdx);
}
};
} // end namespace
......
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