version 3.9
brinkman.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_BRINKMAN_TERM_HH
13#define DUMUX_NAVIERSTOKES_MOMENTUM_BRINKMAN_TERM_HH
14
15#include <type_traits>
16#include <dune/common/typetraits.hh>
17#include <dune/common/exceptions.hh>
18
19#include <dumux/common/math.hh>
24
25namespace Dumux {
26
45template<class NumEqVector, class Problem, class FVElementGeometry, class ElementVolumeVariables>
47 NumEqVector& source,
48 const Problem& problem,
49 const typename FVElementGeometry::Element& element,
50 const FVElementGeometry& fvGeometry,
51 const ElementVolumeVariables& elemVolVars,
52 const typename FVElementGeometry::SubControlVolume& scv
53){
54 if (problem.spatialParams().brinkmanEpsilon(element, fvGeometry, scv) > 0.0)
55 {
56 const auto brinkmanEpsilon = problem.spatialParams().brinkmanEpsilon(element, fvGeometry, scv);
57 const auto invK = problem.spatialParams().inversePermeability(element, fvGeometry, scv);
58 const auto mu = problem.effectiveViscosity(element, fvGeometry, scv);
59
60 using DM = typename FVElementGeometry::GridGeometry::DiscretizationMethod;
61 if constexpr (DiscretizationMethods::isCVFE<DM>)
62 {
63 const auto velocity = elemVolVars[scv].velocity();
64 source -= mu * brinkmanEpsilon * mv(invK, velocity);
65 }
66 else if constexpr (DM{} == DiscretizationMethods::fcstaggered)
67 {
68 if constexpr (Dune::IsNumber<std::decay_t<decltype(invK)>>::value)
69 {
70 const auto velocity = elemVolVars[scv].velocity();
71 source -= mu * brinkmanEpsilon * invK * velocity;
72 }
73 else
74 {
75 // permeability is tensor-valued, use velocity vector reconstruction
76 const auto getVelocitySCV = [&](const auto& scv){ return elemVolVars[scv].velocity(); };
77 const auto velocity = StaggeredVelocityReconstruction::faceVelocity(scv, fvGeometry, getVelocitySCV);
78 source -= mu * brinkmanEpsilon * mv(invK, velocity)[scv.dofAxis()];
79 }
80 }
81 else
82 DUNE_THROW(Dune::NotImplemented, "Brinkman term only implemented for CVFE and Staggered");
83 }
84}
85
86} // end namespace Dumux
87
88#endif
The local element solution class for control-volume finite element methods.
free functions for the evaluation of primary variables inside elements.
Dune::DenseVector< V >::derived_type mv(const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V > &v)
Returns the result of the projection of a vector v with a Matrix M.
Definition: math.hh:829
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
Define some often used mathematical functions.
The available discretization methods in Dumux.
constexpr FCStaggered fcstaggered
Definition: method.hh:151
Definition: adapt.hh:17
void addBrinkmanTerm(NumEqVector &source, const Problem &problem, const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolume &scv)
Definition: brinkman.hh:46
static auto faceVelocity(const SubControlVolume &scv, const FVElementGeometry &fvGeometry, const VelocityHelper &getNormalVelocityDofValue)
Return the velocity vector at dof position given an scv.
Definition: velocityreconstruction.hh:52