24#ifndef DUMUX_PQ1BUBBLE_ELEMENT_SOLUTION_HH
25#define DUMUX_PQ1BUBBLE_ELEMENT_SOLUTION_HH
28#include <dune/common/reservedvector.hh>
37template<
class FVElementGeometry,
class PV>
40 using GridGeometry =
typename FVElementGeometry::GridGeometry;
41 using GridView =
typename GridGeometry::GridView;
42 using Element =
typename GridView::template Codim<0>::Entity;
44 static constexpr int dim = GridView::dimension;
46 static constexpr int numCubeDofs = (1<<dim) + 1;
56 template<
class SolutionVector>
58 const GridGeometry& gridGeometry)
60 update(element, sol, gridGeometry);
64 template<
class ElementVolumeVariables>
66 const FVElementGeometry& fvGeometry)
68 priVars_.resize(fvGeometry.numScv());
69 for (
const auto& scv : scvs(fvGeometry))
70 priVars_[scv.localDofIndex()] = elemVolVars[scv].priVars();
74 template<
class SolutionVector>
75 void update(
const Element& element,
const SolutionVector& sol,
76 const GridGeometry& gridGeometry)
78 const auto numVertices = element.subEntities(GridView::dimension);
79 priVars_.resize(numVertices+1);
80 for (
int vIdx = 0; vIdx < numVertices; ++vIdx)
81 priVars_[vIdx] = sol[gridGeometry.dofMapper().subIndex(element, vIdx, GridView::dimension)];
83 priVars_[numVertices] = sol[gridGeometry.dofMapper().index(element)];
87 template<
class SolutionVector>
88 void update(
const Element& element,
const SolutionVector& sol,
89 const FVElementGeometry& fvGeometry)
91 priVars_.resize(fvGeometry.numScv());
92 for (
const auto& scv : scvs(fvGeometry))
93 priVars_[scv.localDofIndex()] = sol[scv.dofIndex()];
98 {
return priVars_.size(); }
101 template<
typename IndexType>
103 {
return priVars_[i]; }
106 template<
typename IndexType>
108 {
return priVars_[i]; }
111 Dune::ReservedVector<PrimaryVariables, numCubeDofs> priVars_;
118template<
class Element,
class SolutionVector,
class Gr
idGeometry>
119auto elementSolution(
const Element& element,
const SolutionVector& sol,
const GridGeometry& gg)
121 PQ1BubbleElementSolution<
typename GridGeometry::LocalView,
122 std::decay_t<decltype(std::declval<SolutionVector>()[0])>>
125 using PrimaryVariables = std::decay_t<decltype(std::declval<SolutionVector>()[0])>;
126 return PQ1BubbleElementSolution<typename GridGeometry::LocalView, PrimaryVariables>(element, sol, gg);
133template<
class Element,
class ElementVolumeVariables,
class FVElementGeometry>
134auto elementSolution(
const Element& element,
const ElementVolumeVariables& elemVolVars,
const FVElementGeometry& gg)
136 PQ1BubbleElementSolution<FVElementGeometry,
137 typename ElementVolumeVariables::VolumeVariables::PrimaryVariables>>
139 using PrimaryVariables =
typename ElementVolumeVariables::VolumeVariables::PrimaryVariables;
140 return PQ1BubbleElementSolution<FVElementGeometry, PrimaryVariables>(element, elemVolVars, gg);
The available discretization methods in Dumux.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr PQ1Bubble pq1bubble
Definition: method.hh:137
The element solution vector.
Definition: pq1bubble/elementsolution.hh:39
PV PrimaryVariables
export the primary variables type
Definition: pq1bubble/elementsolution.hh:50
PQ1BubbleElementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
Constructor with element and solution and grid geometry.
Definition: pq1bubble/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: pq1bubble/elementsolution.hh:88
PQ1BubbleElementSolution(const Element &element, const ElementVolumeVariables &elemVolVars, const FVElementGeometry &fvGeometry)
Constructor with element and elemVolVars and fvGeometry.
Definition: pq1bubble/elementsolution.hh:65
PQ1BubbleElementSolution()=default
Default constructor.
const PrimaryVariables & operator[](IndexType i) const
bracket operator const access
Definition: pq1bubble/elementsolution.hh:102
void update(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
extract the element solution from the solution vector using a mapper
Definition: pq1bubble/elementsolution.hh:75
std::size_t size() const
return the size of the element solution
Definition: pq1bubble/elementsolution.hh:97