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 OnePTests | ||
10 | * \brief The spatial parameters class for the test problem using the | ||
11 | * 1p box model. | ||
12 | */ | ||
13 | |||
14 | #ifndef DUMUX_1P_TEST_SPATIALPARAMS_HH | ||
15 | #define DUMUX_1P_TEST_SPATIALPARAMS_HH | ||
16 | |||
17 | #include <dumux/porousmediumflow/properties.hh> | ||
18 | #include <dumux/porousmediumflow/fvspatialparams1p.hh> | ||
19 | #include <dumux/material/gstatrandomfield.hh> | ||
20 | |||
21 | namespace Dumux { | ||
22 | |||
23 | /*! | ||
24 | * \ingroup OnePTests | ||
25 | * \brief The spatial parameters class for the test problem using the | ||
26 | * 1p box model. | ||
27 | */ | ||
28 | template<class GridGeometry, class Scalar> | ||
29 | class OnePTestSpatialParams | ||
30 | : public FVPorousMediumFlowSpatialParamsOneP<GridGeometry, Scalar, | ||
31 | OnePTestSpatialParams<GridGeometry, Scalar>> | ||
32 | { | ||
33 | using GridView = typename GridGeometry::GridView; | ||
34 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
35 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
36 | using IndexSet = typename GridView::IndexSet; | ||
37 | |||
38 | using ThisType = OnePTestSpatialParams<GridGeometry, Scalar>; | ||
39 | using ParentType = FVPorousMediumFlowSpatialParamsOneP<GridGeometry, Scalar, ThisType>; | ||
40 | |||
41 | enum { | ||
42 | dim=GridView::dimension, | ||
43 | dimWorld=GridView::dimensionworld | ||
44 | }; | ||
45 | |||
46 | using Element = typename GridView::template Codim<0>::Entity; | ||
47 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
48 | |||
49 | public: | ||
50 | // export permeability type | ||
51 | using PermeabilityType = Scalar; | ||
52 | |||
53 | 5 | OnePTestSpatialParams(std::shared_ptr<const GridGeometry> gridGeometry) | |
54 | : ParentType(gridGeometry), | ||
55 |
1/2✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
|
5 | randomPermeability_(gridGeometry->gridView().size(dim), 0.0), |
56 |
12/28✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 5 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 5 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 5 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 5 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 5 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 5 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 5 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 5 times.
✗ Branch 26 not taken.
✓ Branch 27 taken 5 times.
✗ Branch 28 not taken.
✓ Branch 29 taken 5 times.
✗ Branch 30 not taken.
✓ Branch 31 taken 5 times.
✗ Branch 32 not taken.
✗ Branch 33 not taken.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
|
30 | indexSet_(gridGeometry->gridView().indexSet()) |
57 | { | ||
58 |
1/4✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
5 | randomField_ = getParam<bool>("SpatialParams.RandomField", false); |
59 |
1/2✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
|
5 | permeability_ = getParam<Scalar>("SpatialParams.Permeability"); |
60 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
5 | if(!randomField_) |
61 |
1/2✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
|
5 | permeabilityLens_ = getParam<Scalar>("SpatialParams.PermeabilityLens"); |
62 | else | ||
63 | ✗ | initRandomField(*gridGeometry); | |
64 | |||
65 |
1/2✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
|
5 | lensLowerLeft_ = getParam<GlobalPosition>("SpatialParams.LensLowerLeft"); |
66 |
1/2✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
|
5 | lensUpperRight_ = getParam<GlobalPosition>("SpatialParams.LensUpperRight"); |
67 | 5 | } | |
68 | |||
69 | /*! | ||
70 | * \brief Function for defining the (intrinsic) permeability \f$[m^2]\f$. | ||
71 | * | ||
72 | * \param element The element | ||
73 | * \param scv The sub-control volume | ||
74 | * \param elemSol The element solution vector | ||
75 | * \return The intrinsic permeability | ||
76 | */ | ||
77 | template<class ElementSolution> | ||
78 | 12248 | PermeabilityType permeability(const Element& element, | |
79 | const SubControlVolume& scv, | ||
80 | const ElementSolution& elemSol) const | ||
81 | { | ||
82 |
2/2✓ Branch 0 taken 4152 times.
✓ Branch 1 taken 8096 times.
|
36744 | if (isInLens_(scv.dofPosition())) |
83 | { | ||
84 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4152 times.
|
4152 | if(randomField_) |
85 | ✗ | return randomPermeability_[indexSet_.index(element)]; | |
86 | else | ||
87 | 4152 | return permeabilityLens_; | |
88 | } | ||
89 | else | ||
90 | 8096 | return permeability_; | |
91 | } | ||
92 | |||
93 | /*! \brief Defines the porosity in [-]. | ||
94 | * | ||
95 | * \param globalPos The global position where we evaluate | ||
96 | */ | ||
97 | ✗ | Scalar porosityAtPos(const GlobalPosition& globalPos) const | |
98 | ✗ | { return 0.4; } | |
99 | |||
100 | /*! | ||
101 | * \brief This method allows the generation of a statistical field using gstat | ||
102 | * | ||
103 | * \param gg The finite-volume grid geometry used by the problem | ||
104 | */ | ||
105 | ✗ | void initRandomField(const GridGeometry& gg) | |
106 | { | ||
107 | ✗ | const auto& gridView = gg.gridView(); | |
108 | ✗ | const auto& elementMapper = gg.elementMapper(); | |
109 | ✗ | const auto gStatControlFile = getParam<std::string>("Gstat.ControlFile"); | |
110 | ✗ | const auto gStatInputFile = getParam<std::string>("Gstat.InputFile"); | |
111 | ✗ | const auto outputFilePrefix = getParam<std::string>("Gstat.OutputFilePrefix"); | |
112 | |||
113 | // create random permeability object | ||
114 | using RandomField = GstatRandomField<GridView, Scalar>; | ||
115 | ✗ | RandomField randomPermeabilityField(gridView, elementMapper); | |
116 | ✗ | randomPermeabilityField.create(gStatControlFile, | |
117 | gStatInputFile, | ||
118 | outputFilePrefix + ".dat", | ||
119 | RandomField::FieldType::log10, | ||
120 | true); | ||
121 | ✗ | randomPermeability_.resize(gridView.size(dim), 0.0); | |
122 | |||
123 | // copy vector from the temporary gstat object | ||
124 | ✗ | randomPermeability_ = randomPermeabilityField.data(); | |
125 | ✗ | } | |
126 | |||
127 | //! get the permeability field for output | ||
128 | const std::vector<Scalar>& getPermField() const | ||
129 | ✗ | { return randomPermeability_; } | |
130 | |||
131 | private: | ||
132 | bool isInLens_(const GlobalPosition &globalPos) const | ||
133 | { | ||
134 |
2/2✓ Branch 0 taken 19716 times.
✓ Branch 1 taken 4152 times.
|
23868 | for (int i = 0; i < dimWorld; ++i) { |
135 |
12/12✓ Branch 0 taken 15668 times.
✓ Branch 1 taken 4048 times.
✓ Branch 2 taken 15668 times.
✓ Branch 3 taken 4048 times.
✓ Branch 4 taken 15668 times.
✓ Branch 5 taken 4048 times.
✓ Branch 6 taken 11620 times.
✓ Branch 7 taken 4048 times.
✓ Branch 8 taken 11620 times.
✓ Branch 9 taken 4048 times.
✓ Branch 10 taken 11620 times.
✓ Branch 11 taken 4048 times.
|
59148 | if (globalPos[i] < lensLowerLeft_[i] + eps_ || globalPos[i] > lensUpperRight_[i] - eps_) |
136 | return false; | ||
137 | } | ||
138 | return true; | ||
139 | } | ||
140 | |||
141 | bool randomField_; | ||
142 | GlobalPosition lensLowerLeft_; | ||
143 | GlobalPosition lensUpperRight_; | ||
144 | |||
145 | Scalar permeability_, permeabilityLens_; | ||
146 | std::vector<Scalar> randomPermeability_; | ||
147 | |||
148 | const IndexSet& indexSet_; | ||
149 | static constexpr Scalar eps_ = 1.5e-7; | ||
150 | }; | ||
151 | |||
152 | } // end namespace | ||
153 | |||
154 | #endif | ||
155 |