version 3.8
couplingmanager_pq1bubble.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-FileCopyrightInfo: 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_PQ1BUBBLE_HH
13#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_PQ1BUBBLE_HH
14
15#warning "This file is deprecated and will be removed after 3.7. Use CVFEFreeFlowCouplingManager."
16
17#include <memory>
18#include <tuple>
19#include <vector>
20#include <deque>
21
22#include <dune/common/exceptions.hh>
23#include <dune/common/indices.hh>
24#include <dune/common/float_cmp.hh>
25#include <dune/geometry/referenceelements.hh>
26
29
33
36
39
40namespace Dumux {
41
47template<class Traits>
49: public CouplingManager<Traits>
50{
52public:
53 static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index();
54 static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index();
55
56 // this can be used if the coupling manager is used inside a meta-coupling manager (e.g. multi-binary)
57 // to manager the solution vector storage outside this class
59private:
60 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
61 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
62 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
63 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
64 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
65 template<std::size_t id> using ElementSeed = typename GridView<id>::Grid::template Codim<0>::EntitySeed;
66 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
67 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
68 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
69 template<std::size_t id> using GridVariables = typename Traits::template SubDomain<id>::GridVariables;
70 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
71 template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache;
72 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
73 template<std::size_t id> using VolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::VolumeVariables>;
74
75 using Scalar = typename Traits::Scalar;
76 using SolutionVector = typename Traits::SolutionVector;
77
78 using CouplingStencilType = std::vector<std::size_t>;
79
80 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
81
82 using FluidSystem = typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
83
84 using VelocityVector = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
85 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
86
87 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
88
89 struct MomentumCouplingContext
90 {
91 FVElementGeometry<freeFlowMassIndex> fvGeometry;
92 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
93 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
94 std::size_t eIdx;
95 };
96
97 struct MassAndEnergyCouplingContext
98 {
99 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f, const std::size_t i)
100 : fvGeometry(std::move(f))
101 , eIdx(i)
102 {}
103
104 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
105 std::size_t eIdx;
106 };
107
108 using MomentumDiscretizationMethod = typename GridGeometry<freeFlowMomentumIndex>::DiscretizationMethod;
109 using MassDiscretizationMethod = typename GridGeometry<freeFlowMassIndex>::DiscretizationMethod;
110
111public:
112
113 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
114
118 // \{
119
121 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
122 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
123 GridVariablesTuple&& gridVariables,
124 const SolutionVector& curSol)
125 {
126 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
127 gridVariables_ = gridVariables;
128 this->updateSolution(curSol);
129
130 computeCouplingStencils_();
131 }
132
134 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
135 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
136 GridVariablesTuple&& gridVariables,
137 const SolutionVector& curSol,
138 const SolutionVector& prevSol)
139 {
140 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
141 prevSol_ = &prevSol;
142 isTransient_ = true;
143 }
144
146 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
147 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
148 GridVariablesTuple&& gridVariables,
150 {
151 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
152 gridVariables_ = gridVariables;
153 this->attachSolution(curSol);
154
155 computeCouplingStencils_();
156 }
157
158 // \}
159
163 // \{
164
168 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
169 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
170 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
171 const bool considerPreviousTimeStep = false) const
172 {
173 assert(!(considerPreviousTimeStep && !this->isTransient_));
174 const auto& gg = this->problem(freeFlowMassIndex).gridGeometry();
175 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
176 : this->curSol(freeFlowMassIndex);
177 const auto elemSol = elementSolution(element, sol, gg);
178 return evalSolution(element, element.geometry(), gg, elemSol, scvf.ipGlobal())[pressureIdx];
179 }
180
184 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
185 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
186 const SubControlVolume<freeFlowMomentumIndex>& scv,
187 const bool considerPreviousTimeStep = false) const
188 {
189 assert(!(considerPreviousTimeStep && !this->isTransient_));
190 const auto& gg = this->problem(freeFlowMassIndex).gridGeometry();
191 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
192 : this->curSol(freeFlowMassIndex);
193 const auto elemSol = elementSolution(element, sol, gg);
194 return evalSolution(element, element.geometry(), gg, elemSol, scv.dofPosition())[pressureIdx];
195 }
196
200 Scalar density(const Element<freeFlowMomentumIndex>& element,
201 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
202 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
203 const bool considerPreviousTimeStep = false) const
204 {
205 assert(!(considerPreviousTimeStep && !this->isTransient_));
206 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
207
208 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
209 {
210 const auto eIdx = fvGeometry.elementIndex();
211 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
212
213 const auto& volVars = considerPreviousTimeStep ?
214 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
215 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
216
217 return volVars.density();
218 }
219 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
220 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
221 {
222 // TODO: cache the shape values when Box method is used
223 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
224 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
225 std::vector<ShapeValue> shapeValues;
226 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
227 localBasis.evaluateFunction(ipLocal, shapeValues);
228
229 Scalar rho = 0.0;
230 for (const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
231 {
232 const auto& volVars = considerPreviousTimeStep ?
233 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
234 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
235 rho += volVars.density()*shapeValues[scv.indexInElement()][0];
236 }
237
238 return rho;
239 }
240 else
241 DUNE_THROW(Dune::NotImplemented,
242 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
243 );
244 }
245
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 assert(!(considerPreviousTimeStep && !this->isTransient_));
255 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
256
257 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
258 {
259 const auto eIdx = scv.elementIndex();
260 const auto& scvI = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
261
262 const auto& volVars = considerPreviousTimeStep ?
263 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
264 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
265
266 return volVars.density();
267 }
268 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
269 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
270 {
271 // TODO: cache the shape values when Box method is used
272 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
273 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
274 std::vector<ShapeValue> shapeValues;
275 const auto ipLocal = element.geometry().local(scv.dofPosition());
276 localBasis.evaluateFunction(ipLocal, shapeValues);
277
278 Scalar rho = 0.0;
279 for (const auto& scvI : scvs(this->momentumCouplingContext_()[0].fvGeometry))
280 {
281 const auto& volVars = considerPreviousTimeStep ?
282 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
283 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
284 rho += volVars.density()*shapeValues[scvI.indexInElement()][0];
285 }
286 return rho;
287 }
288 else
289 DUNE_THROW(Dune::NotImplemented,
290 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
291 );
292 }
293
297 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
298 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
299 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
300 const bool considerPreviousTimeStep = false) const
301 {
302 assert(!(considerPreviousTimeStep && !this->isTransient_));
303 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
304
305 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
306 {
307 const auto eIdx = fvGeometry.elementIndex();
308 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
309 const auto& volVars = considerPreviousTimeStep ?
310 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
311 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
312 return volVars.viscosity();
313 }
314 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
315 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
316 {
317 // TODO: cache the shape values when Box method is used
318 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
319 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
320 std::vector<ShapeValue> shapeValues;
321 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
322 localBasis.evaluateFunction(ipLocal, shapeValues);
323
324 Scalar mu = 0.0;
325 for (const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
326 {
327 const auto& volVars = considerPreviousTimeStep ?
328 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
329 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
330 mu += volVars.viscosity()*shapeValues[scv.indexInElement()][0];
331 }
332
333 return mu;
334 }
335 else
336 DUNE_THROW(Dune::NotImplemented,
337 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
338 );
339 }
340
344 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
345 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
346 const SubControlVolume<freeFlowMomentumIndex>& scv,
347 const bool considerPreviousTimeStep = false) const
348 {
349 assert(!(considerPreviousTimeStep && !this->isTransient_));
350 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
351
352 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
353 {
354 const auto eIdx = fvGeometry.elementIndex();
355 const auto& scvI = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
356 const auto& volVars = considerPreviousTimeStep ?
357 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
358 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
359 return volVars.viscosity();
360 }
361 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
362 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
363 {
364 // TODO: cache the shape values when Box method is used
365 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
366 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
367 std::vector<ShapeValue> shapeValues;
368 const auto ipLocal = element.geometry().local(scv.dofPosition());
369 localBasis.evaluateFunction(ipLocal, shapeValues);
370
371 Scalar mu = 0.0;
372 for (const auto& scvI : scvs(this->momentumCouplingContext_()[0].fvGeometry))
373 {
374 const auto& volVars = considerPreviousTimeStep ?
375 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
376 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
377 mu += volVars.viscosity()*shapeValues[scvI.indexInElement()][0];
378 }
379
380 return mu;
381 }
382 else
383 DUNE_THROW(Dune::NotImplemented,
384 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
385 );
386 }
387
391 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
392 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
393 {
394 // TODO: optimize this function for tpfa where the scvf ip coincides with the dof location
395
396 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(element);
397 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, eIdx);
398
399 const auto& fvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
400 const auto& localBasis = fvGeometry.feLocalBasis();
401
402 std::vector<ShapeValue> shapeValues;
403 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
404 localBasis.evaluateFunction(ipLocal, shapeValues);
405
406 // interpolate velocity at scvf
407 VelocityVector velocity(0.0);
408 for (const auto& scv : scvs(fvGeometry))
409 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(freeFlowMomentumIndex)[scv.dofIndex()]);
410
411 return velocity;
412 }
413
417 VelocityVector elementVelocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry) const
418 {
419 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), fvGeometry.element());
420
421 const auto& momentumFvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
422 const auto& localBasis = momentumFvGeometry.feLocalBasis();
423
424 // interpolate velocity at scvf
425 VelocityVector velocity(0.0);
426 std::vector<ShapeValue> shapeValues;
427 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
428
429 for (const auto& scv : scvs(momentumFvGeometry))
430 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(freeFlowMomentumIndex)[scv.dofIndex()]);
431
432 return velocity;
433 }
434
439 template<std::size_t j>
440 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
441 const Element<freeFlowMomentumIndex>& elementI,
442 const SubControlVolume<freeFlowMomentumIndex>& scvI,
443 Dune::index_constant<j> domainJ) const
444 { return emptyStencil_; }
445
460 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
461 const Element<freeFlowMassIndex>& elementI,
462 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
463 {
464 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
465 return massAndEnergyToMomentumStencils_[eIdx];
466 }
467
476 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
477 const Element<freeFlowMomentumIndex>& elementI,
478 Dune::index_constant<freeFlowMassIndex> domainJ) const
479 {
480 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
481 return momentumToMassAndEnergyStencils_[eIdx];
482 }
483
484 // \}
485
489 // \{
490
510 template<std::size_t i, std::size_t j, class LocalAssemblerI>
511 void updateCouplingContext(Dune::index_constant<i> domainI,
512 const LocalAssemblerI& localAssemblerI,
513 Dune::index_constant<j> domainJ,
514 std::size_t dofIdxGlobalJ,
515 const PrimaryVariables<j>& priVarsJ,
516 int pvIdxJ)
517 {
518 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
519
520 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
521 {
522 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
523 {
524 bindCouplingContext_(domainI, localAssemblerI.element());
525
526 const auto& problem = this->problem(domainJ);
527 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
528 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
529 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
530 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
531
532 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
533 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
534 else
535 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
536 }
537 }
538 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
539 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
540 {
541 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
542 {
543 bindCouplingContext_(domainI, localAssemblerI.element());
544
545 const auto& problem = this->problem(domainJ);
546 const auto& deflectedElement = problem.gridGeometry().element(this->momentumCouplingContext_()[0].eIdx);
547 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
548 const auto& fvGeometry = this->momentumCouplingContext_()[0].fvGeometry;
549
550 for (const auto& scv : scvs(fvGeometry))
551 {
552 if(scv.dofIndex() == dofIdxGlobalJ)
553 {
554 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
555 this->gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
556 else
557 this->momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
558 }
559 }
560 }
561 }
562 else
563 DUNE_THROW(Dune::NotImplemented,
564 "Context update for discretization scheme " << MassDiscretizationMethod{}
565 );
566 }
567
568 // \}
569
574 {
575 // use coloring of the momentum discretization for both domains
576 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
577 }
578
585 template<std::size_t i, class AssembleElementFunc>
586 void assembleMultithreaded(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
587 {
588 if (elementSets_.empty())
589 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
590
591 // make this element loop run in parallel
592 // for this we have to color the elements so that we don't get
593 // race conditions when writing into the global matrix
594 // each color can be assembled using multiple threads
595 const auto& grid = this->problem(freeFlowMassIndex).gridGeometry().gridView().grid();
596 for (const auto& elements : elementSets_)
597 {
598 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
599 {
600 const auto element = grid.entity(elements[eIdx]);
601 assembleElement(element);
602 });
603 }
604 }
605
606private:
607 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
608 const Element<freeFlowMomentumIndex>& elementI) const
609 {
610 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
611 if (momentumCouplingContext_().empty())
612 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI));
613 else
614 bindCouplingContext_(domainI, elementI, momentumCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
615 }
616
617 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
618 const Element<freeFlowMomentumIndex>& elementI,
619 const std::size_t eIdx) const
620 {
621 if (momentumCouplingContext_().empty())
622 {
623 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
624 fvGeometry.bind(elementI);
625
626 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
627 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
628
629 auto prevElemVolVars = isTransient_ ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
630 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
631
632 if (isTransient_)
633 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[freeFlowMassIndex]);
634
635 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
636 }
637 else if (eIdx != momentumCouplingContext_()[0].eIdx)
638 {
639 momentumCouplingContext_()[0].eIdx = eIdx;
640 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
641 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->curSol(freeFlowMassIndex));
642
643 if (isTransient_)
644 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[freeFlowMassIndex]);
645 }
646 }
647
648 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
649 const Element<freeFlowMassIndex>& elementI) const
650 {
651 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
652 if (massAndEnergyCouplingContext_().empty())
653 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI));
654 else
655 bindCouplingContext_(domainI, elementI, massAndEnergyCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
656 }
657
658 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
659 const Element<freeFlowMassIndex>& elementI,
660 const std::size_t eIdx) const
661 {
662 if (massAndEnergyCouplingContext_().empty())
663 {
664 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
665 auto fvGeometry = localView(gridGeometry);
666 fvGeometry.bindElement(elementI);
667 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
668 }
669 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
670 {
671 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
672 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
673 }
674 }
675
680 template<std::size_t i>
681 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
682 {
683 if (std::get<i>(gridVariables_))
684 return *std::get<i>(gridVariables_);
685 else
686 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
687 }
688
693 template<std::size_t i>
694 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
695 {
696 if (std::get<i>(gridVariables_))
697 return *std::get<i>(gridVariables_);
698 else
699 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
700 }
701
702
703 void computeCouplingStencils_()
704 {
705 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
706 const auto& massGridGeometry = this->problem(freeFlowMassIndex).gridGeometry();
707 auto momentumFvGeometry = localView(momentumGridGeometry);
708 auto massFvGeometry = localView(massGridGeometry);
709
710 massAndEnergyToMomentumStencils_.clear();
711 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
712
713 momentumToMassAndEnergyStencils_.clear();
714 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
715
716 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
717
718 for (const auto& element : elements(momentumGridGeometry.gridView()))
719 {
720 momentumFvGeometry.bindElement(element);
721 massFvGeometry.bindElement(element);
722 const auto eIdx = momentumFvGeometry.elementIndex();
723
724 for (const auto& scv : scvs(momentumFvGeometry))
725 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
726
727 for (const auto& scv : scvs(massFvGeometry))
728 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
729 }
730 }
731
732 CouplingStencilType emptyStencil_;
733 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
734 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
735
736 // the coupling context exists for each thread
737 // TODO this is a bad pattern, just like mutable caches
738 // we should really construct and pass the context and not store it globally
739 std::vector<MomentumCouplingContext>& momentumCouplingContext_() const
740 {
741 thread_local static std::vector<MomentumCouplingContext> c;
742 return c;
743 }
744
745 // the coupling context exists for each thread
746 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_() const
747 {
748 thread_local static std::vector<MassAndEnergyCouplingContext> c;
749 return c;
750 }
751
753 GridVariablesTuple gridVariables_;
754
755 const SolutionVector* prevSol_;
756 bool isTransient_;
757
758 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
759};
760
762template<class T>
764: public std::false_type {};
765
766} // end namespace Dumux
767
768#endif
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:322
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:287
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:309
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:64
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:338
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:219
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:71
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_pq1bubble.hh:50
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_pq1bubble.hh:53
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_pq1bubble.hh:573
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_pq1bubble.hh:54
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_pq1bubble.hh:344
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, Dune::index_constant< freeFlowMassIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: couplingmanager_pq1bubble.hh:476
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_pq1bubble.hh:168
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_pq1bubble.hh:297
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol, const SolutionVector &prevSol)
use as regular coupling manager in a transient setting
Definition: couplingmanager_pq1bubble.hh:134
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: couplingmanager_pq1bubble.hh:249
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< j > domainJ) const
The coupling stencil of domain I, i.e. which domain J DOFs the given domain I element's residual depe...
Definition: couplingmanager_pq1bubble.hh:440
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMassIndex > domainI, const Element< freeFlowMassIndex > &elementI, Dune::index_constant< freeFlowMomentumIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: couplingmanager_pq1bubble.hh:460
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol)
Methods to be accessed by main.
Definition: couplingmanager_pq1bubble.hh:121
static constexpr auto pressureIdx
Definition: couplingmanager_pq1bubble.hh:113
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume face.
Definition: couplingmanager_pq1bubble.hh:200
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_pq1bubble.hh:184
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, typename ParentType::SolutionVectorStorage &curSol)
use as binary coupling manager in multi model context
Definition: couplingmanager_pq1bubble.hh:146
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: couplingmanager_pq1bubble.hh:391
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_pq1bubble.hh:417
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_pq1bubble.hh:586
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Type traits.
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
PrimaryVariables evalSolution(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:152
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
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: couplingmanager_pq1bubble.hh:511
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
The available discretization methods in Dumux.
The interface of the coupling manager for multi domain problems.
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
Parallel for loop (multithreading)
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78