Skip to content
Snippets Groups Projects
Commit b56ad1ff authored by Timo Koch's avatar Timo Koch
Browse files

[richards] Make analytical robin interface optional

parent 313724a0
No related branches found
No related tags found
1 merge request!2574Fix/richards analytic storage derivative
Pipeline #3781 waiting for manual action
......@@ -34,6 +34,18 @@
#include <dumux/discretization/extrusion.hh>
#include <dumux/flux/referencesystemformulation.hh>
namespace Dumux::Detail {
// helper structs and functions detecting if the user-defined problem class
// implements addRobinFluxDerivatives
template <typename T, typename ...Ts>
using RobinDerivDetector = decltype(std::declval<T>().addRobinFluxDerivatives(std::declval<Ts>()...));
template<class T, typename ...Args>
static constexpr bool hasAddRobinFluxDerivatives()
{ return Dune::Std::is_detected<RobinDerivDetector, T, Args...>::value; }
} // end namespace Dumux::Detail
namespace Dumux {
/*!
......@@ -486,9 +498,11 @@ public:
const ElementFluxVariablesCache& elemFluxVarsCache,
const SubControlVolumeFace& scvf) const
{
/* forward to problem for the user to implement the Robin derivatives*/
problem.addRobinFluxDerivatives(derivativeMatrices, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
if constexpr(Detail::hasAddRobinFluxDerivatives<Problem,
PartialDerivativeMatrices, Element, FVElementGeometry,
ElementVolumeVariables, ElementFluxVariablesCache, SubControlVolumeFace>()
)
problem.addRobinFluxDerivatives(derivativeMatrices, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
}
private:
......
......@@ -73,8 +73,6 @@ class RichardsLensProblem : public PorousMediumFlowProblem<TypeTag>
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
using FVElementGeometry = typename GridGeometry::LocalView;
using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
enum {
// copy some indices for convenience
pressureIdx = Indices::pressureIdx,
......@@ -194,25 +192,6 @@ public:
// \}
/*!
* \brief Adds Robin flux derivatives for wetting phase. This is needed in case of solution dependent Neumann conditions.
*
* \param derivativeMatrices The matrices containing the derivatives
* \param element The element
* \param fvGeometry The finite volume element geometry
* \param curElemVolVars The current element volume variables
* \param elemFluxVarsCache The element flux variables cache
* \param scvf The sub control volume face
*/
template<class PartialDerivativeMatrices, class ElementVolumeVariables, class ElementFluxVariablesCache>
void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
const Element& element,
const FVElementGeometry& fvGeometry,
const ElementVolumeVariables& elemVolVars,
const ElementFluxVariablesCache& elemFluxVarsCache,
const SubControlVolumeFace& scvf) const
{}
private:
PrimaryVariables initial_(const GlobalPosition &globalPos) const
{
......
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