24#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_DIAMOND_HH
25#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_DIAMOND_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();
212 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
213 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
214 std::vector<ShapeValue> shapeValues;
215 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
216 localBasis.evaluateFunction(ipLocal, shapeValues);
219 for (
const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
221 const auto& volVars = considerPreviousTimeStep ?
222 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
223 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
224 rho += volVars.density()*shapeValues[scv.indexInElement()][0];
230 DUNE_THROW(Dune::NotImplemented,
231 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
238 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
239 const SubControlVolume<freeFlowMomentumIndex>& scv,
240 const bool considerPreviousTimeStep =
false)
const
242 assert(!(considerPreviousTimeStep && !this->isTransient_));
243 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
247 const auto eIdx = scv.elementIndex();
248 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
250 const auto& volVars = considerPreviousTimeStep ?
251 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
252 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
254 return volVars.density();
259 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
260 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
261 std::vector<ShapeValue> shapeValues;
262 const auto ipLocal = element.geometry().local(scv.dofPosition());
263 localBasis.evaluateFunction(ipLocal, shapeValues);
266 for (
const auto& scvI : scvs(this->momentumCouplingContext_()[0].fvGeometry))
268 const auto& volVars = considerPreviousTimeStep ?
269 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
270 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
271 rho += volVars.density()*shapeValues[scvI.indexInElement()][0];
276 DUNE_THROW(Dune::NotImplemented,
277 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
285 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
286 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
288 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
292 const auto eIdx = fvGeometry.elementIndex();
293 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
294 const auto& volVars = this->momentumCouplingContext_()[0].curElemVolVars[scv];
295 return volVars.viscosity();
300 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
301 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
302 std::vector<ShapeValue> shapeValues;
303 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
304 localBasis.evaluateFunction(ipLocal, shapeValues);
307 for (
const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
309 const auto& volVars = this->momentumCouplingContext_()[0].curElemVolVars[scv];
310 mu += volVars.viscosity()*shapeValues[scv.indexInElement()][0];
316 DUNE_THROW(Dune::NotImplemented,
317 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
326 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
331 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, eIdx);
333 const auto& fvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
334 const auto& localBasis = fvGeometry.feLocalBasis();
335 std::vector<ShapeValue> shapeValues;
336 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
337 localBasis.evaluateFunction(ipLocal, shapeValues);
340 VelocityVector velocity(0.0);
341 for (
const auto& scv : scvs(fvGeometry))
342 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(
freeFlowMomentumIndex)[scv.dofIndex()]);
350 VelocityVector
elementVelocity(
const FVElementGeometry<freeFlowMassIndex>& fvGeometry)
const
352 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), fvGeometry.element());
354 const auto& momentumFvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
355 const auto& localBasis = momentumFvGeometry.feLocalBasis();
358 VelocityVector velocity(0.0);
359 std::vector<ShapeValue> shapeValues;
360 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
362 for (
const auto& scv : scvs(momentumFvGeometry))
363 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(
freeFlowMomentumIndex)[scv.dofIndex()]);
372 template<std::
size_t j>
374 const Element<freeFlowMomentumIndex>& elementI,
375 const SubControlVolume<freeFlowMomentumIndex>& scvI,
376 Dune::index_constant<j> domainJ)
const
377 {
return emptyStencil_; }
394 const Element<freeFlowMassIndex>& elementI,
395 Dune::index_constant<freeFlowMomentumIndex> domainJ)
const
398 return massAndEnergyToMomentumStencils_[eIdx];
410 const Element<freeFlowMomentumIndex>& elementI,
411 Dune::index_constant<freeFlowMassIndex> domainJ)
const
414 return momentumToMassAndEnergyStencils_[eIdx];
443 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
445 const LocalAssemblerI& localAssemblerI,
446 Dune::index_constant<j> domainJ,
447 std::size_t dofIdxGlobalJ,
448 const PrimaryVariables<j>& priVarsJ,
451 this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
457 bindCouplingContext_(domainI, localAssemblerI.element());
460 const auto& deflectedElement =
problem.gridGeometry().element(dofIdxGlobalJ);
462 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
463 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
465 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
466 gridVars_(
freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol),
problem, deflectedElement, scv);
468 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol),
problem, deflectedElement, scv);
475 bindCouplingContext_(domainI, localAssemblerI.element());
478 const auto& deflectedElement =
problem.gridGeometry().element(this->momentumCouplingContext_()[0].eIdx);
480 const auto& fvGeometry = this->momentumCouplingContext_()[0].fvGeometry;
482 for (
const auto& scv : scvs(fvGeometry))
484 if(scv.dofIndex() == dofIdxGlobalJ)
486 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
487 this->gridVars_(
freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol),
problem, deflectedElement, scv);
489 this->momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol),
problem, deflectedElement, scv);
495 DUNE_THROW(Dune::NotImplemented,
496 "Context update for discretization scheme " << MassDiscretizationMethod{}
518 template<std::
size_t i,
class AssembleElementFunc>
521 if (elementSets_.empty())
522 DUNE_THROW(Dune::InvalidStateException,
"Call computeColorsForAssembly before assembling in parallel!");
529 for (
const auto& elements : elementSets_)
533 const auto element = grid.entity(elements[eIdx]);
534 assembleElement(element);
540 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
541 const Element<freeFlowMomentumIndex>& elementI)
const
544 if (momentumCouplingContext_().empty())
547 bindCouplingContext_(domainI, elementI, momentumCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
550 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
551 const Element<freeFlowMomentumIndex>& elementI,
552 const std::size_t eIdx)
const
554 if (momentumCouplingContext_().empty())
557 fvGeometry.bind(elementI);
566 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[
freeFlowMassIndex]);
568 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
570 else if (eIdx != momentumCouplingContext_()[0].eIdx)
572 momentumCouplingContext_()[0].eIdx = eIdx;
573 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
574 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->
curSol(freeFlowMassIndex));
577 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[
freeFlowMassIndex]);
581 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
582 const Element<freeFlowMassIndex>& elementI)
const
585 if (massAndEnergyCouplingContext_().empty())
586 bindCouplingContext_(domainI, elementI, this->
problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI));
588 bindCouplingContext_(domainI, elementI, massAndEnergyCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
591 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
592 const Element<freeFlowMassIndex>& elementI,
593 const std::size_t eIdx)
const
595 if (massAndEnergyCouplingContext_().empty())
598 auto fvGeometry =
localView(gridGeometry);
599 fvGeometry.bindElement(elementI);
600 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
602 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
604 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
605 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
613 template<std::
size_t i>
614 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
const
616 if (std::get<i>(gridVariables_))
617 return *std::get<i>(gridVariables_);
619 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
626 template<std::
size_t i>
627 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
629 if (std::get<i>(gridVariables_))
630 return *std::get<i>(gridVariables_);
632 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
636 void computeCouplingStencils_()
639 const auto& massGridGeometry = this->
problem(freeFlowMassIndex).gridGeometry();
640 auto momentumFvGeometry =
localView(momentumGridGeometry);
641 auto massFvGeometry =
localView(massGridGeometry);
643 massAndEnergyToMomentumStencils_.clear();
644 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
646 momentumToMassAndEnergyStencils_.clear();
647 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
649 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
651 for (
const auto& element : elements(momentumGridGeometry.gridView()))
653 momentumFvGeometry.bindElement(element);
654 massFvGeometry.bindElement(element);
655 const auto eIdx = momentumFvGeometry.elementIndex();
657 for (
const auto& scv : scvs(momentumFvGeometry))
658 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
660 for (
const auto& scv : scvs(massFvGeometry))
661 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
665 CouplingStencilType emptyStencil_;
666 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
667 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
672 std::vector<MomentumCouplingContext>& momentumCouplingContext_()
const
674 thread_local static std::vector<MomentumCouplingContext> c;
679 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_()
const
681 thread_local static std::vector<MassAndEnergyCouplingContext> c;
686 GridVariablesTuple gridVariables_;
688 const SolutionVector* prevSol_;
691 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
697:
public std::true_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_diamond.hh:444
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 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_diamond.hh:60
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_diamond.hh:64
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_diamond.hh:350
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_diamond.hh:284
Scalar density(const Element< freeFlowMomentumIndex > &element, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: couplingmanager_diamond.hh:238
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_diamond.hh:156
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_diamond.hh:190
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_diamond.hh:505
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_diamond.hh:409
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: couplingmanager_diamond.hh:325
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_diamond.hh:373
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_diamond.hh:144
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_diamond.hh:63
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_diamond.hh:519
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_diamond.hh:178
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_diamond.hh:131
static constexpr auto pressureIdx
Definition: couplingmanager_diamond.hh:123
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_diamond.hh:393
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.