3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
couplingmanager_diamond.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_DIAMOND_HH
25#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_DIAMOND_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 {
211 // TODO: cache the shape values when Box method is used
212 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
213 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
214 std::vector<ShapeValue> shapeValues;
215 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
216 localBasis.evaluateFunction(ipLocal, shapeValues);
217
218 Scalar rho = 0.0;
219 for (const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
220 {
221 const auto& volVars = considerPreviousTimeStep ?
222 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
223 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
224 rho += volVars.density()*shapeValues[scv.indexInElement()][0];
225 }
226
227 return rho;
228 }
229 else
230 DUNE_THROW(Dune::NotImplemented,
231 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
232 );
233 }
234
238 Scalar density(const Element<freeFlowMomentumIndex>& element,
239 const SubControlVolume<freeFlowMomentumIndex>& scv,
240 const bool considerPreviousTimeStep = false) const
241 {
242 assert(!(considerPreviousTimeStep && !this->isTransient_));
243 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
244
245 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
246 {
247 const auto eIdx = scv.elementIndex();
248 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
249
250 const auto& volVars = considerPreviousTimeStep ?
251 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
252 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
253
254 return volVars.density();
255 }
256 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box)
257 {
258 // TODO: cache the shape values when Box method is used
259 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
260 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
261 std::vector<ShapeValue> shapeValues;
262 const auto ipLocal = element.geometry().local(scv.dofPosition());
263 localBasis.evaluateFunction(ipLocal, shapeValues);
264
265 Scalar rho = 0.0;
266 for (const auto& scvI : scvs(this->momentumCouplingContext_()[0].fvGeometry))
267 {
268 const auto& volVars = considerPreviousTimeStep ?
269 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
270 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
271 rho += volVars.density()*shapeValues[scvI.indexInElement()][0];
272 }
273 return rho;
274 }
275 else
276 DUNE_THROW(Dune::NotImplemented,
277 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
278 );
279 }
280
284 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
285 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
286 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
287 {
288 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
289
290 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
291 {
292 const auto eIdx = fvGeometry.elementIndex();
293 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
294 const auto& volVars = this->momentumCouplingContext_()[0].curElemVolVars[scv];
295 return volVars.viscosity();
296 }
297 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box)
298 {
299 // TODO: cache the shape values when Box method is used
300 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
301 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
302 std::vector<ShapeValue> shapeValues;
303 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
304 localBasis.evaluateFunction(ipLocal, shapeValues);
305
306 Scalar mu = 0.0;
307 for (const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
308 {
309 const auto& volVars = this->momentumCouplingContext_()[0].curElemVolVars[scv];
310 mu += volVars.viscosity()*shapeValues[scv.indexInElement()][0];
311 }
312
313 return mu;
314 }
315 else
316 DUNE_THROW(Dune::NotImplemented,
317 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
318 );
319
320 }
321
325 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
326 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
327 {
328 // TODO: optimize this function for tpfa where the scvf ip coincides with the dof location
329
330 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(element);
331 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, eIdx);
332
333 const auto& fvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
334 const auto& localBasis = fvGeometry.feLocalBasis();
335 std::vector<ShapeValue> shapeValues;
336 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
337 localBasis.evaluateFunction(ipLocal, shapeValues);
338
339 // interpolate velocity at scvf
340 VelocityVector velocity(0.0);
341 for (const auto& scv : scvs(fvGeometry))
342 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(freeFlowMomentumIndex)[scv.dofIndex()]);
343
344 return velocity;
345 }
346
350 VelocityVector elementVelocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry) const
351 {
352 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), fvGeometry.element());
353
354 const auto& momentumFvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
355 const auto& localBasis = momentumFvGeometry.feLocalBasis();
356
357 // interpolate velocity at scvf
358 VelocityVector velocity(0.0);
359 std::vector<ShapeValue> shapeValues;
360 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
361
362 for (const auto& scv : scvs(momentumFvGeometry))
363 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(freeFlowMomentumIndex)[scv.dofIndex()]);
364
365 return velocity;
366 }
367
372 template<std::size_t j>
373 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
374 const Element<freeFlowMomentumIndex>& elementI,
375 const SubControlVolume<freeFlowMomentumIndex>& scvI,
376 Dune::index_constant<j> domainJ) const
377 { return emptyStencil_; }
378
393 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
394 const Element<freeFlowMassIndex>& elementI,
395 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
396 {
397 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
398 return massAndEnergyToMomentumStencils_[eIdx];
399 }
400
409 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
410 const Element<freeFlowMomentumIndex>& elementI,
411 Dune::index_constant<freeFlowMassIndex> domainJ) const
412 {
413 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
414 return momentumToMassAndEnergyStencils_[eIdx];
415 }
416
417 // \}
418
422 // \{
423
443 template<std::size_t i, std::size_t j, class LocalAssemblerI>
444 void updateCouplingContext(Dune::index_constant<i> domainI,
445 const LocalAssemblerI& localAssemblerI,
446 Dune::index_constant<j> domainJ,
447 std::size_t dofIdxGlobalJ,
448 const PrimaryVariables<j>& priVarsJ,
449 int pvIdxJ)
450 {
451 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
452
453 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
454 {
455 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
456 {
457 bindCouplingContext_(domainI, localAssemblerI.element());
458
459 const auto& problem = this->problem(domainJ);
460 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
461 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
462 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
463 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
464
465 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
466 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
467 else
468 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
469 }
470 }
471 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box)
472 {
473 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
474 {
475 bindCouplingContext_(domainI, localAssemblerI.element());
476
477 const auto& problem = this->problem(domainJ);
478 const auto& deflectedElement = problem.gridGeometry().element(this->momentumCouplingContext_()[0].eIdx);
479 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
480 const auto& fvGeometry = this->momentumCouplingContext_()[0].fvGeometry;
481
482 for (const auto& scv : scvs(fvGeometry))
483 {
484 if(scv.dofIndex() == dofIdxGlobalJ)
485 {
486 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
487 this->gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
488 else
489 this->momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
490 }
491 }
492 }
493 }
494 else
495 DUNE_THROW(Dune::NotImplemented,
496 "Context update for discretization scheme " << MassDiscretizationMethod{}
497 );
498 }
499
500 // \}
501
506 {
507 // use coloring of the mass discretization for both domains
508 // the diamond coloring is a subset (minimum amount of colors) of cctpfa/box coloring
509 elementSets_ = computeColoring(this->problem(freeFlowMassIndex).gridGeometry()).sets;
510 }
511
518 template<std::size_t i, class AssembleElementFunc>
519 void assembleMultithreaded(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
520 {
521 if (elementSets_.empty())
522 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
523
524 // make this element loop run in parallel
525 // for this we have to color the elements so that we don't get
526 // race conditions when writing into the global matrix
527 // each color can be assembled using multiple threads
528 const auto& grid = this->problem(freeFlowMassIndex).gridGeometry().gridView().grid();
529 for (const auto& elements : elementSets_)
530 {
531 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
532 {
533 const auto element = grid.entity(elements[eIdx]);
534 assembleElement(element);
535 });
536 }
537 }
538
539private:
540 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
541 const Element<freeFlowMomentumIndex>& elementI) const
542 {
543 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
544 if (momentumCouplingContext_().empty())
545 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI));
546 else
547 bindCouplingContext_(domainI, elementI, momentumCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
548 }
549
550 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
551 const Element<freeFlowMomentumIndex>& elementI,
552 const std::size_t eIdx) const
553 {
554 if (momentumCouplingContext_().empty())
555 {
556 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
557 fvGeometry.bind(elementI);
558
559 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
560 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
561
562 auto prevElemVolVars = isTransient_ ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
563 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
564
565 if (isTransient_)
566 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[freeFlowMassIndex]);
567
568 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
569 }
570 else if (eIdx != momentumCouplingContext_()[0].eIdx)
571 {
572 momentumCouplingContext_()[0].eIdx = eIdx;
573 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
574 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->curSol(freeFlowMassIndex));
575
576 if (isTransient_)
577 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[freeFlowMassIndex]);
578 }
579 }
580
581 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
582 const Element<freeFlowMassIndex>& elementI) const
583 {
584 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
585 if (massAndEnergyCouplingContext_().empty())
586 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI));
587 else
588 bindCouplingContext_(domainI, elementI, massAndEnergyCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
589 }
590
591 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
592 const Element<freeFlowMassIndex>& elementI,
593 const std::size_t eIdx) const
594 {
595 if (massAndEnergyCouplingContext_().empty())
596 {
597 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
598 auto fvGeometry = localView(gridGeometry);
599 fvGeometry.bindElement(elementI);
600 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
601 }
602 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
603 {
604 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
605 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
606 }
607 }
608
613 template<std::size_t i>
614 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
615 {
616 if (std::get<i>(gridVariables_))
617 return *std::get<i>(gridVariables_);
618 else
619 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
620 }
621
626 template<std::size_t i>
627 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
628 {
629 if (std::get<i>(gridVariables_))
630 return *std::get<i>(gridVariables_);
631 else
632 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
633 }
634
635
636 void computeCouplingStencils_()
637 {
638 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
639 const auto& massGridGeometry = this->problem(freeFlowMassIndex).gridGeometry();
640 auto momentumFvGeometry = localView(momentumGridGeometry);
641 auto massFvGeometry = localView(massGridGeometry);
642
643 massAndEnergyToMomentumStencils_.clear();
644 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
645
646 momentumToMassAndEnergyStencils_.clear();
647 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
648
649 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
650
651 for (const auto& element : elements(momentumGridGeometry.gridView()))
652 {
653 momentumFvGeometry.bindElement(element);
654 massFvGeometry.bindElement(element);
655 const auto eIdx = momentumFvGeometry.elementIndex();
656
657 for (const auto& scv : scvs(momentumFvGeometry))
658 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
659
660 for (const auto& scv : scvs(massFvGeometry))
661 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
662 }
663 }
664
665 CouplingStencilType emptyStencil_;
666 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
667 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
668
669 // the coupling context exists for each thread
670 // TODO this is a bad pattern, just like mutable caches
671 // we should really construct and pass the context and not store it globally
672 std::vector<MomentumCouplingContext>& momentumCouplingContext_() const
673 {
674 thread_local static std::vector<MomentumCouplingContext> c;
675 return c;
676 }
677
678 // the coupling context exists for each thread
679 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_() const
680 {
681 thread_local static std::vector<MassAndEnergyCouplingContext> c;
682 return c;
683 }
684
686 GridVariablesTuple gridVariables_;
687
688 const SolutionVector* prevSol_;
689 bool isTransient_;
690
691 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
692};
693
695template<class T>
697: public std::true_type {};
698
699} // end namespace Dumux
700
701#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_diamond.hh:444
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 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_diamond.hh:60
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_diamond.hh:64
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_diamond.hh:350
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_diamond.hh:284
Scalar density(const Element< freeFlowMomentumIndex > &element, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: couplingmanager_diamond.hh:238
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_diamond.hh:156
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_diamond.hh:190
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_diamond.hh:505
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_diamond.hh:409
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: couplingmanager_diamond.hh:325
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_diamond.hh:373
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_diamond.hh:144
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_diamond.hh:63
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_diamond.hh:519
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_diamond.hh:178
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_diamond.hh:131
static constexpr auto pressureIdx
Definition: couplingmanager_diamond.hh:123
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_diamond.hh:393
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.