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#include <iostream>
20
21#include <dune/common/exceptions.hh>
22#include <dune/common/indices.hh>
23#include <dune/common/float_cmp.hh>
24#include <dune/geometry/referenceelements.hh>
25
29
34
37
40
41#include "typetraits.hh"
42
43namespace Dumux {
44
50template<class Traits>
52: public CouplingManager<Traits>
53{
55public:
56 static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index();
57 static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index();
58
59 // this can be used if the coupling manager is used inside a meta-coupling manager (e.g. multi-binary)
60 // to manager the solution vector storage outside this class
62private:
63 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
64 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
65 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
66 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
67 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
68 template<std::size_t id> using ElementSeed = typename GridView<id>::Grid::template Codim<0>::EntitySeed;
69 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
70 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
71 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
72 template<std::size_t id> using GridVariables = typename Traits::template SubDomain<id>::GridVariables;
73 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
74 template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache;
75 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
76 template<std::size_t id> using VolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::VolumeVariables>;
77
78 using Scalar = typename Traits::Scalar;
79 using SolutionVector = typename Traits::SolutionVector;
80
81 using CouplingStencilType = std::vector<std::size_t>;
82
83 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
84
85 using FluidSystem = typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
86
87 using GlobalPosition = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
88 using VelocityVector = GlobalPosition;
89 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
90
91 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
92
93 template<class ElementSolution>
94 struct MomentumCouplingContextNoCaching
95 {
96 MomentumCouplingContextNoCaching(ElementSolution&& elemSol)
97 : elemSol_(std::move(elemSol)) {}
98
99 template<class GridVolVars, class FvElementGeometry, class SubControlVolume>
100 auto vars(const GridVolVars& gridVolVars, const FvElementGeometry& fvGeometry, const SubControlVolume& scv) const
101 {
102 const auto& problem = gridVolVars.problem();
103 VolumeVariables<freeFlowMassIndex> volVars;
104 volVars.update(elemSol_, problem, fvGeometry.element(), scv);
105 return volVars;
106 }
107
108 ElementSolution elemSol_;
109 };
110
111 struct MomentumCouplingContextGlobalCaching
112 {
113 template<class GridVolVars, class FvElementGeometry, class SubControlVolume>
114 const auto& vars(const GridVolVars& gridVolVars, const FvElementGeometry& fvGeometry, const SubControlVolume& scv) const
115 {
116 return gridVolVars.volVars(scv);
117 }
118 };
119
120 using MomentumDiscretizationMethod = typename GridGeometry<freeFlowMomentumIndex>::DiscretizationMethod;
121 using MassDiscretizationMethod = typename GridGeometry<freeFlowMassIndex>::DiscretizationMethod;
122
123 template<std::size_t id> using IpData
124 = Dumux::CVFE::InterpolationPointData<typename GridView<id>::template Codim<0>::Entity::Geometry::LocalCoordinate,
125 typename GridView<id>::template Codim<0>::Entity::Geometry::GlobalCoordinate>;
126
127public:
128
129 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
130
134 // \{
135
137 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
138 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
139 GridVariablesTuple&& gridVariables,
140 const SolutionVector& curSol)
141 {
142 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
143 gridVariables_ = gridVariables;
144 this->updateSolution(curSol);
145
146 computeCouplingStencils_();
147 }
148
150 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
151 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
152 GridVariablesTuple&& gridVariables,
153 const SolutionVector& curSol,
154 const SolutionVector& prevSol)
155 {
156 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
157 prevSol_ = &prevSol;
158 isTransient_ = true;
159 }
160
162 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
163 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
164 GridVariablesTuple&& gridVariables,
166 {
167 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
168 gridVariables_ = gridVariables;
169 this->attachSolution(curSol);
170
171 computeCouplingStencils_();
172 }
173
174 // \}
175
179 // \{
180
184 [[deprecated("This method will be removed after release (3.11). Use pressure(..., ipData) instead!")]]
185 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
186 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
187 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
188 const bool considerPreviousTimeStep = false) const
189 {
190 const auto& globalPos = scvf.ipGlobal();
191 const auto& localPos = element.geometry().local(globalPos);
192 return this->pressure(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
193 }
194
198 [[deprecated("This method will be removed after release (3.11). Use pressure(..., ipData) instead!")]]
199 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
200 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
201 const SubControlVolume<freeFlowMomentumIndex>& scv,
202 const bool considerPreviousTimeStep = false) const
203 {
204 return this->pressure(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
205 }
206
210 template <class IpData>
211 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
212 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
213 const IpData& ipData,
214 const bool considerPreviousTimeStep = false) const
215 {
216 assert(!(considerPreviousTimeStep && !this->isTransient_));
217 const auto& gg = this->problem(freeFlowMassIndex).gridGeometry();
218 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
219 : this->curSol(freeFlowMassIndex);
220 const auto elemSol = elementSolution(element, sol, gg);
221 return evalSolutionAtLocalPos(element, element.geometry(), gg, elemSol, ipData.local())[pressureIdx];
222 }
223
227 [[deprecated("This method will be removed after release (3.11). Use density(..., ipData) instead!")]]
228 Scalar density(const Element<freeFlowMomentumIndex>& element,
229 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
230 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
231 const bool considerPreviousTimeStep = false) const
232 {
233 const auto& globalPos = scvf.ipGlobal();
234 const auto& localPos = element.geometry().local(globalPos);
235 return this->density(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
236 }
237
241 [[deprecated("This method will be removed after release (3.11). Use density(..., ipData) instead!")]]
242 Scalar density(const Element<freeFlowMomentumIndex>& element,
243 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
244 const SubControlVolume<freeFlowMomentumIndex>& scv,
245 const bool considerPreviousTimeStep = false) const
246 {
247 return this->density(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
248 }
249
253 template <class IpData>
254 Scalar density(const Element<freeFlowMomentumIndex>& element,
255 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
256 const IpData& ipData,
257 const bool considerPreviousTimeStep = false) const
258 {
259 assert(!(considerPreviousTimeStep && !this->isTransient_));
260
261 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
262 : this->curSol(freeFlowMassIndex);
263
264 auto massFvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
265 massFvGeometry.bind(element);
266 const auto context = makeMomentumCouplingContext_(massFvGeometry, sol);
267 const auto& gridVolVars = considerPreviousTimeStep ? gridVars_(freeFlowMassIndex).curGridVolVars()
268 : gridVars_(freeFlowMassIndex).prevGridVolVars();
269
270 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
271 {
272 const auto eIdx = fvGeometry.elementIndex();
273 const auto& scv = massFvGeometry.scv(eIdx);
274
275 const auto& volVars = context.vars(gridVolVars, massFvGeometry, scv);
276
277 return volVars.density();
278 }
279 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
280 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
281 {
282 // TODO: cache the shape values
283 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
284 const auto& localBasis = massFvGeometry.feLocalBasis();
285 std::vector<ShapeValue> shapeValues;
286 localBasis.evaluateFunction(ipData.local(), shapeValues);
287
288 Scalar rho = 0.0;
289 for (const auto& scv : scvs(massFvGeometry))
290 {
291 const auto& volVars = context.vars(gridVolVars, massFvGeometry, scv);
292 rho += volVars.density()*shapeValues[scv.localDofIndex()][0];
293 }
294
295 return rho;
296 }
297 else
298 DUNE_THROW(Dune::NotImplemented,
299 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
300 );
301 }
302
306 [[deprecated("This method will be removed after release (3.11). Use effectiveViscosity(..., ipData) instead!")]]
307 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
308 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
309 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
310 const bool considerPreviousTimeStep = false) const
311 {
312 const auto& globalPos = scvf.ipGlobal();
313 const auto& localPos = element.geometry().local(globalPos);
314 return this->effectiveViscosity(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
315 }
316
320 [[deprecated("This method will be removed after release (3.11). Use effectiveViscosity(..., ipData) instead!")]]
321 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
322 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
323 const SubControlVolume<freeFlowMomentumIndex>& scv,
324 const bool considerPreviousTimeStep = false) const
325 {
326 return this->effectiveViscosity(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
327 }
328
332 template <class IpData>
333 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
334 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
335 const IpData& ipData,
336 const bool considerPreviousTimeStep = false) const
337 {
338 assert(!(considerPreviousTimeStep && !this->isTransient_));
339
340 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
341 : this->curSol(freeFlowMassIndex);
342
343 auto massFvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
344 massFvGeometry.bind(element);
345 const auto context = makeMomentumCouplingContext_(massFvGeometry, sol);
346 const auto& gridVolVars = considerPreviousTimeStep ? gridVars_(freeFlowMassIndex).curGridVolVars()
347 : gridVars_(freeFlowMassIndex).prevGridVolVars();
348
349 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
350 {
351 const auto eIdx = fvGeometry.elementIndex();
352 const auto& scv = massFvGeometry.scv(eIdx);
353
354 const auto& volVars = context.vars(gridVolVars, massFvGeometry, scv);
355
356 return volVars.viscosity();
357 }
358 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
359 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
360 {
361 // TODO: cache the shape values
362 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
363 const auto& localBasis = massFvGeometry.feLocalBasis();
364 std::vector<ShapeValue> shapeValues;
365 localBasis.evaluateFunction(ipData.local(), shapeValues);
366
367 Scalar mu = 0.0;
368 for (const auto& scv : scvs(massFvGeometry))
369 {
370 const auto& volVars = context.vars(gridVolVars, massFvGeometry, scv);
371 mu += volVars.viscosity()*shapeValues[scv.localDofIndex()][0];
372 }
373
374 return mu;
375 }
376 else
377 DUNE_THROW(Dune::NotImplemented,
378 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
379 );
380 }
381
385 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
386 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
387 {
388 // TODO: optimize this function for tpfa where the scvf ip coincides with the dof location
389 auto fvGeometry = localView(this->problem(freeFlowMomentumIndex).gridGeometry());
390 fvGeometry.bindElement(element);
391
392 const auto& localBasis = fvGeometry.feLocalBasis();
393
394 std::vector<ShapeValue> shapeValues;
395 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
396 localBasis.evaluateFunction(ipLocal, shapeValues);
397
398 // interpolate velocity at scvf
399 VelocityVector velocity(0.0);
400 for (const auto& localDof : localDofs(fvGeometry))
401 velocity.axpy(shapeValues[localDof.index()][0], this->curSol(freeFlowMomentumIndex)[localDof.dofIndex()]);
402
403 return velocity;
404 }
405
409 VelocityVector elementVelocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry) const
410 {
411 auto momentumFvGeometry = localView(this->problem(freeFlowMomentumIndex).gridGeometry());
412 momentumFvGeometry.bindElement(fvGeometry.element());
413
414 const auto& localBasis = momentumFvGeometry.feLocalBasis();
415
416 // interpolate velocity at scvf
417 VelocityVector velocity(0.0);
418 std::vector<ShapeValue> shapeValues;
419 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
420
421 for (const auto& localDof : localDofs(momentumFvGeometry))
422 velocity.axpy(shapeValues[localDof.index()][0], this->curSol(freeFlowMomentumIndex)[localDof.dofIndex()]);
423
424 return velocity;
425 }
426
430 template <class IpData>
431 VelocityVector velocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry,
432 const IpData& ipData,
433 const bool considerPreviousTimeStep = false) const
434 {
435 assert(!(considerPreviousTimeStep && !this->isTransient_));
436
437 const auto& element = fvGeometry.element();
438 const auto& gg = this->problem(freeFlowMomentumIndex).gridGeometry();
439 auto momentumFvGeometry = localView(gg);
440 momentumFvGeometry.bindElement(fvGeometry.element());
441
442 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMomentumIndex]
444
445 const auto elemSol = elementSolution(element, sol, gg);
446 return evalSolutionAtLocalPos(element, element.geometry(), gg, elemSol, ipData.local());
447 }
448
453 template<std::size_t j>
454 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
455 const Element<freeFlowMomentumIndex>& elementI,
456 const SubControlVolume<freeFlowMomentumIndex>& scvI,
457 Dune::index_constant<j> domainJ) const
458 { return emptyStencil_; }
459
474 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
475 const Element<freeFlowMassIndex>& elementI,
476 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
477 {
478 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
479 return massAndEnergyToMomentumStencils_[eIdx];
480 }
481
490 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
491 const Element<freeFlowMomentumIndex>& elementI,
492 Dune::index_constant<freeFlowMassIndex> domainJ) const
493 {
494 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
495 return momentumToMassAndEnergyStencils_[eIdx];
496 }
497
498 // \}
499
503 // \{
504
524 template<std::size_t i, std::size_t j, class LocalAssemblerI>
525 void updateCouplingContext(Dune::index_constant<i> domainI,
526 const LocalAssemblerI& localAssemblerI,
527 Dune::index_constant<j> domainJ,
528 std::size_t dofIdxGlobalJ,
529 const PrimaryVariables<j>& priVarsJ,
530 int pvIdxJ)
531 {
532 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
533
534 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
535 {
536 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
537 {
538 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
539 {
540 const auto& problem = this->problem(domainJ);
541 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
542 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
543 auto fvGeometry = localView(problem.gridGeometry());
544 fvGeometry.bind(deflectedElement);
545 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
546
547 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
548 }
549 }
550 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
551 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
552 {
553 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
554 {
555 const auto& problem = this->problem(domainJ);
556 const auto deflectedElementIdx = problem.gridGeometry().elementMapper().index(localAssemblerI.element());
557 const auto& deflectedElement = problem.gridGeometry().element(deflectedElementIdx);
558 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
559 auto fvGeometry = localView(problem.gridGeometry());
560 fvGeometry.bind(deflectedElement);
561
562 // ToDo: Replace once all mass models are also working with local dofs
563 for (const auto& scv : scvs(fvGeometry))
564 {
565 if (scv.dofIndex() == dofIdxGlobalJ)
566 this->gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
567 }
568 }
569 }
570 else
571 DUNE_THROW(Dune::NotImplemented,
572 "Context update for discretization scheme " << MassDiscretizationMethod{}
573 );
574 }
575 }
576
577 // \}
578
583 {
584 if constexpr (MomentumDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
585 {
586 // use coloring of the mass discretization for both domains
587 // the diamond coloring is a subset (minimum amount of colors) of cctpfa/box coloring
588 elementSets_ = computeColoring(this->problem(freeFlowMassIndex).gridGeometry()).sets;
589 }
590 else
591 {
592 // use coloring of the momentum discretization for both domains
593 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
594 }
595 }
596
603 template<std::size_t i, class AssembleElementFunc>
604 void assembleMultithreaded(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
605 {
606 if (elementSets_.empty())
607 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
608
609 // make this element loop run in parallel
610 // for this we have to color the elements so that we don't get
611 // race conditions when writing into the global matrix
612 // each color can be assembled using multiple threads
613 const auto& grid = this->problem(freeFlowMassIndex).gridGeometry().gridView().grid();
614 for (const auto& elements : elementSets_)
615 {
616 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
617 {
618 const auto element = grid.entity(elements[eIdx]);
619 assembleElement(element);
620 });
621 }
622 }
623
624private:
625 template<class SolutionVector>
626 auto makeMomentumCouplingContext_(const FVElementGeometry<freeFlowMassIndex>& fvGeometry,
627 const SolutionVector& sol) const
628 {
629 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
630 return MomentumCouplingContextGlobalCaching{};
631 else
632 return MomentumCouplingContextNoCaching{elementSolution(fvGeometry.element(), sol, fvGeometry.gridGeometry())};
633 }
634
639 template<std::size_t i>
640 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
641 {
642 if (std::get<i>(gridVariables_))
643 return *std::get<i>(gridVariables_);
644 else
645 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
646 }
647
652 template<std::size_t i>
653 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
654 {
655 if (std::get<i>(gridVariables_))
656 return *std::get<i>(gridVariables_);
657 else
658 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
659 }
660
661 void computeCouplingStencils_()
662 {
663 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
664 const auto& massGridGeometry = this->problem(freeFlowMassIndex).gridGeometry();
665 auto momentumFvGeometry = localView(momentumGridGeometry);
666 auto massFvGeometry = localView(massGridGeometry);
667
668 massAndEnergyToMomentumStencils_.clear();
669 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
670
671 momentumToMassAndEnergyStencils_.clear();
672 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
673
674 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
675 bool hasCubeElements = false;
676
677 for (const auto& element : elements(momentumGridGeometry.gridView()))
678 {
679 if(element.type().isCube())
680 hasCubeElements = true;
681
682 momentumFvGeometry.bindElement(element);
683 massFvGeometry.bindElement(element);
684 const auto eIdx = momentumFvGeometry.elementIndex();
685
686 for (const auto& localDof : localDofs(momentumFvGeometry))
687 massAndEnergyToMomentumStencils_[eIdx].push_back(localDof.dofIndex());
688
689 // ToDo: Replace once all mass models are also working with local dofs
690 for (const auto& scv : scvs(massFvGeometry))
691 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
692 }
693
694 // Print warning for pq1bubble scheme on cube elements if not using the hybrid variant
695 if constexpr ((MomentumDiscretizationMethod{} == DiscretizationMethods::pq1bubble))
696 {
697 if(hasCubeElements && !GridGeometry<freeFlowMomentumIndex>::enableHybridCVFE)
698 {
699 std::cerr << "Warning: Coupled Navier-Stokes problem on cube elements uses non-hybrid pq1bubble. "
700 << "The hybrid variant is recommended because it implements two bubble functions for stability reasons."
701 << std::endl;
702 }
703 }
704 }
705
706 CouplingStencilType emptyStencil_;
707 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
708 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
709
711 GridVariablesTuple gridVariables_;
712
713 const SolutionVector* prevSol_;
714 bool isTransient_;
715
716 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
717};
718
719namespace Detail {
720
721// declaration (specialize for different discretization types)
722template<class Traits, class DiscretizationMethod = typename Detail::MomentumDiscretizationMethod<Traits>::type>
724
725// multi-threading is not supported because we have only one coupling context instance and a mutable cache
726template<class Traits, class D>
727struct CouplingManagerSupportsMultithreadedAssemblySelector<Traits, DiscretizationMethods::CVFE<D>>
728{ using type = std::true_type; };
729
730} // end namespace Detail
731
733template<class T>
736{};
737
738} // end namespace Dumux
739
740#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:40
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_cvfe.hh:53
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:211
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:385
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:137
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:162
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:321
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_cvfe.hh:56
static constexpr auto pressureIdx
Definition: couplingmanager_cvfe.hh:129
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:474
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:490
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:228
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:254
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:307
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_cvfe.hh:409
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:242
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:333
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:199
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_cvfe.hh:604
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:185
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:150
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_cvfe.hh:582
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_cvfe.hh:57
VelocityVector velocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry, const IpData &ipData, const bool considerPreviousTimeStep=false) const
Returns the velocity at an interpolation point.
Definition: couplingmanager_cvfe.hh:431
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:454
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:525
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:162
constexpr CCTpfa cctpfa
Definition: method.hh:154
constexpr Box box
Definition: method.hh:156
constexpr PQ1Bubble pq1bubble
Definition: method.hh:158
Definition: adapt.hh:17
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:241
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