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 properties for the incompressible test. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_ONEP_SUB_TEST_PROBLEM_HH | ||
14 | #define DUMUX_ONEP_SUB_TEST_PROBLEM_HH | ||
15 | |||
16 | #include <dune/common/indices.hh> | ||
17 | |||
18 | #include <dumux/common/boundarytypes.hh> | ||
19 | #include <dumux/common/properties.hh> | ||
20 | #include <dumux/common/parameters.hh> | ||
21 | #include <dumux/common/numeqvector.hh> | ||
22 | |||
23 | #include <dumux/porousmediumflow/problem.hh> | ||
24 | |||
25 | namespace Dumux { | ||
26 | |||
27 | /*! | ||
28 | * \ingroup BoundaryTests | ||
29 | * \brief Multidomain test problem for the incompressible one-phase model. | ||
30 | * | ||
31 | * The circular model domain consists of two subdomains: | ||
32 | * an inner circle and an outer ring. | ||
33 | * Methane is injected in the center and spreads over the | ||
34 | * coupling boundary into the outer domain. | ||
35 | */ | ||
36 | template<class TypeTag, std::size_t tag> | ||
37 | class OnePTestProblem | ||
38 | : public PorousMediumFlowProblem<TypeTag> | ||
39 | { | ||
40 | using ParentType = PorousMediumFlowProblem<TypeTag>; | ||
41 | using CouplingManager = GetPropType<TypeTag, Properties::CouplingManager>; | ||
42 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
43 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
44 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
45 | using GridView = typename GridGeometry::GridView; | ||
46 | using Element = typename GridView::template Codim<0>::Entity; | ||
47 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
48 | using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>; | ||
49 | using NumEqVector = Dumux::NumEqVector<PrimaryVariables>; | ||
50 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
51 | using BoundaryTypes = Dumux::BoundaryTypes<GetPropType<TypeTag, Properties::ModelTraits>::numEq()>; | ||
52 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
53 | static constexpr int dimWorld = GridView::dimensionworld; | ||
54 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
55 | static constexpr auto domainIdx = Dune::index_constant<tag>{}; | ||
56 | |||
57 | public: | ||
58 | 4 | OnePTestProblem(std::shared_ptr<const GridGeometry> gridGeometry, | |
59 | std::shared_ptr<CouplingManager> couplingManager, | ||
60 | const std::string& paramGroup = "") | ||
61 | : ParentType(gridGeometry, paramGroup) | ||
62 |
3/10✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
|
8 | , couplingManager_(couplingManager) |
63 | { | ||
64 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
4 | injectionRate_ = getParam<double>("Problem.InjectionRate"); |
65 |
8/24✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
|
4 | problemName_ = getParam<std::string>("Vtk.OutputName") + "_" + getParamFromGroup<std::string>(paramGroup, "Problem.Name"); |
66 | 4 | } | |
67 | |||
68 | /*! | ||
69 | * \brief The problem name. | ||
70 | */ | ||
71 | const std::string& name() const | ||
72 | { | ||
73 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
2 | return problemName_; |
74 | } | ||
75 | |||
76 | /*! | ||
77 | * \brief Specifies which kind of boundary condition should be | ||
78 | * used for which equation on a given boundary segment. | ||
79 | * | ||
80 | * \param element The finite element | ||
81 | * \param scvf The sub control volume face | ||
82 | */ | ||
83 | 209920 | BoundaryTypes boundaryTypes(const Element &element, | |
84 | const SubControlVolumeFace &scvf) const | ||
85 | { | ||
86 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
209920 | BoundaryTypes values; |
87 | 209920 | values.setAllDirichlet(); | |
88 | |||
89 |
1/7✗ Branch 1 not taken.
✓ Branch 2 taken 104960 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
419840 | if (couplingManager_->isCoupled(domainIdx, scvf)) |
90 | values.setAllCouplingNeumann(); | ||
91 | |||
92 | 209920 | return values; | |
93 | } | ||
94 | |||
95 | /*! | ||
96 | * \brief Evaluates the boundary conditions for a Neumann boundary segment. | ||
97 | * | ||
98 | * This is the method for the case where the Neumann condition is | ||
99 | * potentially solution dependent. | ||
100 | * | ||
101 | * \param element The finite element | ||
102 | * \param fvGeometry The finite-volume geometry | ||
103 | * \param elemVolVars All volume variables for the element | ||
104 | * \param elemFluxVarsCache Flux variables caches for all faces in stencil | ||
105 | * \param scvf The sub control volume face | ||
106 | * | ||
107 | * Negative values mean influx. | ||
108 | * E.g. for the mass balance that would the mass flux in \f$ [ kg / (m^2 \cdot s)] \f$. | ||
109 | */ | ||
110 | template<class ElementVolumeVariables, class ElementFluxVarsCache> | ||
111 | 199424 | NumEqVector neumann(const Element& element, | |
112 | const FVElementGeometry& fvGeometry, | ||
113 | const ElementVolumeVariables& elemVolVars, | ||
114 | const ElementFluxVarsCache& elemFluxVarsCache, | ||
115 | const SubControlVolumeFace& scvf) const | ||
116 | { | ||
117 | 199424 | NumEqVector values(0.0); | |
118 | 199424 | const auto bcTypes = boundaryTypes(element, scvf); | |
119 |
2/4✓ Branch 0 taken 99712 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 99712 times.
✗ Branch 3 not taken.
|
398848 | if (bcTypes.hasCouplingNeumann()) |
120 | 398848 | values[Indices::conti0EqIdx] = couplingManager_->advectiveFluxCoupling(domainIdx, element, fvGeometry, elemVolVars, scvf, FluidSystem::phase0Idx); | |
121 | |||
122 | 199424 | return values; | |
123 | } | ||
124 | |||
125 | /*! | ||
126 | * \brief Evaluates the boundary conditions for a Dirichlet control volume. | ||
127 | * | ||
128 | * \param globalPos The center of the finite volume which ought to be set. | ||
129 | * | ||
130 | * For this method, the \a values parameter stores primary variables. | ||
131 | */ | ||
132 | ✗ | PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const | |
133 | { | ||
134 | 16400 | PrimaryVariables values(0.0); | |
135 | 18376 | values[Indices::pressureIdx] = 1.0e6; | |
136 | if (NumEqVector::dimension > 1) | ||
137 | 484 | values[Indices::pressureIdx + 1] = 1e-10; // fully saturated | |
138 | ✗ | return values; | |
139 | } | ||
140 | |||
141 | /*! | ||
142 | * \brief Applies a vector of point sources which are possibly solution dependent. | ||
143 | * | ||
144 | * \param pointSources A vector of PointSource s that contain | ||
145 | source values for all phases and space positions. | ||
146 | * | ||
147 | * For this method, the values method of the point source | ||
148 | * has to return the absolute rate values in units | ||
149 | * \f$ [ \textnormal{unit of conserved quantity} / s ] \f$. | ||
150 | * Positive values mean that the conserved quantity is created, negative ones mean that it vanishes. | ||
151 | * E.g. for the mass balance that would be a mass rate in \f$ [ kg / s ] \f$. | ||
152 | */ | ||
153 | template<class PointSource> | ||
154 | ✗ | void addPointSources(std::vector<PointSource>& pointSources) const | |
155 | { | ||
156 | ✗ | NumEqVector values(0.0); | |
157 | // if this is the 2p problem inject methane | ||
158 | if (NumEqVector::dimension > 1) | ||
159 | ✗ | values[Indices::conti0EqIdx + 1] = injectionRate_; // kg/s | |
160 | ✗ | pointSources.emplace_back(GlobalPosition(0.0), values); | |
161 | ✗ | } | |
162 | |||
163 | /*! | ||
164 | * \brief Evaluates the initial value for a control volume. | ||
165 | * | ||
166 | * \param globalPos The global position | ||
167 | */ | ||
168 | ✗ | PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const | |
169 | 2460 | { return dirichletAtPos(globalPos); } | |
170 | |||
171 | |||
172 | private: | ||
173 | Scalar injectionRate_; | ||
174 | std::shared_ptr<CouplingManager> couplingManager_; | ||
175 | static constexpr Scalar eps_ = 1e-7; | ||
176 | std::string problemName_; | ||
177 | }; | ||
178 | |||
179 | } // end namespace Dumux | ||
180 | |||
181 | #endif | ||
182 |