version 3.8
facecentered/staggered/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_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_SOLUTION_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_SOLUTION_HH
14
15#include <type_traits>
16#include <array>
17#include <utility>
18
21
22namespace Dumux {
23
28template<class FVElementGeometry, class PV>
30{
31 using GridGeometry = typename FVElementGeometry::GridGeometry;
32 using GridView = typename GridGeometry::GridView;
33 using Element = typename GridView::template Codim<0>::Entity;
34 using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
35
36 static constexpr auto numScvsPerElement = GridView::Grid::dimension * 2;
37
38public:
40 using PrimaryVariables = PV;
41
43
45 template<class SolutionVector>
47 const SolutionVector& sol,
48 const GridGeometry& gridGeometry)
49 {
50 const auto fvGeometry = localView(gridGeometry).bindElement(element);
51
52 for (const auto& scv : scvs(fvGeometry))
53 priVars_[scv.indexInElement()] = sol[scv.dofIndex()];
54 }
55
57 template<class ElementVolumeVariables>
59 const ElementVolumeVariables& elemVolVars,
60 const FVElementGeometry& fvGeometry)
61 {
62 for (const auto& scv : scvs(fvGeometry))
63 priVars_ = elemVolVars[scv].priVars();
64 }
65
68 {
69 priVars_[0] = std::move(priVars);
70 for (int i = 1; i < priVars_.size(); ++i)
71 priVars_[i] = priVars_[0];
72 }
73
76 {
77 priVars_[0] = priVars;
78 for (int i = 1; i < priVars_.size(); ++i)
79 priVars_[i] = priVars_[0];
80 }
81
83 template<class SolutionVector>
84 void update(const Element& element, const SolutionVector& sol,
85 const GridGeometry& gridGeometry)
86 {
87 const auto fvGeometry = localView(gridGeometry).bindElement(element);
88
89 for (const auto& scv : scvs(fvGeometry))
90 priVars_[scv.indexInElement()] = sol[scv.dofIndex()];
91 }
92
94 const PrimaryVariables& operator [](SmallLocalIndexType localScvIdx) const
95 {
96 return priVars_[localScvIdx];
97 }
98
100 PrimaryVariables& operator [](SmallLocalIndexType localScvIdx)
101 {
102 return priVars_[localScvIdx];
103 }
104
106 static constexpr std::size_t size()
107 { return numScvsPerElement; }
108
109private:
110 std::array<PrimaryVariables, numScvsPerElement> priVars_;
111};
112
117template<class Element, class SolutionVector, class GridGeometry>
118auto elementSolution(const Element& element, const SolutionVector& sol, const GridGeometry& gg)
119-> std::enable_if_t<
120 GridGeometry::discMethod == DiscretizationMethods::fcstaggered,
121 FaceCenteredStaggeredElementSolution<
122 typename GridGeometry::LocalView,
123 std::decay_t<decltype(std::declval<SolutionVector>()[0])>
124 >
125>
126{ return { element, sol, gg }; }
127
132template<class Element, class ElementVolumeVariables, class FVElementGeometry>
133auto elementSolution(const Element& element, const ElementVolumeVariables& elemVolVars, const FVElementGeometry& gg)
134-> std::enable_if_t<
135 FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::fcstaggered,
136 FaceCenteredStaggeredElementSolution<
137 FVElementGeometry,
138 typename ElementVolumeVariables::VolumeVariables::PrimaryVariables
139 >
140>
141{ return { element, elemVolVars, gg }; }
142
148template<class FVElementGeometry, class PrimaryVariables>
149auto elementSolution(PrimaryVariables&& priVars)
150-> std::enable_if_t<
151 FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::fcstaggered,
152 FaceCenteredStaggeredElementSolution<
153 FVElementGeometry,
154 std::decay_t<PrimaryVariables>
155 >
156>
157{ return { std::forward<PrimaryVariables>(priVars) }; }
158
159} // end namespace Dumux
160
161#endif
The global face variables class for staggered models.
Definition: facecentered/staggered/elementsolution.hh:30
FaceCenteredStaggeredElementSolution(PrimaryVariables &&priVars)
Constructor with a primary variable object.
Definition: facecentered/staggered/elementsolution.hh:67
FaceCenteredStaggeredElementSolution(const Element &element, const ElementVolumeVariables &elemVolVars, const FVElementGeometry &fvGeometry)
Constructor with element, element volume variables and fv element geometry.
Definition: facecentered/staggered/elementsolution.hh:58
FaceCenteredStaggeredElementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
Constructor with element, solution vector and grid geometry.
Definition: facecentered/staggered/elementsolution.hh:46
PV PrimaryVariables
export the primary variables type
Definition: facecentered/staggered/elementsolution.hh:40
const PrimaryVariables & operator[](SmallLocalIndexType localScvIdx) const
bracket operator const access
Definition: facecentered/staggered/elementsolution.hh:94
static constexpr std::size_t size()
return the size of the element solution
Definition: facecentered/staggered/elementsolution.hh:106
FaceCenteredStaggeredElementSolution(const PrimaryVariables &priVars)
Constructor with a primary variable object.
Definition: facecentered/staggered/elementsolution.hh:75
void update(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
extract the element solution from the solution vector using a mapper
Definition: facecentered/staggered/elementsolution.hh:84
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
constexpr FCStaggered fcstaggered
Definition: method.hh:151
Definition: adapt.hh:17
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:29