12#ifndef DUMUX_LINEAR_HELMHOLTZ_OPERATOR_HH
13#define DUMUX_LINEAR_HELMHOLTZ_OPERATOR_HH
26template <
class Traits>
29 using Scalar =
typename Traits::PrimaryVariables::value_type;
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()]; }
37 Scalar
priVar(
const int pvIdx)
const {
return priVars_[pvIdx]; }
45template<
class TypeTag>
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;
64 using ParentType::ParentType;
67 const SubControlVolume&,
68 const VolumeVariables&)
const
69 {
return NumEqVector(0.0); }
72 const Element&,
const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars,
74 const SubControlVolumeFace& scvf,
75 const ElementFluxVariablesCache& elemFluxVarsCache)
const
77 NumEqVector flux(0.0);
85 const auto& fluxVarCache = elemFluxVarsCache[scvf];
86 Dune::FieldVector<Scalar, dimWorld> gradConcentration(0.0);
87 for (
const auto& scv :
scvs(fvGeometry))
89 const auto& volVars = elemVolVars[scv];
90 gradConcentration.axpy(volVars.concentration(), fluxVarCache.gradN(scv.indexInElement()));
94 flux[0] = -1.0*
vtmv(scvf.unitOuterNormal(), problem.a(), gradConcentration)*scvf.area();
100 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
101 const auto& insideV = elemVolVars[insideScv];
102 const auto& outsideV = elemVolVars[scvf.outsideScvIdx()];
104 fvGeometry, scvf, insideScv, problem.a(), insideV.extrusionFactor()
107 const auto tij = [&]{
109 return scvf.area()*ti;
112 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
114 fvGeometry, scvf, outsideScv, problem.a(), outsideV.extrusionFactor()
117 return scvf.area()*(ti * tj)/(ti + tj);
121 const auto u = insideV.priVar(0);
122 const auto v = outsideV.priVar(0);
123 flux[0] = tij*(u - v);
126 DUNE_THROW(Dune::NotImplemented,
"Discretization method " << GridGeometry::discMethod);
132 const Element&,
const FVElementGeometry&,
133 const ElementVolumeVariables& elemVolVars,
134 const SubControlVolume& scv)
const
136 NumEqVector source(0.0);
137 source[0] = -problem.b()*elemVolVars[scv].priVar(0);
142template<
class TypeTag>
155 template<
class GlobalPosition>
159 Scalar
a()
const {
return a_; }
160 Scalar
b()
const {
return b_; }
165template<
class G,
class D,
class S>
183template<
class TypeTag,
class P>
185{
using type =
typename P::Scalar; };
188template<
class TypeTag,
class P>
191 using type = Dune::FieldVector<
197template<
class TypeTag,
class P>
200 struct Traits {
static constexpr int numEq() {
return 1; } };
204template<
class TypeTag,
class P>
208template<
class TypeTag,
class P>
212template<
class TypeTag,
class P>
219template<
class TypeTag,
class P>
221{
static constexpr bool value =
true; };
223template<
class TypeTag,
class P>
225{
static constexpr bool value =
true; };
227template<
class TypeTag,
class P>
229{
static constexpr bool value =
true; };
231template<
class TypeTag,
class P>
233{
using type =
typename P::Grid; };
252template<
class Discretization,
class Gr
idView,
class Scalar>
259 auto gridGeometry = std::make_shared<GridGeometry>(gridView);
262 auto problem = std::make_shared<Problem>(gridGeometry, a, b);
266 auto gridVariables = std::make_shared<GridVariables>(problem, gridGeometry);
267 SolutionVector x(gridGeometry->numDofs());
268 gridVariables->init(x);
271 auto assembler = std::make_shared<Assembler>(problem, gridGeometry, gridVariables);
272 using A =
typename Assembler::JacobianMatrix;
273 using V =
typename Assembler::ResidualType;
274 auto jacobian = std::make_shared<A>();
275 auto residual = std::make_shared<V>();
276 assembler->setLinearSystem(jacobian, residual);
277 assembler->assembleJacobian(x);
284template<
class Discretization,
class Gr
idView,
class Scalar>
286{
return makeHelmholtzMatrix<Discretization>(gridView, a, 0.0); };
288template<
class LinearOperator,
class Discretization,
class Gr
idView,
class Scalar>
290{
return std::make_shared<LinearOperator>(makeHelmholtzMatrix<Discretization>(gridView, a, b)); }
292template<
class LinearOperator,
class Discretization,
class Gr
idView,
class Comm,
class Scalar>
294{
return std::make_shared<LinearOperator>(makeHelmholtzMatrix<Discretization>(gridView, a, b), comm); }
296template<
class LinearOperator,
class Discretization,
class Gr
idView,
class Scalar>
298{
return std::make_shared<LinearOperator>(makeLaplaceMatrix<Discretization>(gridView, a)); }
300template<
class LinearOperator,
class Discretization,
class Gr
idView,
class Comm,
class Scalar>
302{
return std::make_shared<LinearOperator>(makeLaplaceMatrix<Discretization>(gridView, a), comm); }
A linear system assembler (residual and Jacobian) for finite volume schemes.
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:26
void setAllNeumann()
Set all boundary conditions to Neumann.
Definition: common/boundarytypes.hh:90
Definition: helmholtzoperator.hh:144
Scalar b() const
Definition: helmholtzoperator.hh:160
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Definition: helmholtzoperator.hh:156
Scalar a() const
Definition: helmholtzoperator.hh:159
HelmholtzModelHomogeneousNeumannProblem(std::shared_ptr< const GridGeometry > gridGeometry, Scalar a, Scalar b)
Definition: helmholtzoperator.hh:150
Definition: helmholtzoperator.hh:48
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:131
NumEqVector computeFlux(const Problem &problem, const Element &, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Definition: helmholtzoperator.hh:71
Definition: helmholtzoperator.hh:28
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:103
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:43
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:162
constexpr PQ2 pq2
Definition: method.hh:157
constexpr CCTpfa cctpfa
Definition: method.hh:154
constexpr Box box
Definition: method.hh:156
constexpr PQ1Bubble pq1bubble
Definition: method.hh:158
The energy balance equation for a porous solid.
Definition: common/properties.hh:26
std::shared_ptr< LinearOperator > makeLaplaceLinearOperator(const GridView &gridView, const Scalar a)
Definition: helmholtzoperator.hh:297
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:253
auto makeLaplaceMatrix(const GridView &gridView, const Scalar a=1.0)
make a Laplace matrix operator aΔ
Definition: helmholtzoperator.hh:285
typename Detail::DiscretizationDefaultLocalOperator< TypeTag >::type DiscretizationDefaultLocalOperator
Definition: defaultlocaloperator.hh:27
std::shared_ptr< LinearOperator > makeHelmholtzLinearOperator(const GridView &gridView, const Scalar a, const Scalar b)
Definition: helmholtzoperator.hh:289
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
A helper to deduce a vector with the same size as numbers of equations.
Definition: helmholtzoperator.hh:167
G Grid
Definition: helmholtzoperator.hh:169
D Discretization
Definition: helmholtzoperator.hh:170
S Scalar
Definition: helmholtzoperator.hh:168
Definition: helmholtzoperator.hh:174
std::tuple< typename P::Discretization > InheritsFrom
Definition: helmholtzoperator.hh:175
typename P::Grid type
Definition: helmholtzoperator.hh:233
Traits type
Definition: helmholtzoperator.hh:201
static constexpr int numEq()
Definition: helmholtzoperator.hh:200
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::ModelTraits >::numEq() > type
Definition: helmholtzoperator.hh:194
typename P::Scalar type
Definition: helmholtzoperator.hh:185
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: helmholtzoperator.hh:215
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...