24#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_PQ1BUBBLE_HH
25#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_PQ1BUBBLE_HH
32#include <dune/common/exceptions.hh>
33#include <dune/common/indices.hh>
34#include <dune/common/float_cmp.hh>
35#include <dune/geometry/referenceelements.hh>
70 template<std::
size_t id>
using SubDomainTypeTag =
typename Traits::template SubDomain<id>::TypeTag;
73 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
74 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
75 template<std::
size_t id>
using ElementSeed =
typename GridView<id>::Grid::template Codim<0>::EntitySeed;
76 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
77 template<std::
size_t id>
using SubControlVolume =
typename FVElementGeometry<id>::SubControlVolume;
78 template<std::
size_t id>
using SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
79 template<std::
size_t id>
using GridVariables =
typename Traits::template SubDomain<id>::GridVariables;
80 template<std::
size_t id>
using ElementVolumeVariables =
typename GridVariables<id>::GridVolumeVariables::LocalView;
81 template<std::
size_t id>
using GridFluxVariablesCache =
typename GridVariables<id>::GridFluxVariablesCache;
85 using Scalar =
typename Traits::Scalar;
86 using SolutionVector =
typename Traits::SolutionVector;
90 using GridVariablesTuple =
typename Traits::template TupleOfSharedPtr<GridVariables>;
92 using FluidSystem =
typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
94 using VelocityVector =
typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
95 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
97 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
99 struct MomentumCouplingContext
101 FVElementGeometry<freeFlowMassIndex> fvGeometry;
102 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
103 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
107 struct MassAndEnergyCouplingContext
109 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f,
const std::size_t i)
110 : fvGeometry(std::move(f))
114 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
118 using MomentumDiscretizationMethod =
typename GridGeometry<freeFlowMomentumIndex>::DiscretizationMethod;
119 using MassDiscretizationMethod =
typename GridGeometry<freeFlowMassIndex>::DiscretizationMethod;
123 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
131 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
132 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
133 GridVariablesTuple&& gridVariables,
134 const SolutionVector&
curSol)
136 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
137 gridVariables_ = gridVariables;
140 computeCouplingStencils_();
144 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
145 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
146 GridVariablesTuple&& gridVariables,
147 const SolutionVector&
curSol,
148 const SolutionVector& prevSol)
150 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables),
curSol);
156 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
157 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
158 GridVariablesTuple&& gridVariables,
161 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
162 gridVariables_ = gridVariables;
165 computeCouplingStencils_();
178 Scalar
pressure(
const Element<freeFlowMomentumIndex>& element,
179 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
180 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
190 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
191 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
192 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
193 const bool considerPreviousTimeStep =
false)
const
195 assert(!(considerPreviousTimeStep && !this->isTransient_));
196 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
200 const auto eIdx = fvGeometry.elementIndex();
201 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
203 const auto& volVars = considerPreviousTimeStep ?
204 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
205 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
207 return volVars.density();
213 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
214 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
215 std::vector<ShapeValue> shapeValues;
216 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
217 localBasis.evaluateFunction(ipLocal, shapeValues);
220 for (
const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
222 const auto& volVars = considerPreviousTimeStep ?
223 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
224 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
225 rho += volVars.density()*shapeValues[scv.indexInElement()][0];
231 DUNE_THROW(Dune::NotImplemented,
232 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
239 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
240 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
241 const SubControlVolume<freeFlowMomentumIndex>& scv,
242 const bool considerPreviousTimeStep =
false)
const
244 assert(!(considerPreviousTimeStep && !this->isTransient_));
245 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
249 const auto eIdx = scv.elementIndex();
250 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
252 const auto& volVars = considerPreviousTimeStep ?
253 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
254 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
256 return volVars.density();
262 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
263 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
264 std::vector<ShapeValue> shapeValues;
265 const auto ipLocal = element.geometry().local(scv.dofPosition());
266 localBasis.evaluateFunction(ipLocal, shapeValues);
269 for (
const auto& scvI : scvs(this->momentumCouplingContext_()[0].fvGeometry))
271 const auto& volVars = considerPreviousTimeStep ?
272 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
273 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
274 rho += volVars.density()*shapeValues[scvI.indexInElement()][0];
279 DUNE_THROW(Dune::NotImplemented,
280 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
288 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
289 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
291 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
295 const auto eIdx = fvGeometry.elementIndex();
296 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
297 const auto& volVars = this->momentumCouplingContext_()[0].curElemVolVars[scv];
298 return volVars.viscosity();
304 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
305 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
306 std::vector<ShapeValue> shapeValues;
307 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
308 localBasis.evaluateFunction(ipLocal, shapeValues);
311 for (
const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
313 const auto& volVars = this->momentumCouplingContext_()[0].curElemVolVars[scv];
314 mu += volVars.viscosity()*shapeValues[scv.indexInElement()][0];
320 DUNE_THROW(Dune::NotImplemented,
321 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
330 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
334 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element);
336 const auto& fvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
337 const auto& localBasis = fvGeometry.feLocalBasis();
340 VelocityVector velocity(0.0);
341 std::vector<ShapeValue> shapeValues;
342 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
343 localBasis.evaluateFunction(ipLocal, shapeValues);
345 for (
const auto& scv : scvs(fvGeometry))
346 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(
freeFlowMomentumIndex)[scv.dofIndex()]);
354 VelocityVector
elementVelocity(
const FVElementGeometry<freeFlowMassIndex>& fvGeometry)
const
356 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), fvGeometry.element());
358 const auto& momentumFvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
359 const auto& localBasis = momentumFvGeometry.feLocalBasis();
362 VelocityVector velocity(0.0);
363 std::vector<ShapeValue> shapeValues;
364 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
366 for (
const auto& scv : scvs(momentumFvGeometry))
367 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(
freeFlowMomentumIndex)[scv.dofIndex()]);
376 template<std::
size_t j>
378 const Element<freeFlowMomentumIndex>& elementI,
379 const SubControlVolume<freeFlowMomentumIndex>& scvI,
380 Dune::index_constant<j> domainJ)
const
381 {
return emptyStencil_; }
398 const Element<freeFlowMassIndex>& elementI,
399 Dune::index_constant<freeFlowMomentumIndex> domainJ)
const
402 return massAndEnergyToMomentumStencils_[eIdx];
414 const Element<freeFlowMomentumIndex>& elementI,
415 Dune::index_constant<freeFlowMassIndex> domainJ)
const
418 return momentumToMassAndEnergyStencils_[eIdx];
447 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
449 const LocalAssemblerI& localAssemblerI,
450 Dune::index_constant<j> domainJ,
451 std::size_t dofIdxGlobalJ,
452 const PrimaryVariables<j>& priVarsJ,
455 this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
461 bindCouplingContext_(domainI, localAssemblerI.element());
464 const auto& deflectedElement =
problem.gridGeometry().element(dofIdxGlobalJ);
466 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
467 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
469 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
470 gridVars_(
freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol),
problem, deflectedElement, scv);
472 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol),
problem, deflectedElement, scv);
480 bindCouplingContext_(domainI, localAssemblerI.element());
483 const auto& deflectedElement =
problem.gridGeometry().element(this->momentumCouplingContext_()[0].eIdx);
485 const auto& fvGeometry = this->momentumCouplingContext_()[0].fvGeometry;
487 for (
const auto& scv : scvs(fvGeometry))
489 if(scv.dofIndex() == dofIdxGlobalJ)
491 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
492 this->gridVars_(
freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol),
problem, deflectedElement, scv);
494 this->momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol),
problem, deflectedElement, scv);
500 DUNE_THROW(Dune::NotImplemented,
501 "Context update for discretization scheme " << MassDiscretizationMethod{}
522 template<std::
size_t i,
class AssembleElementFunc>
525 if (elementSets_.empty())
526 DUNE_THROW(Dune::InvalidStateException,
"Call computeColorsForAssembly before assembling in parallel!");
532 const auto& grid = this->
problem(domainId).gridGeometry().gridView().grid();
533 for (
const auto& elements : elementSets_)
537 const auto element = grid.entity(elements[n]);
538 assembleElement(element);
544 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
545 const Element<freeFlowMomentumIndex>& elementI)
const
548 if (momentumCouplingContext_().empty())
551 bindCouplingContext_(domainI, elementI, momentumCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
554 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
555 const Element<freeFlowMomentumIndex>& elementI,
556 const std::size_t eIdx)
const
558 if (momentumCouplingContext_().empty())
561 fvGeometry.bind(elementI);
570 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[
freeFlowMassIndex]);
572 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
574 else if (eIdx != momentumCouplingContext_()[0].eIdx)
576 momentumCouplingContext_()[0].eIdx = eIdx;
577 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
578 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->
curSol(freeFlowMassIndex));
581 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[
freeFlowMassIndex]);
585 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
586 const Element<freeFlowMassIndex>& elementI)
const
589 if (massAndEnergyCouplingContext_().empty())
590 bindCouplingContext_(domainI, elementI, this->
problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI));
592 bindCouplingContext_(domainI, elementI, massAndEnergyCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
595 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
596 const Element<freeFlowMassIndex>& elementI,
597 const std::size_t eIdx)
const
599 if (massAndEnergyCouplingContext_().empty())
602 auto fvGeometry =
localView(gridGeometry);
603 fvGeometry.bindElement(elementI);
604 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
606 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
608 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
609 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
617 template<std::
size_t i>
618 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
const
620 if (std::get<i>(gridVariables_))
621 return *std::get<i>(gridVariables_);
623 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
630 template<std::
size_t i>
631 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
633 if (std::get<i>(gridVariables_))
634 return *std::get<i>(gridVariables_);
636 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
640 void computeCouplingStencils_()
643 const auto& massGridGeometry = this->
problem(freeFlowMassIndex).gridGeometry();
644 auto momentumFvGeometry =
localView(momentumGridGeometry);
645 auto massFvGeometry =
localView(massGridGeometry);
647 massAndEnergyToMomentumStencils_.clear();
648 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
650 momentumToMassAndEnergyStencils_.clear();
651 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
653 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
655 for (
const auto& element : elements(momentumGridGeometry.gridView()))
657 momentumFvGeometry.bindElement(element);
658 massFvGeometry.bindElement(element);
659 const auto eIdx = momentumFvGeometry.elementIndex();
661 for (
const auto& scv : scvs(momentumFvGeometry))
662 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
664 for (
const auto& scv : scvs(massFvGeometry))
665 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
669 CouplingStencilType emptyStencil_;
670 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
671 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
676 std::vector<MomentumCouplingContext>& momentumCouplingContext_()
const
678 thread_local static std::vector<MomentumCouplingContext> c;
683 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_()
const
685 thread_local static std::vector<MassAndEnergyCouplingContext> c;
690 GridVariablesTuple gridVariables_;
692 const SolutionVector* prevSol_;
695 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
701:
public std::false_type {};
Coloring schemes for shared-memory-parallel assembly.
Element solution classes and factory functions.
free functions for the evaluation of primary variables inside elements.
The available discretization methods in Dumux.
Parallel for loop (multithreading)
PrimaryVariables evalSolution(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const BoxElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:165
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:172
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: couplingmanager_pq1bubble.hh:448
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:251
constexpr FCDiamond fcdiamond
Definition: method.hh:141
constexpr CCTpfa cctpfa
Definition: method.hh:134
constexpr Box box
Definition: method.hh:136
A vector of primary variables.
Definition: common/properties.hh:49
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
Definition: common/properties.hh:100
The secondary variables within a sub-control volume.
Definition: common/properties.hh:105
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:334
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:299
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:321
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:76
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:350
void updateSolution(const SolutionVector &curSol)
Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want t...
Definition: multidomain/couplingmanager.hh:231
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:83
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_pq1bubble.hh:60
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_pq1bubble.hh:63
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_pq1bubble.hh:510
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_pq1bubble.hh:64
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given sub control volume face.
Definition: couplingmanager_pq1bubble.hh:178
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, Dune::index_constant< freeFlowMassIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: couplingmanager_pq1bubble.hh:413
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol, const SolutionVector &prevSol)
use as regular coupling manager in a transient setting
Definition: couplingmanager_pq1bubble.hh:144
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: couplingmanager_pq1bubble.hh:239
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< j > domainJ) const
The coupling stencil of domain I, i.e. which domain J DOFs the given domain I element's residual depe...
Definition: couplingmanager_pq1bubble.hh:377
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMassIndex > domainI, const Element< freeFlowMassIndex > &elementI, Dune::index_constant< freeFlowMomentumIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: couplingmanager_pq1bubble.hh:397
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol)
Methods to be accessed by main.
Definition: couplingmanager_pq1bubble.hh:131
static constexpr auto pressureIdx
Definition: couplingmanager_pq1bubble.hh:123
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume face.
Definition: couplingmanager_pq1bubble.hh:190
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given sub control volume face.
Definition: couplingmanager_pq1bubble.hh:287
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, typename ParentType::SolutionVectorStorage &curSol)
use as binary coupling manager in multi model context
Definition: couplingmanager_pq1bubble.hh:156
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: couplingmanager_pq1bubble.hh:329
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_pq1bubble.hh:354
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_pq1bubble.hh:523
trait that is specialized for coupling manager supporting multithreaded assembly
Definition: multidomain/fvassembler.hh:85
A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
Declares all properties used in Dumux.
The interface of the coupling manager for multi domain problems.