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 EmbeddedTests | ||
10 | * \brief The spatial parameters class blood flow problem. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_ROOT_SPATIALPARAMS_HH | ||
14 | #define DUMUX_ROOT_SPATIALPARAMS_HH | ||
15 | |||
16 | #include <dumux/common/parameters.hh> | ||
17 | #include <dumux/io/grid/griddata.hh> | ||
18 | #include <dumux/material/components/simpleh2o.hh> | ||
19 | #include <dumux/porousmediumflow/fvspatialparams1p.hh> | ||
20 | |||
21 | namespace Dumux { | ||
22 | |||
23 | /*! | ||
24 | * \ingroup EmbeddedTests | ||
25 | * \brief Definition of the spatial parameters for the blood flow problem. | ||
26 | */ | ||
27 | template<class GridGeometry, class Scalar> | ||
28 | class RootSpatialParams | ||
29 | : public FVPorousMediumFlowSpatialParamsOneP<GridGeometry, Scalar, RootSpatialParams<GridGeometry, Scalar>> | ||
30 | { | ||
31 | using ThisType = RootSpatialParams<GridGeometry, Scalar>; | ||
32 | using ParentType = FVPorousMediumFlowSpatialParamsOneP<GridGeometry, Scalar, ThisType>; | ||
33 | using Grid = typename GridGeometry::Grid; | ||
34 | using GridView = typename GridGeometry::GridView; | ||
35 | using Element = typename GridView::template Codim<0>::Entity; | ||
36 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
37 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
38 | |||
39 | //! Indices to access the parameters in the dgf file | ||
40 | enum DGFParamIndicesSurfaceModel { | ||
41 | orderIdx = 0, | ||
42 | rootIdIdx = 1, | ||
43 | surfaceIdx = 2, | ||
44 | massIdx = 3, | ||
45 | plantIdx = 5 | ||
46 | }; | ||
47 | |||
48 | public: | ||
49 | // export permeability type | ||
50 | using PermeabilityType = Scalar; | ||
51 | |||
52 | 3 | RootSpatialParams(std::shared_ptr<const GridGeometry> gridGeometry, | |
53 | std::shared_ptr<const GridData<Grid>> gridData) | ||
54 |
3/12✓ Branch 2 taken 3 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 3 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 3 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
|
6 | : ParentType(gridGeometry), gridData_(gridData) |
55 | { | ||
56 |
1/2✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
|
3 | porosity_ = getParam<Scalar>("Root.SpatialParams.Porosity", 0.4); |
57 |
1/2✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
|
3 | constantKx_ = getParam<Scalar>("Root.SpatialParams.Kx", 5.0968e-17); |
58 |
1/4✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
3 | constantKr_ = getParam<Scalar>("Root.SpatialParams.Kr", 2.04e-13); |
59 | |||
60 | 6 | const auto& gv = gridGeometry->gridView(); | |
61 |
1/2✓ Branch 2 taken 3 times.
✗ Branch 3 not taken.
|
3 | radii_.resize(gv.size(0)); |
62 |
5/6✓ Branch 1 taken 3520 times.
✓ Branch 2 taken 3 times.
✓ Branch 3 taken 3520 times.
✓ Branch 4 taken 3 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 3520 times.
|
3523 | for (const auto& element : elements(gv)) |
63 | { | ||
64 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 3520 times.
|
3520 | const auto eIdx = gv.indexSet().index(element); |
65 | 3520 | auto level0element = element; | |
66 |
2/2✓ Branch 0 taken 88 times.
✓ Branch 1 taken 3520 times.
|
7128 | for (auto levelIdx = element.level(); levelIdx != 0; levelIdx--) |
67 | 264 | level0element = level0element.father(); | |
68 | |||
69 |
2/4✓ Branch 1 taken 3520 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 3520 times.
✗ Branch 5 not taken.
|
7040 | const auto& params = gridData_->parameters(level0element); |
70 |
4/4✓ Branch 0 taken 44 times.
✓ Branch 1 taken 3476 times.
✓ Branch 2 taken 44 times.
✓ Branch 3 taken 3476 times.
|
7040 | if (params.size() == 1) // assume only radius is given |
71 | 88 | radii_[eIdx] = params[0]; | |
72 | else // assume DGFParamIndicesSurfaceModel | ||
73 | { | ||
74 |
2/4✓ Branch 1 taken 3476 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 3476 times.
✗ Branch 5 not taken.
|
3476 | const Scalar rootLength = element.geometry().volume(); |
75 |
2/4✓ Branch 1 taken 3476 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 3476 times.
✗ Branch 5 not taken.
|
6952 | const Scalar rootSurface = gridData_->parameters(level0element)[DGFParamIndicesSurfaceModel::surfaceIdx]/(1 << element.level()); |
76 | 6952 | radii_[eIdx] = rootSurface / rootLength / 2.0 / M_PI; | |
77 | } | ||
78 | } | ||
79 | 3 | } | |
80 | |||
81 | /*! | ||
82 | * \brief Returns the intrinsic permeability for the current sub-control volume in [m^2]. | ||
83 | * | ||
84 | * \note Kx has units [m^4/(Pa*s)] so we have to divide by the cross-section area | ||
85 | * and multiply with a characteristic viscosity | ||
86 | */ | ||
87 | template<class ElementSolution> | ||
88 | ✗ | PermeabilityType permeability(const Element& element, | |
89 | const SubControlVolume& scv, | ||
90 | const ElementSolution& elemSol) const | ||
91 | { | ||
92 | ✗ | const Scalar r = radius(this->gridGeometry().elementMapper().index(element)); | |
93 | ✗ | return constantKx_ / (M_PI*r*r) * Components::SimpleH2O<Scalar>::liquidViscosity(285.15, 1e5); | |
94 | } | ||
95 | |||
96 | /*! | ||
97 | * \brief Returns the radius of the circular pipe for the current sub-control volume in [m]. | ||
98 | * | ||
99 | * \param eIdxGlobal the index of the element | ||
100 | */ | ||
101 | Scalar radius(std::size_t eIdxGlobal) const | ||
102 | { | ||
103 | 17360 | return radii_[eIdxGlobal]; | |
104 | } | ||
105 | |||
106 | /*! | ||
107 | * \brief Returns the radial permeability. | ||
108 | * | ||
109 | * \param eIdxGlobal the index of the element | ||
110 | */ | ||
111 | ✗ | Scalar Kr(std::size_t eIdxGlobal) const | |
112 | { | ||
113 | ✗ | return constantKr_; | |
114 | } | ||
115 | |||
116 | const std::vector<Scalar>& getRadii() const | ||
117 |
1/2✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
|
3 | { return radii_; } |
118 | |||
119 | /*! | ||
120 | * \brief Returns the porosity \f$[-]\f$. | ||
121 | * | ||
122 | * \param globalPos the scv center | ||
123 | */ | ||
124 | ✗ | Scalar porosityAtPos(const GlobalPosition& globalPos) const | |
125 | { | ||
126 | ✗ | return porosity_; | |
127 | } | ||
128 | |||
129 | /*! | ||
130 | * \brief Returns how much the domain is extruded at a given sub-control volume. | ||
131 | * | ||
132 | * The extrusion factor here extrudes the 1d line to a circular tube with | ||
133 | * cross-section area pi*r^2. | ||
134 | */ | ||
135 | template<class ElementSolution> | ||
136 | ✗ | Scalar extrusionFactor(const Element &element, | |
137 | const SubControlVolume &scv, | ||
138 | const ElementSolution& elemSol) const | ||
139 | { | ||
140 | ✗ | const auto eIdx = this->gridGeometry().elementMapper().index(element); | |
141 | ✗ | const auto r = radius(eIdx); | |
142 | ✗ | return M_PI*r*r; | |
143 | } | ||
144 | |||
145 | private: | ||
146 | std::shared_ptr<const GridData<Grid>> gridData_; | ||
147 | Scalar porosity_, constantKx_, constantKr_; | ||
148 | std::vector<Scalar> radii_; | ||
149 | }; | ||
150 | |||
151 | } // end namespace Dumux | ||
152 | |||
153 | #endif | ||
154 |