Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup NavierStokesModel | ||
10 | * \brief A helper function to add Brinkman term to the momentum balance | ||
11 | */ | ||
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> | ||
20 | #include <dumux/discretization/method.hh> | ||
21 | #include <dumux/discretization/evalsolution.hh> | ||
22 | #include <dumux/discretization/cvfe/elementsolution.hh> | ||
23 | #include <dumux/freeflow/navierstokes/momentum/velocityreconstruction.hh> | ||
24 | |||
25 | namespace Dumux { | ||
26 | |||
27 | /*! | ||
28 | * \file | ||
29 | * \ingroup NavierStokesModel | ||
30 | * \brief A helper function to add the Brinkman term to the momentum balance | ||
31 | * \addtogroup NavierStokesModel | ||
32 | * @{ | ||
33 | * The Navier-Stokes model can be extended to a Darcy-Brinkman model by adding | ||
34 | * the term: | ||
35 | * \f[ | ||
36 | + \epsilon_B \mu \mathbf{K}^{-1} \mathbf{v} | ||
37 | * \f] | ||
38 | * to the momentum balance. This can be achieved with the helper function \ref addBrinkmanTerm. | ||
39 | * This function relies on the spatial parameters class being based on <code>BrinkmanSpatialParams</code> | ||
40 | * or providing the `brinkmanEpsilon` and `inversePermeability` interfaces. | ||
41 | * These interface functions provide the | ||
42 | * weighting factor \f$ \epsilon_B \f$ and the permeability tensor \f$ \mathbf{K} \f$. | ||
43 | * @} | ||
44 | */ | ||
45 | template<class NumEqVector, class Problem, class FVElementGeometry, class ElementVolumeVariables> | ||
46 | 1178720 | void addBrinkmanTerm( | |
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 |
2/2✓ Branch 2 taken 334844 times.
✓ Branch 3 taken 843876 times.
|
2357440 | if (problem.spatialParams().brinkmanEpsilon(element, fvGeometry, scv) > 0.0) |
55 | { | ||
56 | 669688 | const auto brinkmanEpsilon = problem.spatialParams().brinkmanEpsilon(element, fvGeometry, scv); | |
57 | 669688 | const auto invK = problem.spatialParams().inversePermeability(element, fvGeometry, scv); | |
58 | 334844 | const auto mu = problem.effectiveViscosity(element, fvGeometry, scv); | |
59 | |||
60 | using DM = typename FVElementGeometry::GridGeometry::DiscretizationMethod; | ||
61 | if constexpr (DiscretizationMethods::isCVFE<DM>) | ||
62 | { | ||
63 | 244500 | const auto velocity = elemVolVars[scv].velocity(); | |
64 | 733500 | 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 | 271032 | const auto getVelocitySCV = [&](const auto& scv){ return elemVolVars[scv].velocity(); }; | |
77 | 90344 | const auto velocity = StaggeredVelocityReconstruction::faceVelocity(scv, fvGeometry, getVelocitySCV); | |
78 | 180688 | 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 | 1178720 | } | |
85 | |||
86 | } // end namespace Dumux | ||
87 | |||
88 | #endif | ||
89 |