version 3.10.0
Loading...
Searching...
No Matches
helmholtzoperator.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_LINEAR_HELMHOLTZ_OPERATOR_HH
13#define DUMUX_LINEAR_HELMHOLTZ_OPERATOR_HH
14
15#include <dumux/common/math.hh>
23
25
26template <class Traits>
28{
29 using Scalar = typename Traits::PrimaryVariables::value_type;
30public:
31 using PrimaryVariables = typename Traits::PrimaryVariables;
32
33 template<class ElementSolution, class Problem, class Element, class SubControlVolume>
34 void update(const ElementSolution& elemSol, const Problem&, const Element&, const SubControlVolume& scv)
35 { priVars_ = elemSol[scv.indexInElement()]; }
36
37 Scalar priVar(const int pvIdx) const { return priVars_[pvIdx]; }
38 const PrimaryVariables& priVars() const { return priVars_; }
39 Scalar extrusionFactor() const { return 1.0; }
40
41private:
42 PrimaryVariables priVars_;
43};
44
45template<class TypeTag>
48{
54 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
55 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
57 using FVElementGeometry = typename GridGeometry::LocalView;
58 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
59 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
60 using GridView = typename GridGeometry::GridView;
61 using Element = typename GridView::template Codim<0>::Entity;
62 static constexpr int dimWorld = GridView::dimensionworld;
63public:
64 using ParentType::ParentType;
65
66 NumEqVector computeStorage(const Problem&,
67 const SubControlVolume&,
68 const VolumeVariables&) const
69 { return NumEqVector(0.0); }
70
71 NumEqVector computeFlux(const Problem& problem,
72 const Element&, const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars,
74 const SubControlVolumeFace& scvf,
75 const ElementFluxVariablesCache& elemFluxVarsCache) const
76 {
77 NumEqVector flux(0.0);
78
79 // CVFE schemes
80 if constexpr (GridGeometry::discMethod == DiscretizationMethods::box
81 || GridGeometry::discMethod == DiscretizationMethods::fcdiamond
82 || GridGeometry::discMethod == DiscretizationMethods::pq1bubble)
83 {
84 const auto& fluxVarCache = elemFluxVarsCache[scvf];
85 Dune::FieldVector<Scalar, dimWorld> gradConcentration(0.0);
86 for (const auto& scv : scvs(fvGeometry))
87 {
88 const auto& volVars = elemVolVars[scv];
89 gradConcentration.axpy(volVars.concentration(), fluxVarCache.gradN(scv.indexInElement()));
90 }
91
92 NumEqVector flux;
93 flux[0] = -1.0*vtmv(scvf.unitOuterNormal(), problem.a(), gradConcentration)*scvf.area();
94 }
95
96 // CCTpfa
97 else if constexpr (GridGeometry::discMethod == DiscretizationMethods::cctpfa)
98 {
99 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
100 const auto& insideV = elemVolVars[insideScv];
101 const auto& outsideV = elemVolVars[scvf.outsideScvIdx()];
102 const auto ti = computeTpfaTransmissibility(
103 fvGeometry, scvf, insideScv, problem.a(), insideV.extrusionFactor()
104 );
105
106 const auto tij = [&]{
107 if (scvf.boundary())
108 return scvf.area()*ti;
109 else
110 {
111 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
112 const Scalar tj = -computeTpfaTransmissibility(
113 fvGeometry, scvf, outsideScv, problem.a(), outsideV.extrusionFactor()
114 );
115
116 return scvf.area()*(ti * tj)/(ti + tj);
117 }
118 }();
119
120 const auto u = insideV.priVar(0);
121 const auto v = outsideV.priVar(0);
122 flux[0] = tij*(u - v);
123 }
124 else
125 DUNE_THROW(Dune::NotImplemented, "Discretization method " << GridGeometry::discMethod);
126
127 return flux;
128 }
129
130 NumEqVector computeSource(const Problem& problem,
131 const Element&, const FVElementGeometry&,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolume& scv) const
134 {
135 NumEqVector source(0.0);
136 source[0] = -problem.b()*elemVolVars[scv].priVar(0);
137 return source;
138 }
139};
140
141template<class TypeTag>
143{
144 using ParentType = FVProblem<TypeTag>;
148public:
149 HelmholtzModelHomogeneousNeumannProblem(std::shared_ptr<const GridGeometry> gridGeometry, Scalar a, Scalar b)
150 : ParentType(gridGeometry)
151 , a_(a), b_(b)
152 {}
153
154 template<class GlobalPosition>
155 BoundaryTypes boundaryTypesAtPos(const GlobalPosition& globalPos) const
156 { BoundaryTypes values; values.setAllNeumann(); return values; }
157
158 Scalar a() const { return a_; }
159 Scalar b() const { return b_; }
160private:
161 Scalar a_, b_;
162};
163
164template<class G, class D, class S>
165struct Policy
166{
167 using Scalar = S;
168 using Grid = G;
169 using Discretization = D;
170};
171
172template<class P>
173struct TTag {
174 using InheritsFrom = std::tuple<typename P::Discretization>;
175};
176
177} // end namespace Dumux::Detail::HelmholtzOperator
178
179namespace Dumux::Properties {
180
182template<class TypeTag, class P>
183struct Scalar<TypeTag, Dumux::Detail::HelmholtzOperator::TTag<P>>
184{ using type = typename P::Scalar; };
185
187template<class TypeTag, class P>
188struct PrimaryVariables<TypeTag, Dumux::Detail::HelmholtzOperator::TTag<P>>
189{
190 using type = Dune::FieldVector<
193 >;
194};
195
196template<class TypeTag, class P>
197struct ModelTraits<TypeTag, Dumux::Detail::HelmholtzOperator::TTag<P>>
198{
199 struct Traits { static constexpr int numEq() { return 1; } };
200 using type = Traits;
201};
202
203template<class TypeTag, class P>
206
207template<class TypeTag, class P>
210
211template<class TypeTag, class P>
217
218template<class TypeTag, class P>
219struct EnableGridVolumeVariablesCache<TypeTag, Dumux::Detail::HelmholtzOperator::TTag<P>>
220{ static constexpr bool value = true; };
221
222template<class TypeTag, class P>
223struct EnableGridFluxVariablesCache<TypeTag, Dumux::Detail::HelmholtzOperator::TTag<P>>
224{ static constexpr bool value = true; };
225
226template<class TypeTag, class P>
227struct EnableGridGeometryCache<TypeTag, Dumux::Detail::HelmholtzOperator::TTag<P>>
228{ static constexpr bool value = true; };
229
230template<class TypeTag, class P>
232{ using type = typename P::Grid; };
233
234} // end namespace Dumux::Properties
235
236namespace Dumux {
237
251template<class Discretization, class GridView, class Scalar>
252auto makeHelmholtzMatrix(const GridView& gridView, const Scalar a = 1.0, const Scalar b = 1.0)
253{
254 using TypeTag = Detail::HelmholtzOperator::TTag<
256 >;
258 auto gridGeometry = std::make_shared<GridGeometry>(gridView);
259
261 auto problem = std::make_shared<Problem>(gridGeometry, a, b);
262
265 auto gridVariables = std::make_shared<GridVariables>(problem, gridGeometry);
266 SolutionVector x(gridGeometry->numDofs());
267 gridVariables->init(x);
268
270 auto assembler = std::make_shared<Assembler>(problem, gridGeometry, gridVariables);
271 using A = typename Assembler::JacobianMatrix;
272 using V = typename Assembler::ResidualType;
273 auto jacobian = std::make_shared<A>();
274 auto residual = std::make_shared<V>();
275 assembler->setLinearSystem(jacobian, residual);
276 assembler->assembleJacobian(x);
277 return jacobian;
278};
279
283template<class Discretization, class GridView, class Scalar>
284auto makeLaplaceMatrix(const GridView& gridView, const Scalar a = 1.0)
285{ return makeHelmholtzMatrix<Discretization>(gridView, a, 0.0); };
286
287template<class LinearOperator, class Discretization, class GridView, class Scalar>
288std::shared_ptr<LinearOperator> makeHelmholtzLinearOperator(const GridView& gridView, const Scalar a, const Scalar b)
289{ return std::make_shared<LinearOperator>(makeHelmholtzMatrix<Discretization>(gridView, a, b)); }
290
291template<class LinearOperator, class Discretization, class GridView, class Comm, class Scalar>
292std::shared_ptr<LinearOperator> makeHelmholtzLinearOperator(const GridView& gridView, const Comm& comm, const Scalar a, const Scalar b)
293{ return std::make_shared<LinearOperator>(makeHelmholtzMatrix<Discretization>(gridView, a, b), comm); }
294
295template<class LinearOperator, class Discretization, class GridView, class Scalar>
296std::shared_ptr<LinearOperator> makeLaplaceLinearOperator(const GridView& gridView, const Scalar a)
297{ return std::make_shared<LinearOperator>(makeLaplaceMatrix<Discretization>(gridView, a)); }
298
299template<class LinearOperator, class Discretization, class GridView, class Comm, class Scalar>
300std::shared_ptr<LinearOperator> makeLaplaceLinearOperator(const GridView& gridView, const Comm& comm, const Scalar a)
301{ return std::make_shared<LinearOperator>(makeLaplaceMatrix<Discretization>(gridView, a), comm); }
302
303} // end namespace Dumux
304
305#endif
A linear system assembler (residual and Jacobian) for finite volume schemes.
Class to specify the type of a boundary.
Definition common/boundarytypes.hh:26
Scalar b() const
Definition helmholtzoperator.hh:159
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Definition helmholtzoperator.hh:155
Scalar a() const
Definition helmholtzoperator.hh:158
HelmholtzModelHomogeneousNeumannProblem(std::shared_ptr< const GridGeometry > gridGeometry, Scalar a, Scalar b)
Definition helmholtzoperator.hh:149
NumEqVector computeStorage(const Problem &, const SubControlVolume &, const VolumeVariables &) const
Definition helmholtzoperator.hh:66
NumEqVector computeSource(const Problem &problem, const Element &, const FVElementGeometry &, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Definition helmholtzoperator.hh:130
NumEqVector computeFlux(const Problem &problem, const Element &, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Definition helmholtzoperator.hh:71
typename Traits::PrimaryVariables PrimaryVariables
Definition helmholtzoperator.hh:31
void update(const ElementSolution &elemSol, const Problem &, const Element &, const SubControlVolume &scv)
Definition helmholtzoperator.hh:34
const PrimaryVariables & priVars() const
Definition helmholtzoperator.hh:38
Scalar priVar(const int pvIdx) const
Definition helmholtzoperator.hh:37
Scalar extrusionFactor() const
Definition helmholtzoperator.hh:39
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition assembly/fvassembler.hh:95
FVProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
Constructor.
Definition common/fvproblem.hh:84
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition common/fvproblem.hh:520
Class to specify the type of a boundary.
Base class for all finite volume problems.
Defines all properties used in Dumux.
The default local operator than can be specialized for each discretization scheme.
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition tpfa/computetransmissibility.hh:36
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition math.hh:880
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition numeqvector.hh:34
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296
Define some often used mathematical functions.
Definition helmholtzoperator.hh:24
constexpr FCDiamond fcdiamond
Definition method.hh:152
constexpr CCTpfa cctpfa
Definition method.hh:145
constexpr Box box
Definition method.hh:147
constexpr PQ1Bubble pq1bubble
Definition method.hh:148
Definition gridcapabilities.hh:57
The energy balance equation for a porous solid.
Definition common/properties.hh:26
Definition adapt.hh:17
std::shared_ptr< LinearOperator > makeLaplaceLinearOperator(const GridView &gridView, const Scalar a)
Definition helmholtzoperator.hh:296
auto makeHelmholtzMatrix(const GridView &gridView, const Scalar a=1.0, const Scalar b=1.0)
make a Helmholtz matrix operator (aΔ + bI)
Definition helmholtzoperator.hh:252
auto makeLaplaceMatrix(const GridView &gridView, const Scalar a=1.0)
make a Laplace matrix operator aΔ
Definition helmholtzoperator.hh:284
const Scalar PengRobinsonMixture< Scalar, StaticParameters >::u
Definition pengrobinsonmixture.hh:138
typename Detail::DiscretizationDefaultLocalOperator< TypeTag >::type DiscretizationDefaultLocalOperator
Definition defaultlocaloperator.hh:26
std::shared_ptr< LinearOperator > makeHelmholtzLinearOperator(const GridView &gridView, const Scalar a, const Scalar b)
Definition helmholtzoperator.hh:288
A helper to deduce a vector with the same size as numbers of equations.
Definition helmholtzoperator.hh:166
G Grid
Definition helmholtzoperator.hh:168
D Discretization
Definition helmholtzoperator.hh:169
S Scalar
Definition helmholtzoperator.hh:167
Definition helmholtzoperator.hh:173
std::tuple< typename P::Discretization > InheritsFrom
Definition helmholtzoperator.hh:174
Dumux::Detail::HelmholtzOperator::HelmholtzModelLocalResidual< TypeTag > type
Definition helmholtzoperator.hh:205
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::ModelTraits >::numEq() > type
Definition helmholtzoperator.hh:190
Dumux::Detail::HelmholtzOperator::HelmholtzModelHomogeneousNeumannProblem< TypeTag > type
Definition helmholtzoperator.hh:209
typename P::Scalar type
Definition helmholtzoperator.hh:184
Dumux::Detail::HelmholtzOperator::HelmholtzModelVolumeVariables< Traits > type
Definition helmholtzoperator.hh:215
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition helmholtzoperator.hh:214
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...