version 3.11-dev
couplingmanager_cvfe.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_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_CVFE_HH
13#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_CVFE_HH
14
15#include <memory>
16#include <tuple>
17#include <vector>
18#include <deque>
19
20#include <dune/common/exceptions.hh>
21#include <dune/common/indices.hh>
22#include <dune/common/float_cmp.hh>
23#include <dune/geometry/referenceelements.hh>
24
28
33
36
39
40#include "typetraits.hh"
41
42namespace Dumux {
43
49template<class Traits>
51: public CouplingManager<Traits>
52{
54public:
55 static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index();
56 static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index();
57
58 // this can be used if the coupling manager is used inside a meta-coupling manager (e.g. multi-binary)
59 // to manager the solution vector storage outside this class
61private:
62 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
63 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
64 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
65 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
66 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
67 template<std::size_t id> using ElementSeed = typename GridView<id>::Grid::template Codim<0>::EntitySeed;
68 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
69 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
70 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
71 template<std::size_t id> using GridVariables = typename Traits::template SubDomain<id>::GridVariables;
72 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
73 template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache;
74 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
75 template<std::size_t id> using VolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::VolumeVariables>;
76
77 using Scalar = typename Traits::Scalar;
78 using SolutionVector = typename Traits::SolutionVector;
79
80 using CouplingStencilType = std::vector<std::size_t>;
81
82 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
83
84 using FluidSystem = typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
85
86 using GlobalPosition = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
87 using VelocityVector = GlobalPosition;
88 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
89
90 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
91
92 struct MomentumCouplingContext
93 {
94 FVElementGeometry<freeFlowMassIndex> fvGeometry;
95 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
96 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
97 std::size_t eIdx;
98 };
99
100 struct MassAndEnergyCouplingContext
101 {
102 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f, const std::size_t i)
103 : fvGeometry(std::move(f))
104 , eIdx(i)
105 {}
106
107 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
108 std::size_t eIdx;
109 };
110
111 using MomentumDiscretizationMethod = typename GridGeometry<freeFlowMomentumIndex>::DiscretizationMethod;
112 using MassDiscretizationMethod = typename GridGeometry<freeFlowMassIndex>::DiscretizationMethod;
113
114 template<std::size_t id> using IpData
115 = Dumux::CVFE::InterpolationPointData<typename GridView<id>::template Codim<0>::Entity::Geometry::LocalCoordinate,
116 typename GridView<id>::template Codim<0>::Entity::Geometry::GlobalCoordinate>;
117
118public:
119
120 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
121
125 // \{
126
128 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
129 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
130 GridVariablesTuple&& gridVariables,
131 const SolutionVector& curSol)
132 {
133 this->momentumCouplingContext_().clear();
134 this->massAndEnergyCouplingContext_().clear();
135
136 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
137 gridVariables_ = gridVariables;
138 this->updateSolution(curSol);
139
140 computeCouplingStencils_();
141 }
142
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)
149 {
150 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
151 prevSol_ = &prevSol;
152 isTransient_ = true;
153 }
154
156 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
157 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
158 GridVariablesTuple&& gridVariables,
160 {
161 this->momentumCouplingContext_().clear();
162 this->massAndEnergyCouplingContext_().clear();
163
164 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
165 gridVariables_ = gridVariables;
166 this->attachSolution(curSol);
167
168 computeCouplingStencils_();
169 }
170
171 // \}
172
176 // \{
177
181 [[deprecated("This method will be removed after release (3.11). Use pressure(..., ipData) instead!")]]
182 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
183 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
184 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
185 const bool considerPreviousTimeStep = false) const
186 {
187 const auto& globalPos = scvf.ipGlobal();
188 const auto& localPos = element.geometry().local(globalPos);
189 return this->pressure(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
190 }
191
195 [[deprecated("This method will be removed after release (3.11). Use pressure(..., ipData) instead!")]]
196 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
197 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
198 const SubControlVolume<freeFlowMomentumIndex>& scv,
199 const bool considerPreviousTimeStep = false) const
200 {
201 return this->pressure(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
202 }
203
207 template <class IpData>
208 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
209 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
210 const IpData& ipData,
211 const bool considerPreviousTimeStep = false) const
212 {
213 assert(!(considerPreviousTimeStep && !this->isTransient_));
214 const auto& gg = this->problem(freeFlowMassIndex).gridGeometry();
215 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
216 : this->curSol(freeFlowMassIndex);
217 const auto elemSol = elementSolution(element, sol, gg);
218 return evalSolutionAtLocalPos(element, element.geometry(), gg, elemSol, ipData.local())[pressureIdx];
219 }
220
224 [[deprecated("This method will be removed after release (3.11). Use density(..., ipData) instead!")]]
225 Scalar density(const Element<freeFlowMomentumIndex>& element,
226 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
227 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
228 const bool considerPreviousTimeStep = false) const
229 {
230 const auto& globalPos = scvf.ipGlobal();
231 const auto& localPos = element.geometry().local(globalPos);
232 return this->density(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
233 }
234
238 [[deprecated("This method will be removed after release (3.11). Use density(..., ipData) instead!")]]
239 Scalar density(const Element<freeFlowMomentumIndex>& element,
240 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
241 const SubControlVolume<freeFlowMomentumIndex>& scv,
242 const bool considerPreviousTimeStep = false) const
243 {
244 return this->density(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
245 }
246
250 template <class IpData>
251 Scalar density(const Element<freeFlowMomentumIndex>& element,
252 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
253 const IpData& ipData,
254 const bool considerPreviousTimeStep = false) const
255 {
256 assert(!(considerPreviousTimeStep && !this->isTransient_));
257 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
258
259 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
260 {
261 const auto eIdx = fvGeometry.elementIndex();
262 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
263
264 const auto& volVars = considerPreviousTimeStep ?
265 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
266 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
267
268 return volVars.density();
269 }
270 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
271 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
272 {
273 // TODO: cache the shape values when Box method is used
274 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
275 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
276 std::vector<ShapeValue> shapeValues;
277 localBasis.evaluateFunction(ipData.local(), shapeValues);
278
279 Scalar rho = 0.0;
280 for (const auto& localDof : localDofs(this->momentumCouplingContext_()[0].fvGeometry))
281 {
282 const auto& volVars = considerPreviousTimeStep ?
283 this->momentumCouplingContext_()[0].prevElemVolVars[localDof.index()]
284 : this->momentumCouplingContext_()[0].curElemVolVars[localDof.index()];
285 rho += volVars.density()*shapeValues[localDof.index()][0];
286 }
287
288 return rho;
289 }
290 else
291 DUNE_THROW(Dune::NotImplemented,
292 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
293 );
294 }
295
299 [[deprecated("This method will be removed after release (3.11). Use effectiveViscosity(..., ipData) instead!")]]
300 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
301 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
302 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
303 const bool considerPreviousTimeStep = false) const
304 {
305 const auto& globalPos = scvf.ipGlobal();
306 const auto& localPos = element.geometry().local(globalPos);
307 return this->effectiveViscosity(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
308 }
309
313 [[deprecated("This method will be removed after release (3.11). Use effectiveViscosity(..., ipData) instead!")]]
314 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
315 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
316 const SubControlVolume<freeFlowMomentumIndex>& scv,
317 const bool considerPreviousTimeStep = false) const
318 {
319 return this->effectiveViscosity(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
320 }
321
325 template <class IpData>
326 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
327 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
328 const IpData& ipData,
329 const bool considerPreviousTimeStep = false) const
330 {
331 assert(!(considerPreviousTimeStep && !this->isTransient_));
332 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
333
334 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
335 {
336 const auto eIdx = fvGeometry.elementIndex();
337 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
338
339 const auto& volVars = considerPreviousTimeStep ?
340 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
341 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
342
343 return volVars.viscosity();
344 }
345 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
346 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
347 {
348 // TODO: cache the shape values
349 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
350 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
351 std::vector<ShapeValue> shapeValues;
352 localBasis.evaluateFunction(ipData.local(), shapeValues);
353
354 Scalar mu = 0.0;
355 for (const auto& localDof : localDofs(this->momentumCouplingContext_()[0].fvGeometry))
356 {
357 const auto& volVars = considerPreviousTimeStep ?
358 this->momentumCouplingContext_()[0].prevElemVolVars[localDof.index()]
359 : this->momentumCouplingContext_()[0].curElemVolVars[localDof.index()];
360 mu += volVars.viscosity()*shapeValues[localDof.index()][0];
361 }
362
363 return mu;
364 }
365 else
366 DUNE_THROW(Dune::NotImplemented,
367 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
368 );
369 }
370
374 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
375 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
376 {
377 // TODO: optimize this function for tpfa where the scvf ip coincides with the dof location
378
379 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(element);
380 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, eIdx);
381
382 const auto& fvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
383 const auto& localBasis = fvGeometry.feLocalBasis();
384
385 std::vector<ShapeValue> shapeValues;
386 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
387 localBasis.evaluateFunction(ipLocal, shapeValues);
388
389 // interpolate velocity at scvf
390 VelocityVector velocity(0.0);
391 for (const auto& localDof : localDofs(fvGeometry))
392 velocity.axpy(shapeValues[localDof.index()][0], this->curSol(freeFlowMomentumIndex)[localDof.dofIndex()]);
393
394 return velocity;
395 }
396
400 VelocityVector elementVelocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry) const
401 {
402 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), fvGeometry.element());
403
404 const auto& momentumFvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
405 const auto& localBasis = momentumFvGeometry.feLocalBasis();
406
407 // interpolate velocity at scvf
408 VelocityVector velocity(0.0);
409 std::vector<ShapeValue> shapeValues;
410 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
411
412 for (const auto& localDof : localDofs(momentumFvGeometry))
413 velocity.axpy(shapeValues[localDof.index()][0], this->curSol(freeFlowMomentumIndex)[localDof.dofIndex()]);
414
415 return velocity;
416 }
417
422 template<std::size_t j>
423 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
424 const Element<freeFlowMomentumIndex>& elementI,
425 const SubControlVolume<freeFlowMomentumIndex>& scvI,
426 Dune::index_constant<j> domainJ) const
427 { return emptyStencil_; }
428
443 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
444 const Element<freeFlowMassIndex>& elementI,
445 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
446 {
447 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
448 return massAndEnergyToMomentumStencils_[eIdx];
449 }
450
459 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
460 const Element<freeFlowMomentumIndex>& elementI,
461 Dune::index_constant<freeFlowMassIndex> domainJ) const
462 {
463 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
464 return momentumToMassAndEnergyStencils_[eIdx];
465 }
466
467 // \}
468
472 // \{
473
493 template<std::size_t i, std::size_t j, class LocalAssemblerI>
494 void updateCouplingContext(Dune::index_constant<i> domainI,
495 const LocalAssemblerI& localAssemblerI,
496 Dune::index_constant<j> domainJ,
497 std::size_t dofIdxGlobalJ,
498 const PrimaryVariables<j>& priVarsJ,
499 int pvIdxJ)
500 {
501 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
502
503 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
504 {
505 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
506 {
507 bindCouplingContext_(domainI, localAssemblerI.element());
508
509 const auto& problem = this->problem(domainJ);
510 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
511 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
512 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
513 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
514
515 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
516 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
517 else
518 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
519 }
520 }
521 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
522 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
523 {
524 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
525 {
526 bindCouplingContext_(domainI, localAssemblerI.element());
527
528 const auto& problem = this->problem(domainJ);
529 const auto& deflectedElement = problem.gridGeometry().element(this->momentumCouplingContext_()[0].eIdx);
530 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
531 const auto& fvGeometry = this->momentumCouplingContext_()[0].fvGeometry;
532
533 // ToDo: Replace once all mass models are also working with local dofs
534 for (const auto& scv : scvs(fvGeometry))
535 {
536 if(scv.dofIndex() == dofIdxGlobalJ)
537 {
538 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
539 this->gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
540 else
541 this->momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
542 }
543 }
544 }
545 }
546 else
547 DUNE_THROW(Dune::NotImplemented,
548 "Context update for discretization scheme " << MassDiscretizationMethod{}
549 );
550 }
551
552 // \}
553
558 {
559 if constexpr (MomentumDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
560 {
561 // use coloring of the mass discretization for both domains
562 // the diamond coloring is a subset (minimum amount of colors) of cctpfa/box coloring
563 elementSets_ = computeColoring(this->problem(freeFlowMassIndex).gridGeometry()).sets;
564 }
565 else
566 {
567 // use coloring of the momentum discretization for both domains
568 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
569 }
570 }
571
578 template<std::size_t i, class AssembleElementFunc>
579 void assembleMultithreaded(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
580 {
581 if (elementSets_.empty())
582 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
583
584 // make this element loop run in parallel
585 // for this we have to color the elements so that we don't get
586 // race conditions when writing into the global matrix
587 // each color can be assembled using multiple threads
588 const auto& grid = this->problem(freeFlowMassIndex).gridGeometry().gridView().grid();
589 for (const auto& elements : elementSets_)
590 {
591 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
592 {
593 const auto element = grid.entity(elements[eIdx]);
594 assembleElement(element);
595 });
596 }
597 }
598
599private:
600 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
601 const Element<freeFlowMomentumIndex>& elementI) const
602 {
603 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
604 if (momentumCouplingContext_().empty())
605 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI));
606 else
607 bindCouplingContext_(domainI, elementI, momentumCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
608 }
609
610 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
611 const Element<freeFlowMomentumIndex>& elementI,
612 const std::size_t eIdx) const
613 {
614 if (momentumCouplingContext_().empty())
615 {
616 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
617 fvGeometry.bind(elementI);
618
619 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
620 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
621
622 auto prevElemVolVars = isTransient_ ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
623 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
624
625 if (isTransient_)
626 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[freeFlowMassIndex]);
627
628 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
629 }
630 else if (eIdx != momentumCouplingContext_()[0].eIdx)
631 {
632 momentumCouplingContext_()[0].eIdx = eIdx;
633 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
634 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->curSol(freeFlowMassIndex));
635
636 if (isTransient_)
637 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[freeFlowMassIndex]);
638 }
639 }
640
641 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
642 const Element<freeFlowMassIndex>& elementI) const
643 {
644 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
645 if (massAndEnergyCouplingContext_().empty())
646 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI));
647 else
648 bindCouplingContext_(domainI, elementI, massAndEnergyCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
649 }
650
651 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
652 const Element<freeFlowMassIndex>& elementI,
653 const std::size_t eIdx) const
654 {
655 if (massAndEnergyCouplingContext_().empty())
656 {
657 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
658 auto fvGeometry = localView(gridGeometry);
659 fvGeometry.bindElement(elementI);
660 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
661 }
662 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
663 {
664 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
665 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
666 }
667 }
668
673 template<std::size_t i>
674 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
675 {
676 if (std::get<i>(gridVariables_))
677 return *std::get<i>(gridVariables_);
678 else
679 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
680 }
681
686 template<std::size_t i>
687 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
688 {
689 if (std::get<i>(gridVariables_))
690 return *std::get<i>(gridVariables_);
691 else
692 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
693 }
694
695
696 void computeCouplingStencils_()
697 {
698 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
699 const auto& massGridGeometry = this->problem(freeFlowMassIndex).gridGeometry();
700 auto momentumFvGeometry = localView(momentumGridGeometry);
701 auto massFvGeometry = localView(massGridGeometry);
702
703 massAndEnergyToMomentumStencils_.clear();
704 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
705
706 momentumToMassAndEnergyStencils_.clear();
707 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
708
709 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
710
711 for (const auto& element : elements(momentumGridGeometry.gridView()))
712 {
713 momentumFvGeometry.bindElement(element);
714 massFvGeometry.bindElement(element);
715 const auto eIdx = momentumFvGeometry.elementIndex();
716
717 for (const auto& localDof : localDofs(momentumFvGeometry))
718 massAndEnergyToMomentumStencils_[eIdx].push_back(localDof.dofIndex());
719
720 // ToDo: Replace once all mass models are also working with local dofs
721 for (const auto& scv : scvs(massFvGeometry))
722 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
723 }
724 }
725
726 CouplingStencilType emptyStencil_;
727 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
728 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
729
730 std::vector<MomentumCouplingContext>& momentumCouplingContext_() const
731 { return momentumCouplingContextImpl_; }
732
733 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_() const
734 { return massAndEnergyCouplingContextImpl_; }
735
736 mutable std::vector<MassAndEnergyCouplingContext> massAndEnergyCouplingContextImpl_;
737 mutable std::vector<MomentumCouplingContext> momentumCouplingContextImpl_;
738
740 GridVariablesTuple gridVariables_;
741
742 const SolutionVector* prevSol_;
743 bool isTransient_;
744
745 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
746};
747
748namespace Detail {
749
750// declaration (specialize for different discretization types)
751template<class Traits, class DiscretizationMethod = typename Detail::MomentumDiscretizationMethod<Traits>::type>
753
754// multi-threading is not supported because we have only one coupling context instance and a mutable cache
755template<class Traits, class D>
756struct CouplingManagerSupportsMultithreadedAssemblySelector<Traits, DiscretizationMethods::CVFE<D>>
757{ using type = std::false_type; };
758
759} // end namespace Detail
760
762template<class T>
765{};
766
767} // end namespace Dumux
768
769#endif
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:25
const LocalPosition & local() const
The local position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:35
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_cvfe.hh:52
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const IpData &ipData, const bool considerPreviousTimeStep=false) const
Returns the pressure at a given interpolation point.
Definition: couplingmanager_cvfe.hh:208
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: couplingmanager_cvfe.hh:374
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_cvfe.hh:128
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const typename ParentType::SolutionVectorStorage &curSol)
use as binary coupling manager in multi model context
Definition: couplingmanager_cvfe.hh:156
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the effective viscosity at a given sub control volume.
Definition: couplingmanager_cvfe.hh:314
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_cvfe.hh:55
static constexpr auto pressureIdx
Definition: couplingmanager_cvfe.hh:120
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_cvfe.hh:443
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_cvfe.hh:459
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_cvfe.hh:225
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const IpData &ipData, const bool considerPreviousTimeStep=false) const
Returns the density at a given position.
Definition: couplingmanager_cvfe.hh:251
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the effective viscosity at a given sub control volume face.
Definition: couplingmanager_cvfe.hh:300
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_cvfe.hh:400
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_cvfe.hh:239
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const IpData &ipData, const bool considerPreviousTimeStep=false) const
Returns the effective viscosity at a given position.
Definition: couplingmanager_cvfe.hh:326
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the pressure at a given sub control volume.
Definition: couplingmanager_cvfe.hh:196
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_cvfe.hh:579
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the pressure at a given sub control volume face.
Definition: couplingmanager_cvfe.hh:182
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_cvfe.hh:144
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_cvfe.hh:557
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_cvfe.hh:56
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_cvfe.hh:423
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
void attachSolution(const SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:311
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:276
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:298
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:53
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:327
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:208
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:60
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Type traits.
Classes representing interpolation point data for control-volume finite element schemes.
Element solution classes and factory functions.
free functions for the evaluation of primary variables inside elements.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
PrimaryVariables evalSolutionAtLocalPos(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::LocalCoordinate &localPos, bool ignoreState=false)
Interpolates a given cvfe element solution at a given local position. Uses the finite element cache o...
Definition: evalsolution.hh:173
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_cvfe.hh:494
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Class representing dofs on elements for control-volume finite element schemes.
The available discretization methods in Dumux.
The interface of the coupling manager for multi domain problems.
Some useful type traits.
A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
constexpr FCDiamond fcdiamond
Definition: method.hh:152
constexpr CCTpfa cctpfa
Definition: method.hh:145
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:239
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50
Parallel for loop (multithreading)
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78