version 3.10-dev
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-FileCopyrightInfo: 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>
19
20namespace Dumux {
21
26template<class FVElementGeometry, class PV>
28{
29 using GridGeometry = typename FVElementGeometry::GridGeometry;
30 using GridView = typename GridGeometry::GridView;
31 using Element = typename GridView::template Codim<0>::Entity;
32
33 static constexpr int dim = GridView::dimension;
34
35 // Dofs are located at corners and in the element center
36 static constexpr int numCubeDofs = Detail::LocalDofTraits<
37 GridView, typename GridGeometry::DiscretizationMethod
38 >::numCubeElementDofs;
39
40public:
42 using PrimaryVariables = PV;
43
46
48 template<class SolutionVector>
49 CVFEElementSolution(const Element& element, const SolutionVector& sol,
50 const GridGeometry& gridGeometry)
51 {
52 update(element, sol, gridGeometry);
53 }
54
56 template<class ElementVolumeVariables>
57 CVFEElementSolution(const Element& element, const ElementVolumeVariables& elemVolVars,
58 const FVElementGeometry& fvGeometry)
59 {
60 priVars_.resize(fvGeometry.numScv());
61 for (const auto& scv : scvs(fvGeometry))
62 priVars_[scv.localDofIndex()] = elemVolVars[scv].priVars();
63 }
64
66 template<class SolutionVector>
67 void update(const Element& element, const SolutionVector& sol,
68 const GridGeometry& gridGeometry)
69 {
70 // TODO: this implementation only works if there is only one dof per codim/entity
71 // As local-global mappings are provided by the grid geometry
72 // this logic should maybe become part of the grid geometry.
73 const auto& localCoeff = gridGeometry.feCache().get(element.type()).localCoefficients();
74 priVars_.resize(localCoeff.size());
75 for (int localDofIdx = 0; localDofIdx < localCoeff.size(); ++localDofIdx)
76 {
77 const auto& localKey = localCoeff.localKey(localDofIdx);
78 priVars_[localDofIdx] = sol[
79 gridGeometry.dofMapper().subIndex(
80 element, localKey.subEntity(), localKey.codim()
81 )
82 ];
83 }
84 }
85
87 template<class SolutionVector>
88 void update(const Element& element, const SolutionVector& sol,
89 const FVElementGeometry& fvGeometry)
90 {
91 priVars_.resize(fvGeometry.numScv());
92 for (const auto& scv : scvs(fvGeometry))
93 priVars_[scv.localDofIndex()] = sol[scv.dofIndex()];
94 }
95
97 std::size_t size() const
98 { return priVars_.size(); }
99
101 template<typename IndexType>
102 const PrimaryVariables& operator [](IndexType i) const
103 { return priVars_[i]; }
104
106 template<typename IndexType>
108 { return priVars_[i]; }
109
110private:
111 Dune::ReservedVector<PrimaryVariables, numCubeDofs> priVars_;
112};
113
115template<class Element, class SolutionVector, class GridGeometry>
116auto elementSolution(const Element& element, const SolutionVector& sol, const GridGeometry& gg)
117-> std::enable_if_t<DiscretizationMethods::isCVFE<typename GridGeometry::DiscretizationMethod>,
118 CVFEElementSolution<typename GridGeometry::LocalView,
119 std::decay_t<decltype(std::declval<SolutionVector>()[0])>>
120 >
121{
122 using PrimaryVariables = std::decay_t<decltype(std::declval<SolutionVector>()[0])>;
123 return CVFEElementSolution<typename GridGeometry::LocalView, PrimaryVariables>(element, sol, gg);
124}
125
127template<class Element, class ElementVolumeVariables, class FVElementGeometry>
128auto elementSolution(const Element& element, const ElementVolumeVariables& elemVolVars, const FVElementGeometry& gg)
129-> std::enable_if_t<DiscretizationMethods::isCVFE<typename FVElementGeometry::GridGeometry::DiscretizationMethod>,
130 CVFEElementSolution<FVElementGeometry,
131 typename ElementVolumeVariables::VolumeVariables::PrimaryVariables>>
132{
133 using PrimaryVariables = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables;
134 return CVFEElementSolution<FVElementGeometry, PrimaryVariables>(element, elemVolVars, gg);
135}
136
137} // end namespace Dumux
138
139#endif
The element solution vector.
Definition: cvfe/elementsolution.hh:28
const PrimaryVariables & operator[](IndexType i) const
bracket operator const access
Definition: cvfe/elementsolution.hh:102
PV PrimaryVariables
export the primary variables type
Definition: cvfe/elementsolution.hh:42
CVFEElementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
Constructor with element and solution and grid geometry.
Definition: cvfe/elementsolution.hh:49
CVFEElementSolution(const Element &element, const ElementVolumeVariables &elemVolVars, const FVElementGeometry &fvGeometry)
Constructor with element and elemVolVars and fvGeometry.
Definition: cvfe/elementsolution.hh:57
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:88
std::size_t size() const
return the size of the element solution
Definition: cvfe/elementsolution.hh:97
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:67
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
Element-specific traits of grid geometries / discretization schemes.
The available discretization methods in Dumux.
Definition: adapt.hh:17
Definition: localdoftraits.hh:18