version 3.11-dev
freeflow/navierstokes/mass/problem.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_NAVIERSTOKES_MASS_PROBLEM_HH
13#define DUMUX_NAVIERSTOKES_MASS_PROBLEM_HH
14
15#include <dune/common/exceptions.hh>
16#include <dune/common/typetraits.hh>
22
23namespace Dumux {
24
25// default implementation
26template<class TypeTag, class DiscretizationMethod>
28{
30 using Implementation = GetPropType<TypeTag, Properties::Problem>;
31
33 using GridView = typename GridGeometry::GridView;
34 using FVElementGeometry = typename GridGeometry::LocalView;
35 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
36 using Element = typename GridView::template Codim<0>::Entity;
37 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
38 using VelocityVector = GlobalPosition;
42
43 static constexpr bool isCoupled_ = !std::is_empty_v<CouplingManager>;
44
45public:
46
51 using InitialValues = Dune::FieldVector<Scalar, ModelTraits::numEq()>;
52 using Sources = Dune::FieldVector<Scalar, ModelTraits::numEq()>;
53 using DirichletValues = Dune::FieldVector<Scalar, ModelTraits::numEq()>;
54 using BoundaryFluxes = Dune::FieldVector<Scalar, ModelTraits::numEq()>;
55
57 using BoundaryTypes = Dumux::BoundaryTypes<ModelTraits::numEq()>;
58
60 static constexpr bool isMomentumProblem() { return false; }
61
68 NavierStokesMassProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
69 std::shared_ptr<CouplingManager> couplingManager,
70 const std::string& paramGroup = "")
72 , couplingManager_(couplingManager)
73 {}
74
80 NavierStokesMassProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
81 const std::string& paramGroup = "")
83 {}
84
88 VelocityVector faceVelocity(const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const SubControlVolumeFace& scvf) const
91 {
92 if constexpr (isCoupled_)
93 return couplingManager_->faceVelocity(element, scvf);
94 else
95 return asImp_().velocityAtPos(scvf.ipGlobal());
96 }
97
101 template<class IpData>
102 VelocityVector velocity(const FVElementGeometry& fvGeometry, const IpData& ipData) const
103 {
104 if constexpr (isCoupled_)
105 return couplingManager_->velocity(fvGeometry, ipData);
106 else
107 return asImp_().velocityAtPos(ipData.global());
108 }
109
113 VelocityVector elementVelocity(const FVElementGeometry& fvGeometry) const
114 {
115 if constexpr (isCoupled_)
116 return couplingManager_->elementVelocity(fvGeometry);
117 else
118 return asImp_().velocityAtPos(fvGeometry.element().geometry().center());
119 }
120
124 VelocityVector velocityAtPos(const GlobalPosition&) const
125 {
126 DUNE_THROW(Dune::NotImplemented, "velocityAtPos not implemented");
127 }
128
129 const CouplingManager& couplingManager() const
130 {
131 if constexpr (isCoupled_)
132 return *couplingManager_;
133 else
134 DUNE_THROW(Dune::InvalidStateException,
135 "Accessing coupling manager of an uncoupled problem is not possible."
136 );
137 }
138
139private:
141 Implementation &asImp_()
142 { return *static_cast<Implementation *>(this); }
143
145 const Implementation &asImp_() const
146 { return *static_cast<const Implementation *>(this); }
147
148 std::shared_ptr<CouplingManager> couplingManager_;
149};
150
151template<class TypeTag>
153{
155 using Implementation = GetPropType<TypeTag, Properties::Problem>;
156
158 using GridView = typename GridGeometry::GridView;
159 using FVElementGeometry = typename GridGeometry::LocalView;
160 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
161 using Element = typename GridView::template Codim<0>::Entity;
162 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
163 using VelocityVector = GlobalPosition;
167
168 static constexpr bool isCoupled_ = !std::is_empty_v<CouplingManager>;
169
170public:
171
176 using InitialValues = Dune::FieldVector<Scalar, ModelTraits::numEq()>;
177 using Sources = Dune::FieldVector<Scalar, ModelTraits::numEq()>;
178 using DirichletValues = Dune::FieldVector<Scalar, ModelTraits::numEq()>;
179 using BoundaryFluxes = Dune::FieldVector<Scalar, ModelTraits::numEq()>;
180
182 using BoundaryTypes = Dumux::BoundaryTypes<ModelTraits::numEq()>;
183
185 static constexpr bool isMomentumProblem() { return false; }
186 static constexpr bool providesIntegralInterface() { return true; }
187
194 CVFENavierStokesMassProblem(std::shared_ptr<const GridGeometry> gridGeometry,
195 std::shared_ptr<CouplingManager> couplingManager,
196 const std::string& paramGroup = "")
198 , couplingManager_(couplingManager)
199 {}
200
206 CVFENavierStokesMassProblem(std::shared_ptr<const GridGeometry> gridGeometry,
207 const std::string& paramGroup = "")
209 {}
210
211 using ParentType::source;
227 template<class ElementVariables, class IpData>
228 Sources source(const FVElementGeometry& fvGeometry,
229 const ElementVariables& elemVars,
230 const IpData& ipData) const
231 {
232 return asImp_().sourceAtPos(ipData.global());
233 }
234
242 template<class Intersection>
243 BoundaryTypes boundaryTypes(const FVElementGeometry& fvGeometry,
244 const Intersection& intersection) const
245 {
246 // forward it to the method which only takes the global coordinate
247 return asImp_().boundaryTypesAtPos(intersection.geometry().center());
248 }
249
258 template<class FaceIpData>
259 DirichletValues dirichlet(const FVElementGeometry& fvGeometry,
260 const FaceIpData& faceIpData) const
261 {
262 // forward it to the method which only takes the global coordinate
263 return asImp_().dirichletAtPos(faceIpData.global());
264 }
265
274 template<class ElementVariables, class ElementFluxVariablesCache>
275 BoundaryFluxes boundaryFluxIntegral(const FVElementGeometry& fvGeometry,
276 const ElementVariables& elemVars,
277 const ElementFluxVariablesCache& elemFluxVarsCache,
278 const SubControlVolumeFace& scvf) const
279 {
280 BoundaryFluxes flux(0.0);
281 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
282 flux += qpData.weight() * asImp_().boundaryFlux(fvGeometry, elemVars, elemFluxVarsCache, qpData.ipData());
283
284 return flux * elemVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
285 }
286
295 template<class ElementVariables, class ElementFluxVariablesCache, class FaceIpData>
296 BoundaryFluxes boundaryFlux(const FVElementGeometry& fvGeometry,
297 const ElementVariables& elemVars,
298 const ElementFluxVariablesCache& elemFluxVarsCache,
299 const FaceIpData& faceIpData) const
300 {
301 return asImp_().boundaryFluxAtPos(faceIpData.global());
302 }
303
307 BoundaryFluxes boundaryFluxAtPos(const GlobalPosition& globalPos) const
308 { return BoundaryFluxes(0.0); }
309
313 VelocityVector faceVelocity(const Element& element,
314 const FVElementGeometry& fvGeometry,
315 const SubControlVolumeFace& scvf) const
316 {
317 if constexpr (isCoupled_)
318 return couplingManager_->faceVelocity(element, scvf);
319 else
320 return asImp_().velocityAtPos(scvf.ipGlobal());
321 }
322
326 template<class IpData>
327 VelocityVector velocity(const FVElementGeometry& fvGeometry, const IpData& ipData) const
328 {
329 if constexpr (isCoupled_)
330 return couplingManager_->velocity(fvGeometry, ipData);
331 else
332 return asImp_().velocityAtPos(ipData.global());
333 }
334
338 VelocityVector elementVelocity(const FVElementGeometry& fvGeometry) const
339 {
340 if constexpr (isCoupled_)
341 return couplingManager_->elementVelocity(fvGeometry);
342 else
343 return asImp_().velocityAtPos(fvGeometry.element().geometry().center());
344 }
345
349 VelocityVector velocityAtPos(const GlobalPosition&) const
350 {
351 DUNE_THROW(Dune::NotImplemented, "velocityAtPos not implemented");
352 }
353
354 const CouplingManager& couplingManager() const
355 {
356 if constexpr (isCoupled_)
357 return *couplingManager_;
358 else
359 DUNE_THROW(Dune::InvalidStateException,
360 "Accessing coupling manager of an uncoupled problem is not possible."
361 );
362 }
363
364private:
366 Implementation &asImp_()
367 { return *static_cast<Implementation *>(this); }
368
370 const Implementation &asImp_() const
371 { return *static_cast<const Implementation *>(this); }
372
373 std::shared_ptr<CouplingManager> couplingManager_;
374};
375
382template<class TypeTag>
385>;
386
387} // end namespace Dumux
388
389#endif
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:26
Definition: freeflow/navierstokes/mass/problem.hh:153
CVFENavierStokesMassProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor for usage without a coupling manager.
Definition: freeflow/navierstokes/mass/problem.hh:206
BoundaryFluxes boundaryFluxIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluates the boundary flux integral for a scvf.
Definition: freeflow/navierstokes/mass/problem.hh:275
CVFENavierStokesMassProblem(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< CouplingManager > couplingManager, const std::string &paramGroup="")
The constructor.
Definition: freeflow/navierstokes/mass/problem.hh:194
Dune::FieldVector< Scalar, ModelTraits::numEq()> DirichletValues
Definition: freeflow/navierstokes/mass/problem.hh:178
VelocityVector faceVelocity(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
A default, i.e. if the user's does not overload any boundaryFlux method.
Definition: freeflow/navierstokes/mass/problem.hh:313
Dune::FieldVector< Scalar, ModelTraits::numEq()> Sources
Definition: freeflow/navierstokes/mass/problem.hh:177
Dune::FieldVector< Scalar, ModelTraits::numEq()> BoundaryFluxes
Definition: freeflow/navierstokes/mass/problem.hh:179
static constexpr bool isMomentumProblem()
this problem is used for the mass balance model
Definition: freeflow/navierstokes/mass/problem.hh:185
VelocityVector velocity(const FVElementGeometry &fvGeometry, const IpData &ipData) const
Returns the velocity at a given interpolation point.
Definition: freeflow/navierstokes/mass/problem.hh:327
BoundaryFluxes boundaryFlux(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const FaceIpData &faceIpData) const
Evaluates the boundary flux related to a localDof at a given interpolation point.
Definition: freeflow/navierstokes/mass/problem.hh:296
DirichletValues dirichlet(const FVElementGeometry &fvGeometry, const FaceIpData &faceIpData) const
Evaluate the boundary conditions for a Dirichlet control volume.
Definition: freeflow/navierstokes/mass/problem.hh:259
Sources source(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const IpData &ipData) const
Evaluate the source term at a given interpolation point, related to the residual of a local dof.
Definition: freeflow/navierstokes/mass/problem.hh:228
static constexpr bool providesIntegralInterface()
Definition: freeflow/navierstokes/mass/problem.hh:186
Dune::FieldVector< Scalar, ModelTraits::numEq()> InitialValues
Definition: freeflow/navierstokes/mass/problem.hh:176
BoundaryFluxes boundaryFluxAtPos(const GlobalPosition &globalPos) const
Returns the boundary flux at a given position.
Definition: freeflow/navierstokes/mass/problem.hh:307
VelocityVector velocityAtPos(const GlobalPosition &) const
Returns the velocity at a given position.
Definition: freeflow/navierstokes/mass/problem.hh:349
const CouplingManager & couplingManager() const
Definition: freeflow/navierstokes/mass/problem.hh:354
BoundaryTypes boundaryTypes(const FVElementGeometry &fvGeometry, const Intersection &intersection) const
Specifies which kind of boundary condition should be used for which equation on a given boundary face...
Definition: freeflow/navierstokes/mass/problem.hh:243
VelocityVector elementVelocity(const FVElementGeometry &fvGeometry) const
Returns the velocity at the element center.
Definition: freeflow/navierstokes/mass/problem.hh:338
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:43
const std::string & paramGroup() const
The parameter group in which to retrieve runtime parameters.
Definition: common/fvproblem.hh:524
PrimaryVariables dirichlet(const Element &element, const SubControlVolumeFace &scvf) const
Evaluate the boundary conditions for a dirichlet control volume face.
Definition: common/fvproblem.hh:178
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:520
NumEqVector source(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: common/fvproblem.hh:307
Base class for all finite-volume problems using spatial parameters.
Definition: fvproblemwithspatialparams.hh:29
Definition: freeflow/navierstokes/mass/problem.hh:28
VelocityVector elementVelocity(const FVElementGeometry &fvGeometry) const
Returns the velocity at the element center.
Definition: freeflow/navierstokes/mass/problem.hh:113
Dune::FieldVector< Scalar, ModelTraits::numEq()> BoundaryFluxes
Definition: freeflow/navierstokes/mass/problem.hh:54
const CouplingManager & couplingManager() const
Definition: freeflow/navierstokes/mass/problem.hh:129
VelocityVector faceVelocity(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the normal velocity at a given sub control volume face.
Definition: freeflow/navierstokes/mass/problem.hh:88
VelocityVector velocity(const FVElementGeometry &fvGeometry, const IpData &ipData) const
Returns the velocity at a given interpolation point.
Definition: freeflow/navierstokes/mass/problem.hh:102
Dune::FieldVector< Scalar, ModelTraits::numEq()> DirichletValues
Definition: freeflow/navierstokes/mass/problem.hh:53
NavierStokesMassProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor for usage without a coupling manager.
Definition: freeflow/navierstokes/mass/problem.hh:80
Dune::FieldVector< Scalar, ModelTraits::numEq()> Sources
Definition: freeflow/navierstokes/mass/problem.hh:52
NavierStokesMassProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< CouplingManager > couplingManager, const std::string &paramGroup="")
The constructor.
Definition: freeflow/navierstokes/mass/problem.hh:68
VelocityVector velocityAtPos(const GlobalPosition &) const
Returns the velocity at a given position.
Definition: freeflow/navierstokes/mass/problem.hh:124
Dune::FieldVector< Scalar, ModelTraits::numEq()> InitialValues
Definition: freeflow/navierstokes/mass/problem.hh:51
static constexpr bool isMomentumProblem()
this problem is used for the mass balance model
Definition: freeflow/navierstokes/mass/problem.hh:60
Defines all properties used in Dumux.
Base class for all finite volume problems that are parameterized.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition: quadraturerules.hh:148
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
Quadrature rules over sub-control volumes and sub-control volume faces.