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 | /*! | ||
9 | * \file | ||
10 | * \ingroup NavierStokesTests | ||
11 | * \brief Darcy Brinkman model for a single-domain evaluation of coupled freeflow and porous medium flows | ||
12 | */ | ||
13 | |||
14 | #ifndef DUMUX_BRINKMAN_SPATIAL_PARAMS_HH | ||
15 | #define DUMUX_BRINKMAN_SPATIAL_PARAMS_HH | ||
16 | |||
17 | #include <cmath> | ||
18 | #include <numeric> | ||
19 | #include <algorithm> | ||
20 | #include <dune/common/fmatrix.hh> | ||
21 | #include <dumux/freeflow/spatialparams.hh> | ||
22 | |||
23 | namespace Dumux { | ||
24 | |||
25 | /*! | ||
26 | * \ingroup NavierStokesTests | ||
27 | * \brief The spatial parameters class for the Darcy-Brinkman model test | ||
28 | */ | ||
29 | template<class GridGeometry, class Scalar> | ||
30 |
2/12✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
|
4 | class BrinkmanTestSpatialParams |
31 | : public BrinkmanSpatialParams<GridGeometry, Scalar, BrinkmanTestSpatialParams<GridGeometry, Scalar>> | ||
32 | { | ||
33 | using ParentType = BrinkmanSpatialParams<GridGeometry, Scalar, BrinkmanTestSpatialParams<GridGeometry, Scalar>>; | ||
34 | using GridView = typename GridGeometry::GridView; | ||
35 | using Element = typename GridView::template Codim<0>::Entity; | ||
36 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
37 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
38 | |||
39 | static constexpr int dimWorld = GridView::dimensionworld; | ||
40 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
41 | using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>; | ||
42 | |||
43 | public: | ||
44 | using PermeabilityType = DimWorldMatrix; | ||
45 | |||
46 | 8 | BrinkmanTestSpatialParams(std::shared_ptr<const GridGeometry> gridGeometry) | |
47 | : ParentType(gridGeometry) | ||
48 | , permeability_(0.0) | ||
49 |
7/16✓ Branch 2 taken 4 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 4 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 4 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 4 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 4 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 4 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 4 times.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
|
64 | , inversePermeability_(0.0) |
50 | { | ||
51 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
8 | initPermeability_(); |
52 | |||
53 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
8 | pmLowerLeft_ = getParam<GlobalPosition>("SpatialParams.PorousMediumLowerLeft"); |
54 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
8 | pmUpperRight_ = getParam<GlobalPosition>("SpatialParams.PorousMediumUpperRight"); |
55 | 16 | pmCenter_ = 0.5*(pmLowerLeft_ + pmUpperRight_); | |
56 | 8 | pmDimensions_ = pmUpperRight_ - pmLowerLeft_; | |
57 |
1/4✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
8 | pmRoundedness_ = getParam<Scalar>("SpatialParams.PorousMediumRoundedness", 0.0); |
58 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
8 | pmTransitionHalfLength_ = 0.5*getParam<Scalar>("SpatialParams.PorousMediumTransitionLength"); |
59 | 8 | } | |
60 | |||
61 | 10800 | PermeabilityType permeabilityAtPos(const GlobalPosition& globalPos) const | |
62 | 10800 | { return permeability_ * brinkmanEpsilonAtPos(globalPos); } | |
63 | |||
64 | ✗ | PermeabilityType inversePermeability(const Element& element, const FVElementGeometry& fvGeometry, const SubControlVolume& scv) const | |
65 | 334844 | { return inversePermeability_; } | |
66 | |||
67 | ✗ | Scalar brinkmanEpsilon(const Element& element, const FVElementGeometry& fvGeometry, const SubControlVolume& scv) const | |
68 | 2692284 | { return brinkmanEpsilonAtPos(scv.dofPosition()); } | |
69 | |||
70 | 3070328 | Scalar brinkmanEpsilonAtPos(const GlobalPosition& globalPos) const | |
71 | { | ||
72 | 3070328 | const auto p = globalPos - pmCenter_; | |
73 | 3070328 | GlobalPosition q(0.0); | |
74 |
2/2✓ Branch 0 taken 3070328 times.
✓ Branch 1 taken 1535164 times.
|
9210984 | for (int i = 0; i < dimWorld; ++i) |
75 | 30703280 | q[i] += std::abs(p[i]) - 0.5*pmDimensions_[i] + pmRoundedness_; | |
76 | |||
77 |
8/8✓ Branch 0 taken 1222078 times.
✓ Branch 1 taken 313086 times.
✓ Branch 2 taken 1222078 times.
✓ Branch 3 taken 313086 times.
✓ Branch 4 taken 649296 times.
✓ Branch 5 taken 885868 times.
✓ Branch 6 taken 649296 times.
✓ Branch 7 taken 885868 times.
|
11028968 | const auto length = std::hypot(std::max(q[0], 0.0), std::max(q[1], 0.0)); |
78 |
8/8✓ Branch 0 taken 299384 times.
✓ Branch 1 taken 1235780 times.
✓ Branch 2 taken 299384 times.
✓ Branch 3 taken 1235780 times.
✓ Branch 4 taken 299384 times.
✓ Branch 5 taken 1235780 times.
✓ Branch 6 taken 993328 times.
✓ Branch 7 taken 541836 times.
|
9809752 | const auto sdf = std::min(std::max(q[0], q[1]), 0.0) + length - pmRoundedness_; |
79 |
2/2✓ Branch 0 taken 941384 times.
✓ Branch 1 taken 593780 times.
|
3070328 | const auto clampedSDF = -std::clamp(sdf, -pmTransitionHalfLength_, pmTransitionHalfLength_) + pmTransitionHalfLength_; |
80 | 3070328 | return 0.5*clampedSDF/pmTransitionHalfLength_; | |
81 | } | ||
82 | |||
83 | private: | ||
84 | 8 | void initPermeability_() | |
85 | { | ||
86 | 8 | Scalar k = getParam<Scalar>("SpatialParams.Permeability"); | |
87 | 8 | Scalar anisotropyRatio = getParam<Scalar>("SpatialParams.AnisotropyRatio", 0.0); | |
88 | 16 | permeability_[0][0] = k; | |
89 | 16 | permeability_[1][1] = k * anisotropyRatio; | |
90 | |||
91 | // rotate the tensor by theta | ||
92 | 8 | Scalar theta_ = getParam<Scalar>("SpatialParams.PermeabilityRotation", 0.0); | |
93 | |||
94 | // degrees to Radians for the rotation angle, store rotation entries | ||
95 | 8 | Scalar radTheta = theta_ * M_PI / 180.0; | |
96 | 8 | Scalar cosTheta = std::cos(radTheta); | |
97 | 8 | Scalar sinTheta = std::sin(radTheta); | |
98 | |||
99 | // create a rotation matrix according to the rotation angle | ||
100 | PermeabilityType rotationMatrix; | ||
101 | 16 | rotationMatrix[0][0] = cosTheta; | |
102 | 16 | rotationMatrix[0][1] = sinTheta * -1.0; | |
103 | 16 | rotationMatrix[1][0] = sinTheta; | |
104 | 16 | rotationMatrix[1][1] = cosTheta; | |
105 | |||
106 | // rotate the permeability tensor | ||
107 | 8 | PermeabilityType originalPermeability = permeability_; | |
108 | 24 | permeability_ = rotationMatrix * originalPermeability * getTransposed(rotationMatrix); | |
109 | |||
110 | 8 | inversePermeability_ = permeability_; | |
111 | 16 | inversePermeability_.invert(); | |
112 | 8 | } | |
113 | |||
114 | PermeabilityType permeability_; | ||
115 | PermeabilityType inversePermeability_; | ||
116 | |||
117 | PermeabilityType ffPermeability_; | ||
118 | |||
119 | GlobalPosition pmLowerLeft_; | ||
120 | GlobalPosition pmUpperRight_; | ||
121 | GlobalPosition pmCenter_; | ||
122 | GlobalPosition pmDimensions_; | ||
123 | Scalar pmRoundedness_; | ||
124 | Scalar pmTransitionHalfLength_; | ||
125 | |||
126 | static constexpr Scalar eps_ = 1e-7; | ||
127 | }; | ||
128 | |||
129 | } // end namespace Dumux | ||
130 | |||
131 | #endif | ||
132 |