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
27#include <dumux/common/concepts/variables_.hh>
30
35
38
41
42#include "typetraits.hh"
43
44namespace Dumux {
45
51template<class Traits>
53: public CouplingManager<Traits>
54{
56public:
57 static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index();
58 static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index();
59
60 // this can be used if the coupling manager is used inside a meta-coupling manager (e.g. multi-binary)
61 // to manager the solution vector storage outside this class
63private:
64 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
65 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
66 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
67 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
68 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
69 template<std::size_t id> using ElementSeed = typename GridView<id>::Grid::template Codim<0>::EntitySeed;
70 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
71 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
72 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
73 template<std::size_t id> using GridVariables = typename Traits::template SubDomain<id>::GridVariables;
74 template<std::size_t id> using GridVariablesCache = Concept::GridVariablesCache_t<GridVariables<id>>;
75 template<std::size_t id> using ElementVariables = typename GridVariablesCache<id>::LocalView;
76 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
77 template<std::size_t id> using Variables = Concept::Variables_t<GridVariables<id>>;
78
79 using Scalar = typename Traits::Scalar;
80 using SolutionVector = typename Traits::SolutionVector;
81
82 using CouplingStencilType = std::vector<std::size_t>;
83
84 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
85
86 using FluidSystem = typename Variables<freeFlowMassIndex>::FluidSystem;
87
88 using GlobalPosition = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
89 using VelocityVector = GlobalPosition;
90 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
91
92 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
93
94 template<class ElementSolution>
95 struct MomentumCouplingContextNoCaching
96 {
97 MomentumCouplingContextNoCaching(ElementSolution&& elemSol)
98 : elemSol_(std::move(elemSol)) {}
99
100 template<class GridVarsCache, class FvElementGeometry, class SubControlVolume>
101 auto vars(const GridVarsCache& gridVarsCache, const FvElementGeometry& fvGeometry, const SubControlVolume& scv) const
102 {
103 const auto& problem = gridVarsCache.problem();
104 Variables<freeFlowMassIndex> variables;
105 if constexpr (Concept::FVGridVariables<GridVariables<freeFlowMassIndex>>)
106 variables.update(elemSol_, problem, fvGeometry.element(), scv);
107 else
108 variables.update(elemSol_, problem, fvGeometry, ipData(fvGeometry, scv));
109 return variables;
110 }
111
112 ElementSolution elemSol_;
113 };
114
115 struct MomentumCouplingContextGlobalCaching
116 {
117 template<class GridVarsCache, class FvElementGeometry, class SubControlVolume>
118 const auto& vars(const GridVarsCache& gridVarsCache, const FvElementGeometry& fvGeometry, const SubControlVolume& scv) const
119 {
120 if constexpr (requires { gridVarsCache.volVars(scv); })
121 return gridVarsCache.volVars(scv);
122 else
123 return gridVarsCache.variables(scv);
124 }
125 };
126
127 using MomentumDiscretizationMethod = typename GridGeometry<freeFlowMomentumIndex>::DiscretizationMethod;
128 using MassDiscretizationMethod = typename GridGeometry<freeFlowMassIndex>::DiscretizationMethod;
129
130 template<std::size_t id> using IpData
131 = Dumux::CVFE::InterpolationPointData<typename GridView<id>::template Codim<0>::Entity::Geometry::LocalCoordinate,
132 typename GridView<id>::template Codim<0>::Entity::Geometry::GlobalCoordinate>;
133
134public:
135
136 static constexpr auto pressureIdx = Variables<freeFlowMassIndex>::Indices::pressureIdx;
137
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 {
149 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
150 gridVariables_ = gridVariables;
151 this->updateSolution(curSol);
152
153 computeCouplingStencils_();
154 }
155
157 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
158 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
159 GridVariablesTuple&& gridVariables,
160 const SolutionVector& curSol,
161 const SolutionVector& prevSol)
162 {
163 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
164 prevSol_ = &prevSol;
165 isTransient_ = true;
166 }
167
169 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
170 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
171 GridVariablesTuple&& gridVariables,
173 {
174 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
175 gridVariables_ = gridVariables;
176 this->attachSolution(curSol);
177
178 computeCouplingStencils_();
179 }
180
181 // \}
182
186 // \{
187
191 [[deprecated("This method will be removed after release (3.11). Use pressure(..., ipData) instead!")]]
192 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
193 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
194 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
195 const bool considerPreviousTimeStep = false) const
196 {
197 const auto& globalPos = scvf.ipGlobal();
198 const auto& localPos = element.geometry().local(globalPos);
199 return this->pressure(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
200 }
201
205 [[deprecated("This method will be removed after release (3.11). Use pressure(..., ipData) instead!")]]
206 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
207 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
208 const SubControlVolume<freeFlowMomentumIndex>& scv,
209 const bool considerPreviousTimeStep = false) const
210 {
211 return this->pressure(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
212 }
213
217 template <class IpData>
218 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
219 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
220 const IpData& ipData,
221 const bool considerPreviousTimeStep = false) const
222 {
223 assert(!(considerPreviousTimeStep && !this->isTransient_));
224 const auto& gg = this->problem(freeFlowMassIndex).gridGeometry();
225 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
226 : this->curSol(freeFlowMassIndex);
227 const auto elemSol = elementSolution(element, sol, gg);
228 return evalSolutionAtLocalPos(element, element.geometry(), gg, elemSol, ipData.local())[pressureIdx];
229 }
230
234 [[deprecated("This method will be removed after release (3.11). Use density(..., ipData) instead!")]]
235 Scalar density(const Element<freeFlowMomentumIndex>& element,
236 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
237 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
238 const bool considerPreviousTimeStep = false) const
239 {
240 const auto& globalPos = scvf.ipGlobal();
241 const auto& localPos = element.geometry().local(globalPos);
242 return this->density(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
243 }
244
248 [[deprecated("This method will be removed after release (3.11). Use density(..., ipData) instead!")]]
249 Scalar density(const Element<freeFlowMomentumIndex>& element,
250 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
251 const SubControlVolume<freeFlowMomentumIndex>& scv,
252 const bool considerPreviousTimeStep = false) const
253 {
254 return this->density(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
255 }
256
260 template <class IpData>
261 Scalar density(const Element<freeFlowMomentumIndex>& element,
262 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
263 const IpData& ipData,
264 const bool considerPreviousTimeStep = false) const
265 {
266 assert(!(considerPreviousTimeStep && !this->isTransient_));
267
268 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
269 : this->curSol(freeFlowMassIndex);
270
271 auto massFvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
272 massFvGeometry.bind(element);
273 const auto context = makeMomentumCouplingContext_(massFvGeometry, sol);
274 const auto& gridVarsCache = subDomainGridVars_(Dune::index_constant<freeFlowMassIndex>{}, considerPreviousTimeStep);
275
276 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
277 {
278 const auto eIdx = fvGeometry.elementIndex();
279 const auto& scv = massFvGeometry.scv(eIdx);
280
281 const auto& variables = context.vars(gridVarsCache, massFvGeometry, scv);
282
283 return variables.density();
284 }
285 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
286 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
287 {
288 // TODO: cache the shape values
289 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
290 const auto& localBasis = massFvGeometry.feLocalBasis();
291 std::vector<ShapeValue> shapeValues;
292 localBasis.evaluateFunction(ipData.local(), shapeValues);
293
294 Scalar rho = 0.0;
295 for (const auto& scv : scvs(massFvGeometry))
296 {
297 const auto& variables = context.vars(gridVarsCache, massFvGeometry, scv);
298 rho += variables.density()*shapeValues[scv.localDofIndex()][0];
299 }
300
301 return rho;
302 }
303 else
304 DUNE_THROW(Dune::NotImplemented,
305 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
306 );
307 }
308
312 [[deprecated("This method will be removed after release (3.11). Use effectiveViscosity(..., ipData) instead!")]]
313 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
314 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
315 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
316 const bool considerPreviousTimeStep = false) const
317 {
318 const auto& globalPos = scvf.ipGlobal();
319 const auto& localPos = element.geometry().local(globalPos);
320 return this->effectiveViscosity(element, fvGeometry, IpData<freeFlowMassIndex>(localPos, globalPos), considerPreviousTimeStep);
321 }
322
326 [[deprecated("This method will be removed after release (3.11). Use effectiveViscosity(..., ipData) instead!")]]
327 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
328 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
329 const SubControlVolume<freeFlowMomentumIndex>& scv,
330 const bool considerPreviousTimeStep = false) const
331 {
332 return this->effectiveViscosity(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
333 }
334
338 template <class IpData>
339 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
340 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
341 const IpData& ipData,
342 const bool considerPreviousTimeStep = false) const
343 {
344 assert(!(considerPreviousTimeStep && !this->isTransient_));
345
346 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
347 : this->curSol(freeFlowMassIndex);
348
349 auto massFvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
350 massFvGeometry.bind(element);
351 const auto context = makeMomentumCouplingContext_(massFvGeometry, sol);
352 const auto& gridVarsCache = subDomainGridVars_(Dune::index_constant<freeFlowMassIndex>{}, considerPreviousTimeStep);
353
354 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
355 {
356 const auto eIdx = fvGeometry.elementIndex();
357 const auto& scv = massFvGeometry.scv(eIdx);
358
359 const auto& variables = context.vars(gridVarsCache, massFvGeometry, scv);
360
361 return variables.viscosity();
362 }
363 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
364 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
365 {
366 // TODO: cache the shape values
367 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
368 const auto& localBasis = massFvGeometry.feLocalBasis();
369 std::vector<ShapeValue> shapeValues;
370 localBasis.evaluateFunction(ipData.local(), shapeValues);
371
372 Scalar mu = 0.0;
373 for (const auto& scv : scvs(massFvGeometry))
374 {
375 const auto& variables = context.vars(gridVarsCache, massFvGeometry, scv);
376 mu += variables.viscosity()*shapeValues[scv.localDofIndex()][0];
377 }
378
379 return mu;
380 }
381 else
382 DUNE_THROW(Dune::NotImplemented,
383 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
384 );
385 }
386
390 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
391 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
392 {
393 // TODO: optimize this function for tpfa where the scvf ip coincides with the dof location
394 auto fvGeometry = localView(this->problem(freeFlowMomentumIndex).gridGeometry());
395 fvGeometry.bindElement(element);
396
397 const auto& localBasis = fvGeometry.feLocalBasis();
398
399 std::vector<ShapeValue> shapeValues;
400 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
401 localBasis.evaluateFunction(ipLocal, shapeValues);
402
403 // interpolate velocity at scvf
404 VelocityVector velocity(0.0);
405 for (const auto& localDof : localDofs(fvGeometry))
406 velocity.axpy(shapeValues[localDof.index()][0], this->curSol(freeFlowMomentumIndex)[localDof.dofIndex()]);
407
408 return velocity;
409 }
410
414 VelocityVector elementVelocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry) const
415 {
416 auto momentumFvGeometry = localView(this->problem(freeFlowMomentumIndex).gridGeometry());
417 momentumFvGeometry.bindElement(fvGeometry.element());
418
419 const auto& localBasis = momentumFvGeometry.feLocalBasis();
420
421 // interpolate velocity at scvf
422 VelocityVector velocity(0.0);
423 std::vector<ShapeValue> shapeValues;
424 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
425
426 for (const auto& localDof : localDofs(momentumFvGeometry))
427 velocity.axpy(shapeValues[localDof.index()][0], this->curSol(freeFlowMomentumIndex)[localDof.dofIndex()]);
428
429 return velocity;
430 }
431
435 template <class IpData>
436 VelocityVector velocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry,
437 const IpData& ipData,
438 const bool considerPreviousTimeStep = false) const
439 {
440 assert(!(considerPreviousTimeStep && !this->isTransient_));
441
442 const auto& element = fvGeometry.element();
443 const auto& gg = this->problem(freeFlowMomentumIndex).gridGeometry();
444 auto momentumFvGeometry = localView(gg);
445 momentumFvGeometry.bindElement(fvGeometry.element());
446
447 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMomentumIndex]
449
450 const auto elemSol = elementSolution(element, sol, gg);
451 return evalSolutionAtLocalPos(element, element.geometry(), gg, elemSol, ipData.local());
452 }
453
458 template<std::size_t j>
459 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
460 const Element<freeFlowMomentumIndex>& elementI,
461 const SubControlVolume<freeFlowMomentumIndex>& scvI,
462 Dune::index_constant<j> domainJ) const
463 { return emptyStencil_; }
464
479 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
480 const Element<freeFlowMassIndex>& elementI,
481 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
482 {
483 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
484 return massAndEnergyToMomentumStencils_[eIdx];
485 }
486
495 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
496 const Element<freeFlowMomentumIndex>& elementI,
497 Dune::index_constant<freeFlowMassIndex> domainJ) const
498 {
499 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
500 return momentumToMassAndEnergyStencils_[eIdx];
501 }
502
503 // \}
504
508 // \{
509
529 template<std::size_t i, std::size_t j, class LocalAssemblerI>
530 void updateCouplingContext(Dune::index_constant<i> domainI,
531 const LocalAssemblerI& localAssemblerI,
532 Dune::index_constant<j> domainJ,
533 std::size_t dofIdxGlobalJ,
534 const PrimaryVariables<j>& priVarsJ,
535 int pvIdxJ)
536 {
537 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
538
539 if constexpr (GridVariablesCache<freeFlowMassIndex>::cachingEnabled)
540 {
541 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
542 {
543 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
544 {
545 const auto& problem = this->problem(domainJ);
546 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
547 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
548 auto fvGeometry = localView(problem.gridGeometry());
549 fvGeometry.bind(deflectedElement);
550 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
551
552 if constexpr (Concept::FVGridVariables<GridVariables<freeFlowMassIndex>>)
553 subDomainVariables_(Dune::index_constant<freeFlowMassIndex>{}, /*current*/true, scv).update(std::move(elemSol), problem, deflectedElement, scv);
554 else
555 subDomainVariables_(Dune::index_constant<freeFlowMassIndex>{}, /*current*/true, ipData(fvGeometry, scv)).update(std::move(elemSol), problem, deflectedElement, scv);
556 }
557 }
558 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
559 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
560 {
561 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
562 {
563 const auto& problem = this->problem(domainJ);
564 const auto deflectedElementIdx = problem.gridGeometry().elementMapper().index(localAssemblerI.element());
565 const auto& deflectedElement = problem.gridGeometry().element(deflectedElementIdx);
566 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
567 auto fvGeometry = localView(problem.gridGeometry());
568 fvGeometry.bind(deflectedElement);
569
570 // ToDo: Replace once all mass models are also working with local dofs
571 for (const auto& scv : scvs(fvGeometry))
572 {
573 if (scv.dofIndex() == dofIdxGlobalJ)
574 {
575 if constexpr (Concept::FVGridVariables<GridVariables<freeFlowMassIndex>>)
576 this->subDomainVariables_(Dune::index_constant<freeFlowMassIndex>{}, /*current*/true, scv).update(std::move(elemSol), problem, deflectedElement, scv);
577 else
578 this->subDomainVariables_(Dune::index_constant<freeFlowMassIndex>{}, /*current*/true, scv).update(std::move(elemSol), problem, fvGeometry, ipData(fvGeometry, scv));
579 }
580 }
581 }
582 }
583 else
584 DUNE_THROW(Dune::NotImplemented,
585 "Context update for discretization scheme " << MassDiscretizationMethod{}
586 );
587 }
588 }
589
590 // \}
591
596 {
597 if constexpr (MomentumDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
598 {
599 // use coloring of the mass discretization for both domains
600 // the diamond coloring is a subset (minimum amount of colors) of cctpfa/box coloring
601 elementSets_ = computeColoring(this->problem(freeFlowMassIndex).gridGeometry()).sets;
602 }
603 else
604 {
605 // use coloring of the momentum discretization for both domains
606 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
607 }
608 }
609
616 template<std::size_t i, class AssembleElementFunc>
617 void assembleMultithreaded(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
618 {
619 if (elementSets_.empty())
620 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
621
622 // make this element loop run in parallel
623 // for this we have to color the elements so that we don't get
624 // race conditions when writing into the global matrix
625 // each color can be assembled using multiple threads
626 const auto& grid = this->problem(freeFlowMassIndex).gridGeometry().gridView().grid();
627 for (const auto& elements : elementSets_)
628 {
629 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
630 {
631 const auto element = grid.entity(elements[eIdx]);
632 assembleElement(element);
633 });
634 }
635 }
636
637private:
638 template<std::size_t i>
639 const auto& subDomainGridVars_(Dune::index_constant<i> domainIdx, bool current) const
640 {
641 if constexpr (Concept::FVGridVariables<GridVariables<i>>)
642 return current ? gridVars_(domainIdx).curGridVolVars()
643 : gridVars_(domainIdx).prevGridVolVars();
644 else
645 return current ? gridVars_(domainIdx).curGridVars()
646 : gridVars_(domainIdx).prevGridVars();
647 }
648
649 template<std::size_t i>
650 auto& subDomainGridVars_(Dune::index_constant<i> domainIdx, bool current)
651 {
652 if constexpr (Concept::FVGridVariables<GridVariables<i>>)
653 return current ? gridVars_(domainIdx).curGridVolVars()
654 : gridVars_(domainIdx).prevGridVolVars();
655 else
656 return current ? gridVars_(domainIdx).curGridVars()
657 : gridVars_(domainIdx).prevGridVars();
658 }
659
660 template<std::size_t i, class Scv>
661 auto& subDomainVariables_(Dune::index_constant<i> domainIdx, bool current, const Scv& scv)
662 {
663 if constexpr (Concept::FVGridVariables<GridVariables<i>>)
664 return subDomainGridVars_(domainIdx, current).volVars(scv);
665 else
666 return subDomainGridVars_(domainIdx, current).variables(scv);
667 }
668
669 template<std::size_t i, class Scv>
670 const auto& subDomainVariables_(Dune::index_constant<i> domainIdx, bool current, const Scv& scv) const
671 {
672 if constexpr (Concept::FVGridVariables<GridVariables<i>>)
673 return subDomainGridVars_(domainIdx, current).volVars(scv);
674 else
675 return subDomainGridVars_(domainIdx, current).variables(scv);
676 }
677
678 template<class SolutionVector>
679 auto makeMomentumCouplingContext_(const FVElementGeometry<freeFlowMassIndex>& fvGeometry,
680 const SolutionVector& sol) const
681 {
682 if constexpr (GridVariablesCache<freeFlowMassIndex>::cachingEnabled)
683 return MomentumCouplingContextGlobalCaching{};
684 else
685 return MomentumCouplingContextNoCaching{elementSolution(fvGeometry.element(), sol, fvGeometry.gridGeometry())};
686 }
687
692 template<std::size_t i>
693 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
694 {
695 if (std::get<i>(gridVariables_))
696 return *std::get<i>(gridVariables_);
697 else
698 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
699 }
700
705 template<std::size_t i>
706 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
707 {
708 if (std::get<i>(gridVariables_))
709 return *std::get<i>(gridVariables_);
710 else
711 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
712 }
713
714 void computeCouplingStencils_()
715 {
716 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
717 const auto& massGridGeometry = this->problem(freeFlowMassIndex).gridGeometry();
718 auto momentumFvGeometry = localView(momentumGridGeometry);
719 auto massFvGeometry = localView(massGridGeometry);
720
721 massAndEnergyToMomentumStencils_.clear();
722 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
723
724 momentumToMassAndEnergyStencils_.clear();
725 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
726
727 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
728 bool hasCubeElements = false;
729
730 for (const auto& element : elements(momentumGridGeometry.gridView()))
731 {
732 if(element.type().isCube())
733 hasCubeElements = true;
734
735 momentumFvGeometry.bindElement(element);
736 massFvGeometry.bindElement(element);
737 const auto eIdx = momentumFvGeometry.elementIndex();
738
739 for (const auto& localDof : localDofs(momentumFvGeometry))
740 massAndEnergyToMomentumStencils_[eIdx].push_back(localDof.dofIndex());
741
742 // ToDo: Replace once all mass models are also working with local dofs
743 for (const auto& scv : scvs(massFvGeometry))
744 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
745 }
746
747 // Print warning for pq1bubble scheme on cube elements if not using the hybrid variant
748 if constexpr ((MomentumDiscretizationMethod{} == DiscretizationMethods::pq1bubble))
749 {
750 if(hasCubeElements && !GridGeometry<freeFlowMomentumIndex>::enableHybridCVFE)
751 {
752 std::cerr << "Warning: Coupled Navier-Stokes problem on cube elements uses non-hybrid pq1bubble. "
753 << "The hybrid variant is recommended because it implements two bubble functions for stability reasons."
754 << std::endl;
755 }
756 }
757 }
758
759 CouplingStencilType emptyStencil_;
760 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
761 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
762
764 GridVariablesTuple gridVariables_;
765
766 const SolutionVector* prevSol_;
767 bool isTransient_;
768
769 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
770};
771
772namespace Detail {
773
774// declaration (specialize for different discretization types)
775template<class Traits, class DiscretizationMethod = typename Detail::MomentumDiscretizationMethod<Traits>::type>
777
778// multi-threading is not supported because we have only one coupling context instance and a mutable cache
779template<class Traits, class D>
780struct CouplingManagerSupportsMultithreadedAssemblySelector<Traits, DiscretizationMethods::CVFE<D>>
781{ using type = std::true_type; };
782
783} // end namespace Detail
784
786template<class T>
789{};
790
791} // end namespace Dumux
792
793#endif
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:31
const LocalPosition & local() const
The local position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:46
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_cvfe.hh:54
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:218
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:390
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:144
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:169
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:327
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_cvfe.hh:57
static constexpr auto pressureIdx
Definition: couplingmanager_cvfe.hh:136
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:479
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:495
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:235
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:261
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:313
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_cvfe.hh:414
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:249
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:339
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:206
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_cvfe.hh:617
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:192
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:157
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_cvfe.hh:595
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_cvfe.hh:58
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:436
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:459
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:332
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:297
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:319
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:348
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:229
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:530
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