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 NonEquilibriumModel | ||
10 | * \brief Class storing scv and scvf variables. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_NONEQUILIBRIUM_GRID_VARIABLES_HH | ||
14 | #define DUMUX_NONEQUILIBRIUM_GRID_VARIABLES_HH | ||
15 | |||
16 | #include <memory> | ||
17 | #include <dune/common/fvector.hh> | ||
18 | #include <dune/grid/common/partitionset.hh> | ||
19 | |||
20 | #include <dumux/common/properties.hh> | ||
21 | #include <dumux/discretization/method.hh> | ||
22 | #include <dumux/discretization/fvgridvariables.hh> | ||
23 | #include <dumux/porousmediumflow/velocity.hh> | ||
24 | #include <dumux/porousmediumflow/fluxvariables.hh> | ||
25 | |||
26 | namespace Dumux { | ||
27 | |||
28 | /*! | ||
29 | * \ingroup NonEquilibriumModel | ||
30 | * \brief This class stores the velocities which are used to compute Reynolds | ||
31 | * numbers for the source terms of non-equilibrium models. | ||
32 | */ | ||
33 | template<class TypeTag> | ||
34 | class NonEquilibriumGridVariables | ||
35 | : public FVGridVariables<GetPropType<TypeTag, Properties::GridGeometry>, | ||
36 | GetPropType<TypeTag, Properties::GridVolumeVariables>, | ||
37 | GetPropType<TypeTag, Properties::GridFluxVariablesCache>> | ||
38 | { | ||
39 | using ThisType = NonEquilibriumGridVariables<TypeTag>; | ||
40 | using ParentType = FVGridVariables<GetPropType<TypeTag, Properties::GridGeometry>, | ||
41 | GetPropType<TypeTag, Properties::GridVolumeVariables>, | ||
42 | GetPropType<TypeTag, Properties::GridFluxVariablesCache>>; | ||
43 | |||
44 | using VelocityBackend = PorousMediumFlowVelocity<ThisType, PorousMediumFluxVariables<TypeTag>>; | ||
45 | |||
46 | static constexpr auto dim = ParentType::GridGeometry::GridView::dimension; // Grid and world dimension | ||
47 | static constexpr auto dimWorld = ParentType::GridGeometry::GridView::dimensionworld; | ||
48 | static constexpr int numPhases = ParentType::VolumeVariables::numFluidPhases(); | ||
49 | static constexpr bool isBox = ParentType::GridGeometry::discMethod == DiscretizationMethods::box; | ||
50 | |||
51 | public: | ||
52 | //! Export the type used for scalar values | ||
53 | using typename ParentType::Scalar; | ||
54 | using typename ParentType::GridGeometry; | ||
55 | |||
56 | //! Constructor | ||
57 | template<class Problem> | ||
58 | 5 | NonEquilibriumGridVariables(std::shared_ptr<Problem> problem, | |
59 | std::shared_ptr<const GridGeometry> gridGeometry) | ||
60 |
2/8✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 5 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
|
20 | : ParentType(problem, gridGeometry) |
61 | { | ||
62 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 5 times.
|
13 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
63 |
5/12✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 8 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 8 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 8 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 8 times.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
|
37 | velocityNorm_[phaseIdx].assign(gridGeometry->numDofs(), 0.0); |
64 | |||
65 |
3/6✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 5 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 5 times.
|
5 | velocityBackend_ = std::make_unique<VelocityBackend>(*this); |
66 | 5 | } | |
67 | |||
68 | template<class SolutionVector> | ||
69 | 445 | void calcVelocityAverage(const SolutionVector& curSol) | |
70 | { | ||
71 | using Scalar = typename SolutionVector::field_type; | ||
72 | using VelocityVector = typename Dune::FieldVector<Scalar, dimWorld>; | ||
73 | |||
74 | 445 | std::array<std::vector<VelocityVector>, numPhases> velocity; | |
75 | |||
76 |
2/2✓ Branch 0 taken 552 times.
✓ Branch 1 taken 445 times.
|
997 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
77 | { | ||
78 | if(isBox && dim == 1) | ||
79 | ✗ | velocity[phaseIdx].resize(this->gridGeometry_->gridView().size(0)); | |
80 | else | ||
81 |
4/8✓ Branch 1 taken 552 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 552 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 552 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 552 times.
✗ Branch 11 not taken.
|
1974 | velocity[phaseIdx].resize(this->gridGeometry_->numDofs()); |
82 | } | ||
83 | |||
84 |
2/4✓ Branch 1 taken 445 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 445 times.
✗ Branch 5 not taken.
|
1090 | auto fvGeometry = localView(*this->gridGeometry_); |
85 |
2/7✓ Branch 1 taken 445 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 445 times.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
1335 | auto elemVolVars = localView(this->curGridVolVars()); |
86 |
3/7✓ Branch 1 taken 445 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 445 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 245 times.
✗ Branch 7 not taken.
|
1135 | auto elemFluxVarsCache = localView(this->gridFluxVarsCache()); |
87 | |||
88 |
4/10✓ Branch 1 taken 445 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 445 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 445 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 97805 times.
✗ Branch 10 not taken.
|
196500 | for (const auto& element : elements(this->gridGeometry_->gridView(), Dune::Partitions::interior)) |
89 | { | ||
90 | 292080 | const auto eIdxGlobal = this->gridGeometry_->elementMapper().index(element); | |
91 |
1/2✓ Branch 1 taken 97360 times.
✗ Branch 2 not taken.
|
97360 | fvGeometry.bind(element); |
92 |
1/2✓ Branch 1 taken 97360 times.
✗ Branch 2 not taken.
|
97360 | elemVolVars.bind(element, fvGeometry, curSol); |
93 |
1/2✓ Branch 1 taken 97360 times.
✗ Branch 2 not taken.
|
97360 | elemFluxVarsCache.bind(element, fvGeometry, elemVolVars); |
94 | |||
95 |
2/2✓ Branch 0 taken 127120 times.
✓ Branch 1 taken 97360 times.
|
224480 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
96 | { | ||
97 |
3/6✓ Branch 1 taken 127120 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 127120 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 127120 times.
✗ Branch 8 not taken.
|
381360 | velocityBackend_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx); |
98 | |||
99 |
4/4✓ Branch 0 taken 264880 times.
✓ Branch 1 taken 230320 times.
✓ Branch 2 taken 218080 times.
✓ Branch 3 taken 183520 times.
|
896800 | for (auto&& scv : scvs(fvGeometry)) |
100 | { | ||
101 | 368080 | const auto dofIdxGlobal = scv.dofIndex(); | |
102 | if (isBox && dim == 1) | ||
103 | ✗ | velocityNorm_[phaseIdx][dofIdxGlobal] = velocity[phaseIdx][eIdxGlobal].two_norm(); | |
104 | else | ||
105 | 2576560 | velocityNorm_[phaseIdx][dofIdxGlobal] = velocity[phaseIdx][dofIdxGlobal].two_norm(); | |
106 | } | ||
107 | } //end phases | ||
108 | } //end elements | ||
109 | 445 | } // end calcVelocity | |
110 | |||
111 | /*! | ||
112 | * \brief Access to the averaged (magnitude of) velocity for each vertex. | ||
113 | * | ||
114 | * \param phaseIdx The index of the fluid phase | ||
115 | * \param dofIdxGlobal The global index of the degree of freedom | ||
116 | */ | ||
117 | const Scalar volumeDarcyMagVelocity(const unsigned int phaseIdx, | ||
118 | const unsigned int dofIdxGlobal) const | ||
119 |
0/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
11258517 | { return velocityNorm_[phaseIdx][dofIdxGlobal]; } |
120 | |||
121 | private: | ||
122 | std::array<std::vector<Dune::FieldVector<Scalar, 1> > , numPhases> velocityNorm_; | ||
123 | std::unique_ptr<VelocityBackend> velocityBackend_; | ||
124 | }; | ||
125 | |||
126 | } // end namespace Dumux | ||
127 | |||
128 | #endif | ||
129 |