version 3.11-dev
Loading...
Searching...
No Matches
cvfe/elementsolution.hh
Go to the documentation of this file.
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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_CVFE_ELEMENT_SOLUTION_HH
13#define DUMUX_CVFE_ELEMENT_SOLUTION_HH
14
15#include <type_traits>
16#include <dune/common/reservedvector.hh>
17
18#include <dumux/common/typetraits/localdofs_.hh>
21
22namespace Dumux {
23
24namespace Detail {
25
26template<class ElementVariables>
27constexpr auto primaryVariablesType()
28{
29 if constexpr (requires { typename ElementVariables::Variables::PrimaryVariables; })
30 return std::type_identity<typename ElementVariables::Variables::PrimaryVariables>{};
31 else if constexpr (requires { typename ElementVariables::VolumeVariables::PrimaryVariables; })
32 return std::type_identity<typename ElementVariables::VolumeVariables::PrimaryVariables>{};
33 else
34 return std::type_identity<void>{};
35}
36
37template<class ElementVariables>
39
40} // end namespace Detail
41
46template<class FVElementGeometry, class PV>
48{
49 using GridGeometry = typename FVElementGeometry::GridGeometry;
50 using GridView = typename GridGeometry::GridView;
51 using Element = typename GridView::template Codim<0>::Entity;
52
53 static constexpr int dim = GridView::dimension;
54
55 // Maximum number of local element dofs
56 static constexpr int maxNumLocalDofs = Detail::LocalDofs::maxNumLocalDofs<FVElementGeometry>();
57
58public:
60 using PrimaryVariables = PV;
61
64
66 template<class SolutionVector>
67 CVFEElementSolution(const Element& element, const SolutionVector& sol,
68 const GridGeometry& gridGeometry)
69 {
70 update(element, sol, gridGeometry);
71 }
72
74 template<class ElementVariables>
75 CVFEElementSolution(const Element& element, const ElementVariables& elemVars,
76 const FVElementGeometry& fvGeometry)
77 {
78 priVars_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
79 for (const auto& localDof : localDofs(fvGeometry))
80 priVars_[localDof.index()] = elemVars[localDof.index()].priVars();
81 }
82
84 template<class SolutionVector>
85 void update(const Element& element, const SolutionVector& sol,
86 const GridGeometry& gridGeometry)
87 {
88 // For schemes with multiple DOFs per entity (e.g. PQ3 with 2 DOFs per edge)
89 // the formula subIndex + localKey.index() may need an orientation flip for
90 // simplex elements (see HybridPQ3GeometryHelper::dofIndex for details).
91 // If the grid geometry exposes a flip-aware dofIndex(element, localKey)
92 // method, use it; otherwise fall back to the generic formula.
93 // As local-global mappings are provided by the grid geometry
94 // this logic should maybe become part of the grid geometry.
95 const auto& localCoeff = gridGeometry.feCache().get(element.type()).localCoefficients();
96 priVars_.resize(localCoeff.size());
97 for (int localDofIdx = 0; localDofIdx < localCoeff.size(); ++localDofIdx)
98 {
99 const auto& localKey = localCoeff.localKey(localDofIdx);
100 if constexpr (requires { gridGeometry.dofIndex(element, localKey); })
101 priVars_[localDofIdx] = sol[gridGeometry.dofIndex(element, localKey)];
102 else
103 priVars_[localDofIdx] = sol[
104 gridGeometry.dofMapper().subIndex(
105 element, localKey.subEntity(), localKey.codim()
106 ) + localKey.index()
107 ];
108 }
109 }
110
112 template<class SolutionVector>
113 void update(const Element& element, const SolutionVector& sol,
114 const FVElementGeometry& fvGeometry)
115 {
116 priVars_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
117 for (const auto& localDof : localDofs(fvGeometry))
118 priVars_[localDof.index()] = sol[localDof.dofIndex()];
119 }
120
122 std::size_t size() const
123 { return priVars_.size(); }
124
126 template<typename IndexType>
127 const PrimaryVariables& operator [](IndexType i) const
128 { return priVars_[i]; }
129
131 template<typename IndexType>
133 { return priVars_[i]; }
134
135private:
136 Dune::ReservedVector<PrimaryVariables, maxNumLocalDofs> priVars_;
137};
138
140template<class Element, class SolutionVector, class GridGeometry>
141auto elementSolution(const Element& element, const SolutionVector& sol, const GridGeometry& gg)
142-> std::enable_if_t<DiscretizationMethods::isCVFE<typename GridGeometry::DiscretizationMethod>,
143 CVFEElementSolution<typename GridGeometry::LocalView,
144 std::decay_t<decltype(std::declval<SolutionVector>()[0])>>
145 >
146{
147 using PrimaryVariables = std::decay_t<decltype(std::declval<SolutionVector>()[0])>;
149}
150
152template<class Element, class ElementVariables, class FVElementGeometry>
153auto elementSolution(const Element& element, const ElementVariables& elemVars, const FVElementGeometry& gg)
154-> std::enable_if_t<DiscretizationMethods::isCVFE<typename FVElementGeometry::GridGeometry::DiscretizationMethod>
155 && !std::is_void_v<Detail::PrimaryVariables_t<ElementVariables>>,
156 CVFEElementSolution<FVElementGeometry,
158 >
159{
160 using PrimaryVariables = Detail::PrimaryVariables_t<ElementVariables>;
162}
163
164} // end namespace Dumux
165
166#endif
The element solution vector.
Definition cvfe/elementsolution.hh:48
CVFEElementSolution(const Element &element, const ElementVariables &elemVars, const FVElementGeometry &fvGeometry)
Constructor with element and elemVolVars and fvGeometry.
Definition cvfe/elementsolution.hh:75
const PrimaryVariables & operator[](IndexType i) const
bracket operator const access
Definition cvfe/elementsolution.hh:127
PV PrimaryVariables
export the primary variables type
Definition cvfe/elementsolution.hh:60
CVFEElementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
Constructor with element and solution and grid geometry.
Definition cvfe/elementsolution.hh:67
void update(const Element &element, const SolutionVector &sol, const FVElementGeometry &fvGeometry)
extract the element solution from the solution vector using a local fv geometry
Definition cvfe/elementsolution.hh:113
std::size_t size() const
return the size of the element solution
Definition cvfe/elementsolution.hh:122
void update(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
extract the element solution from the solution vector using a mapper
Definition cvfe/elementsolution.hh:85
CVFEElementSolution()=default
Default constructor.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition cellcentered/elementsolution.hh:101
Class representing dofs on elements for control-volume finite element schemes.
The available discretization methods in Dumux.
Definition cvfelocalresidual.hh:25
constexpr auto primaryVariablesType()
Definition cvfe/elementsolution.hh:27
typename decltype(primaryVariablesType< ElementVariables >())::type PrimaryVariables_t
Definition cvfe/elementsolution.hh:38
Definition adapt.hh:17
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition localdof.hh:50