3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_PQ1BUBBLE_HH
25#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_PQ1BUBBLE_HH
26
27#include <memory>
28#include <tuple>
29#include <vector>
30#include <deque>
31
32#include <dune/common/exceptions.hh>
33#include <dune/common/indices.hh>
34#include <dune/common/float_cmp.hh>
35#include <dune/geometry/referenceelements.hh>
36
39
43
46
49
50namespace Dumux {
51
57template<class Traits>
59: public CouplingManager<Traits>
60{
62public:
63 static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index();
64 static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index();
65
66 // this can be used if the coupling manager is used inside a meta-coupling manager (e.g. multi-binary)
67 // to manager the solution vector storage outside this class
69private:
70 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
71 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
72 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
73 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
74 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
75 template<std::size_t id> using ElementSeed = typename GridView<id>::Grid::template Codim<0>::EntitySeed;
76 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
77 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
78 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
79 template<std::size_t id> using GridVariables = typename Traits::template SubDomain<id>::GridVariables;
80 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
81 template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache;
82 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
83 template<std::size_t id> using VolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::VolumeVariables>;
84
85 using Scalar = typename Traits::Scalar;
86 using SolutionVector = typename Traits::SolutionVector;
87
88 using CouplingStencilType = std::vector<std::size_t>;
89
90 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
91
92 using FluidSystem = typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
93
94 using VelocityVector = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
95 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
96
97 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
98
99 struct MomentumCouplingContext
100 {
101 FVElementGeometry<freeFlowMassIndex> fvGeometry;
102 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
103 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
104 std::size_t eIdx;
105 };
106
107 struct MassAndEnergyCouplingContext
108 {
109 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f, const std::size_t i)
110 : fvGeometry(std::move(f))
111 , eIdx(i)
112 {}
113
114 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
115 std::size_t eIdx;
116 };
117
118 using MomentumDiscretizationMethod = typename GridGeometry<freeFlowMomentumIndex>::DiscretizationMethod;
119 using MassDiscretizationMethod = typename GridGeometry<freeFlowMassIndex>::DiscretizationMethod;
120
121public:
122
123 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
124
128 // \{
129
131 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
132 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
133 GridVariablesTuple&& gridVariables,
134 const SolutionVector& curSol)
135 {
136 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
137 gridVariables_ = gridVariables;
138 this->updateSolution(curSol);
139
140 computeCouplingStencils_();
141 }
142
144 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
145 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
146 GridVariablesTuple&& gridVariables,
147 const SolutionVector& curSol,
148 const SolutionVector& prevSol)
149 {
150 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
151 prevSol_ = &prevSol;
152 isTransient_ = true;
153 }
154
156 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
157 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
158 GridVariablesTuple&& gridVariables,
160 {
161 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
162 gridVariables_ = gridVariables;
163 this->attachSolution(curSol);
164
165 computeCouplingStencils_();
166 }
167
168 // \}
169
173 // \{
174
178 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
179 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
180 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
181 {
182 const auto& gg = this->problem(freeFlowMassIndex).gridGeometry();
183 const auto elemSol = elementSolution(element, this->curSol(freeFlowMassIndex), gg);
184 return evalSolution(element, element.geometry(), gg, elemSol, scvf.ipGlobal())[pressureIdx];
185 }
186
190 Scalar density(const Element<freeFlowMomentumIndex>& element,
191 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
192 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
193 const bool considerPreviousTimeStep = false) const
194 {
195 assert(!(considerPreviousTimeStep && !this->isTransient_));
196 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
197
198 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
199 {
200 const auto eIdx = fvGeometry.elementIndex();
201 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
202
203 const auto& volVars = considerPreviousTimeStep ?
204 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
205 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
206
207 return volVars.density();
208 }
209 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
210 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
211 {
212 // TODO: cache the shape values when Box method is used
213 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
214 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
215 std::vector<ShapeValue> shapeValues;
216 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
217 localBasis.evaluateFunction(ipLocal, shapeValues);
218
219 Scalar rho = 0.0;
220 for (const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
221 {
222 const auto& volVars = considerPreviousTimeStep ?
223 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
224 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
225 rho += volVars.density()*shapeValues[scv.indexInElement()][0];
226 }
227
228 return rho;
229 }
230 else
231 DUNE_THROW(Dune::NotImplemented,
232 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
233 );
234 }
235
239 Scalar density(const Element<freeFlowMomentumIndex>& element,
240 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
241 const SubControlVolume<freeFlowMomentumIndex>& scv,
242 const bool considerPreviousTimeStep = false) const
243 {
244 assert(!(considerPreviousTimeStep && !this->isTransient_));
245 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
246
247 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
248 {
249 const auto eIdx = scv.elementIndex();
250 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
251
252 const auto& volVars = considerPreviousTimeStep ?
253 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
254 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
255
256 return volVars.density();
257 }
258 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
259 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
260 {
261 // TODO: cache the shape values when Box method is used
262 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
263 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
264 std::vector<ShapeValue> shapeValues;
265 const auto ipLocal = element.geometry().local(scv.dofPosition());
266 localBasis.evaluateFunction(ipLocal, shapeValues);
267
268 Scalar rho = 0.0;
269 for (const auto& scvI : scvs(this->momentumCouplingContext_()[0].fvGeometry))
270 {
271 const auto& volVars = considerPreviousTimeStep ?
272 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
273 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
274 rho += volVars.density()*shapeValues[scvI.indexInElement()][0];
275 }
276 return rho;
277 }
278 else
279 DUNE_THROW(Dune::NotImplemented,
280 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
281 );
282 }
283
287 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
288 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
289 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
290 {
291 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
292
293 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
294 {
295 const auto eIdx = fvGeometry.elementIndex();
296 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
297 const auto& volVars = this->momentumCouplingContext_()[0].curElemVolVars[scv];
298 return volVars.viscosity();
299 }
300 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
301 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
302 {
303 // TODO: cache the shape values when Box method is used
304 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
305 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
306 std::vector<ShapeValue> shapeValues;
307 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
308 localBasis.evaluateFunction(ipLocal, shapeValues);
309
310 Scalar mu = 0.0;
311 for (const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
312 {
313 const auto& volVars = this->momentumCouplingContext_()[0].curElemVolVars[scv];
314 mu += volVars.viscosity()*shapeValues[scv.indexInElement()][0];
315 }
316
317 return mu;
318 }
319 else
320 DUNE_THROW(Dune::NotImplemented,
321 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
322 );
323
324 }
325
329 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
330 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
331 {
332 // TODO: optimize this function for tpfa where the scvf ip coincides with the dof location
333
334 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element);
335
336 const auto& fvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
337 const auto& localBasis = fvGeometry.feLocalBasis();
338
339 // interpolate velocity at scvf
340 VelocityVector velocity(0.0);
341 std::vector<ShapeValue> shapeValues;
342 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
343 localBasis.evaluateFunction(ipLocal, shapeValues);
344
345 for (const auto& scv : scvs(fvGeometry))
346 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(freeFlowMomentumIndex)[scv.dofIndex()]);
347
348 return velocity;
349 }
350
354 VelocityVector elementVelocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry) const
355 {
356 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), fvGeometry.element());
357
358 const auto& momentumFvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
359 const auto& localBasis = momentumFvGeometry.feLocalBasis();
360
361 // interpolate velocity at scvf
362 VelocityVector velocity(0.0);
363 std::vector<ShapeValue> shapeValues;
364 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
365
366 for (const auto& scv : scvs(momentumFvGeometry))
367 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(freeFlowMomentumIndex)[scv.dofIndex()]);
368
369 return velocity;
370 }
371
376 template<std::size_t j>
377 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
378 const Element<freeFlowMomentumIndex>& elementI,
379 const SubControlVolume<freeFlowMomentumIndex>& scvI,
380 Dune::index_constant<j> domainJ) const
381 { return emptyStencil_; }
382
397 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
398 const Element<freeFlowMassIndex>& elementI,
399 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
400 {
401 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
402 return massAndEnergyToMomentumStencils_[eIdx];
403 }
404
413 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
414 const Element<freeFlowMomentumIndex>& elementI,
415 Dune::index_constant<freeFlowMassIndex> domainJ) const
416 {
417 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
418 return momentumToMassAndEnergyStencils_[eIdx];
419 }
420
421 // \}
422
426 // \{
427
447 template<std::size_t i, std::size_t j, class LocalAssemblerI>
448 void updateCouplingContext(Dune::index_constant<i> domainI,
449 const LocalAssemblerI& localAssemblerI,
450 Dune::index_constant<j> domainJ,
451 std::size_t dofIdxGlobalJ,
452 const PrimaryVariables<j>& priVarsJ,
453 int pvIdxJ)
454 {
455 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
456
457 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
458 {
459 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
460 {
461 bindCouplingContext_(domainI, localAssemblerI.element());
462
463 const auto& problem = this->problem(domainJ);
464 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
465 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
466 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
467 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
468
469 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
470 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
471 else
472 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
473 }
474 }
475 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
476 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
477 {
478 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
479 {
480 bindCouplingContext_(domainI, localAssemblerI.element());
481
482 const auto& problem = this->problem(domainJ);
483 const auto& deflectedElement = problem.gridGeometry().element(this->momentumCouplingContext_()[0].eIdx);
484 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
485 const auto& fvGeometry = this->momentumCouplingContext_()[0].fvGeometry;
486
487 for (const auto& scv : scvs(fvGeometry))
488 {
489 if(scv.dofIndex() == dofIdxGlobalJ)
490 {
491 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
492 this->gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
493 else
494 this->momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
495 }
496 }
497 }
498 }
499 else
500 DUNE_THROW(Dune::NotImplemented,
501 "Context update for discretization scheme " << MassDiscretizationMethod{}
502 );
503 }
504
505 // \}
506
511 {
512 // use coloring of the momentum discretization for both domains
513 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
514 }
515
522 template<std::size_t i, class AssembleElementFunc>
523 void assembleMultithreaded(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
524 {
525 if (elementSets_.empty())
526 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
527
528 // make this element loop run in parallel
529 // for this we have to color the elements so that we don't get
530 // race conditions when writing into the global matrix
531 // each color can be assembled using multiple threads
532 const auto& grid = this->problem(domainId).gridGeometry().gridView().grid();
533 for (const auto& elements : elementSets_)
534 {
535 Dumux::parallelFor(elements.size(), [&](const std::size_t n)
536 {
537 const auto element = grid.entity(elements[n]);
538 assembleElement(element);
539 });
540 }
541 }
542
543private:
544 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
545 const Element<freeFlowMomentumIndex>& elementI) const
546 {
547 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
548 if (momentumCouplingContext_().empty())
549 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI));
550 else
551 bindCouplingContext_(domainI, elementI, momentumCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
552 }
553
554 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
555 const Element<freeFlowMomentumIndex>& elementI,
556 const std::size_t eIdx) const
557 {
558 if (momentumCouplingContext_().empty())
559 {
560 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
561 fvGeometry.bind(elementI);
562
563 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
564 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
565
566 auto prevElemVolVars = isTransient_ ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
567 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
568
569 if (isTransient_)
570 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[freeFlowMassIndex]);
571
572 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
573 }
574 else if (eIdx != momentumCouplingContext_()[0].eIdx)
575 {
576 momentumCouplingContext_()[0].eIdx = eIdx;
577 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
578 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->curSol(freeFlowMassIndex));
579
580 if (isTransient_)
581 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[freeFlowMassIndex]);
582 }
583 }
584
585 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
586 const Element<freeFlowMassIndex>& elementI) const
587 {
588 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
589 if (massAndEnergyCouplingContext_().empty())
590 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI));
591 else
592 bindCouplingContext_(domainI, elementI, massAndEnergyCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
593 }
594
595 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
596 const Element<freeFlowMassIndex>& elementI,
597 const std::size_t eIdx) const
598 {
599 if (massAndEnergyCouplingContext_().empty())
600 {
601 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
602 auto fvGeometry = localView(gridGeometry);
603 fvGeometry.bindElement(elementI);
604 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
605 }
606 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
607 {
608 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
609 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
610 }
611 }
612
617 template<std::size_t i>
618 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
619 {
620 if (std::get<i>(gridVariables_))
621 return *std::get<i>(gridVariables_);
622 else
623 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
624 }
625
630 template<std::size_t i>
631 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
632 {
633 if (std::get<i>(gridVariables_))
634 return *std::get<i>(gridVariables_);
635 else
636 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
637 }
638
639
640 void computeCouplingStencils_()
641 {
642 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
643 const auto& massGridGeometry = this->problem(freeFlowMassIndex).gridGeometry();
644 auto momentumFvGeometry = localView(momentumGridGeometry);
645 auto massFvGeometry = localView(massGridGeometry);
646
647 massAndEnergyToMomentumStencils_.clear();
648 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
649
650 momentumToMassAndEnergyStencils_.clear();
651 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
652
653 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
654
655 for (const auto& element : elements(momentumGridGeometry.gridView()))
656 {
657 momentumFvGeometry.bindElement(element);
658 massFvGeometry.bindElement(element);
659 const auto eIdx = momentumFvGeometry.elementIndex();
660
661 for (const auto& scv : scvs(momentumFvGeometry))
662 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
663
664 for (const auto& scv : scvs(massFvGeometry))
665 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
666 }
667 }
668
669 CouplingStencilType emptyStencil_;
670 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
671 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
672
673 // the coupling context exists for each thread
674 // TODO this is a bad pattern, just like mutable caches
675 // we should really construct and pass the context and not store it globally
676 std::vector<MomentumCouplingContext>& momentumCouplingContext_() const
677 {
678 thread_local static std::vector<MomentumCouplingContext> c;
679 return c;
680 }
681
682 // the coupling context exists for each thread
683 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_() const
684 {
685 thread_local static std::vector<MassAndEnergyCouplingContext> c;
686 return c;
687 }
688
690 GridVariablesTuple gridVariables_;
691
692 const SolutionVector* prevSol_;
693 bool isTransient_;
694
695 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
696};
697
699template<class T>
701: public std::false_type {};
702
703} // end namespace Dumux
704
705#endif
Coloring schemes for shared-memory-parallel assembly.
Type traits.
Element solution classes and factory functions.
free functions for the evaluation of primary variables inside elements.
The available discretization methods in Dumux.
Parallel for loop (multithreading)
PrimaryVariables evalSolution(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const BoxElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:165
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:172
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: couplingmanager_pq1bubble.hh:448
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:251
constexpr FCDiamond fcdiamond
Definition: method.hh:141
constexpr CCTpfa cctpfa
Definition: method.hh:134
constexpr Box box
Definition: method.hh:136
A vector of primary variables.
Definition: common/properties.hh:49
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
Definition: common/properties.hh:100
The secondary variables within a sub-control volume.
Definition: common/properties.hh:105
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:334
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:299
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:321
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:76
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:350
void updateSolution(const SolutionVector &curSol)
Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want t...
Definition: multidomain/couplingmanager.hh:231
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:83
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_pq1bubble.hh:60
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_pq1bubble.hh:63
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_pq1bubble.hh:510
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_pq1bubble.hh:64
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given sub control volume face.
Definition: couplingmanager_pq1bubble.hh:178
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:413
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:144
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:239
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:377
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:397
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:131
static constexpr auto pressureIdx
Definition: couplingmanager_pq1bubble.hh:123
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:190
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given sub control volume face.
Definition: couplingmanager_pq1bubble.hh:287
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:156
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:329
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_pq1bubble.hh:354
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_pq1bubble.hh:523
trait that is specialized for coupling manager supporting multithreaded assembly
Definition: multidomain/fvassembler.hh:85
A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
Declares all properties used in Dumux.
The interface of the coupling manager for multi domain problems.