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 BoundaryTests | ||
10 | * \brief The spatial parameters class for the test problem using the 1p cc model. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_1P_TEST_SPATIALPARAMS_HH | ||
14 | #define DUMUX_1P_TEST_SPATIALPARAMS_HH | ||
15 | |||
16 | #include <dumux/porousmediumflow/fvspatialparams1p.hh> | ||
17 | #include "testcase.hh" | ||
18 | |||
19 | namespace Dumux { | ||
20 | |||
21 | /*! | ||
22 | * \ingroup BoundaryTests | ||
23 | * \brief The spatial parameters class for the test problem using the | ||
24 | * 1p cc model. | ||
25 | */ | ||
26 | template<class GridGeometry, class Scalar> | ||
27 |
1/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
|
12 | class ConvergenceTestSpatialParams |
28 | : public FVPorousMediumFlowSpatialParamsOneP<GridGeometry, Scalar, | ||
29 | ConvergenceTestSpatialParams<GridGeometry, Scalar>> | ||
30 | { | ||
31 | using GridView = typename GridGeometry::GridView; | ||
32 | using ParentType = FVPorousMediumFlowSpatialParamsOneP< | ||
33 | GridGeometry, Scalar, ConvergenceTestSpatialParams<GridGeometry, Scalar> | ||
34 | >; | ||
35 | |||
36 | using Element = typename GridView::template Codim<0>::Entity; | ||
37 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
38 | |||
39 | static constexpr int dimWorld = GridView::dimensionworld; | ||
40 | using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>; | ||
41 | |||
42 | |||
43 | public: | ||
44 | // export permeability type | ||
45 | using PermeabilityType = DimWorldMatrix; | ||
46 | |||
47 | 12 | ConvergenceTestSpatialParams(std::shared_ptr<const GridGeometry> gridGeometry, | |
48 | const DarcyStokesTestCase testCase) | ||
49 | : ParentType(gridGeometry) | ||
50 |
2/8✓ Branch 2 taken 12 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
|
12 | , testCase_(testCase) |
51 | { | ||
52 |
1/2✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
|
12 | alphaBJ_ = getParam<Scalar>("Darcy.SpatialParams.AlphaBeaversJoseph"); |
53 | 12 | } | |
54 | |||
55 | /*! | ||
56 | * \brief Function for defining the (intrinsic) permeability \f$[m^2]\f$. | ||
57 | * | ||
58 | * \param element The element | ||
59 | * \param scv The sub control volume | ||
60 | * \param elemSol The element solution vector | ||
61 | * \return the intrinsic permeability | ||
62 | */ | ||
63 | template<class SubControlVolume, class ElementSolution> | ||
64 | 3742230 | PermeabilityType permeability(const Element& element, | |
65 | const SubControlVolume& scv, | ||
66 | const ElementSolution& elemSol) const | ||
67 | { | ||
68 | 3742230 | PermeabilityType K(0.0); | |
69 | |||
70 |
2/2✓ Branch 0 taken 1427530 times.
✓ Branch 1 taken 2314700 times.
|
3742230 | if (testCase_ == DarcyStokesTestCase::Schneider) |
71 | { | ||
72 | using std::cos; | ||
73 | using std::sin; | ||
74 | using std::exp; | ||
75 | |||
76 | static constexpr Scalar c = 0.0; | ||
77 | static constexpr Scalar omega = M_PI; | ||
78 | |||
79 | 2855060 | const Scalar x = scv.center()[0]; | |
80 | 2855060 | K[0][0] = 1.0; | |
81 | 2855060 | K[0][1] = -c/(2*omega) * sin(omega*x); | |
82 | 5710120 | K[1][0] = K[0][1]; | |
83 | 4282590 | K[1][1] = exp(-2)*(1 + c*cos(omega*x)); | |
84 | } | ||
85 | else | ||
86 | { | ||
87 |
4/6✓ Branch 0 taken 9 times.
✓ Branch 1 taken 2314691 times.
✓ Branch 3 taken 9 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 9 times.
✗ Branch 7 not taken.
|
2314700 | const static Scalar permeability = getParam<Scalar>("Darcy.SpatialParams.Permeability", 1.0); |
88 | 4629400 | K[0][0] = permeability; | |
89 | 6944100 | K[1][1] = permeability; | |
90 | } | ||
91 | |||
92 | 3742230 | return K; | |
93 | } | ||
94 | |||
95 | /*! \brief Defines the porosity in [-]. | ||
96 | * | ||
97 | * \param globalPos The global position | ||
98 | */ | ||
99 | ✗ | Scalar porosityAtPos(const GlobalPosition& globalPos) const | |
100 | ✗ | { return 0.4; } | |
101 | |||
102 | /*! \brief Defines the Beavers-Joseph coefficient in [-]. | ||
103 | * | ||
104 | * \param globalPos The global position | ||
105 | */ | ||
106 | ✗ | Scalar beaversJosephCoeffAtPos(const GlobalPosition& globalPos) const | |
107 | ✗ | { return alphaBJ_; } | |
108 | |||
109 | private: | ||
110 | DarcyStokesTestCase testCase_; | ||
111 | Scalar permeability_; | ||
112 | Scalar alphaBJ_; | ||
113 | }; | ||
114 | |||
115 | } // end namespace Dumux | ||
116 | |||
117 | #endif | ||
118 |