3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
couplingmanager_staggered.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_STAGGERED_HH
25#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_STAGGERED_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
38
42
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 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
96
97 struct MomentumCouplingContext
98 {
99 FVElementGeometry<freeFlowMassIndex> fvGeometry;
100 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
101 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
102 std::size_t eIdx;
103 };
104
105 struct MassAndEnergyCouplingContext
106 {
107 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f, const std::size_t i)
108 : fvGeometry(std::move(f))
109 , eIdx(i)
110 {}
111
112 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
113 std::size_t eIdx;
114 };
115
116public:
117
118 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
119
123 // \{
124
126 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
127 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
128 GridVariablesTuple&& gridVariables,
129 const SolutionVector& curSol)
130 {
131 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
132 gridVariables_ = gridVariables;
133 this->updateSolution(curSol);
134
135 computeCouplingStencils_();
136 }
137
139 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
140 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
141 GridVariablesTuple&& gridVariables,
142 const SolutionVector& curSol,
143 const SolutionVector& prevSol)
144 {
145 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
146 prevSol_ = &prevSol;
147 isTransient_ = true;
148 }
149
151 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
152 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
153 GridVariablesTuple&& gridVariables,
155 {
156 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
157 gridVariables_ = gridVariables;
158 this->attachSolution(curSol);
159
160 computeCouplingStencils_();
161 }
162
163 // \}
164
166
185 template<std::size_t j, class LocalAssemblerI>
186 decltype(auto) evalCouplingResidual(Dune::index_constant<freeFlowMomentumIndex> domainI,
187 const LocalAssemblerI& localAssemblerI,
188 const SubControlVolume<freeFlowMomentumIndex>& scvI,
189 Dune::index_constant<j> domainJ,
190 std::size_t dofIdxGlobalJ) const
191 {
192 const auto& problem = localAssemblerI.problem();
193 const auto& element = localAssemblerI.element();
194 const auto& fvGeometry = localAssemblerI.fvGeometry();
195 const auto& curElemVolVars = localAssemblerI.curElemVolVars();
196 const auto& prevElemVolVars = localAssemblerI.prevElemVolVars();
197 typename LocalAssemblerI::ElementResidualVector residual(localAssemblerI.element().subEntities(1));
198 const auto& localResidual = localAssemblerI.localResidual();
199
200 localResidual.evalSource(residual, problem, element, fvGeometry, curElemVolVars, scvI);
201
202 for (const auto& scvf : scvfs(fvGeometry, scvI))
203 localResidual.evalFlux(residual, problem, element, fvGeometry, curElemVolVars, localAssemblerI.elemBcTypes(), localAssemblerI.elemFluxVarsCache(), scvf);
204
205 if (!localAssemblerI.assembler().isStationaryProblem())
206 {
207 assert(isTransient_);
208 localResidual.evalStorage(residual, problem, element, fvGeometry, prevElemVolVars, curElemVolVars, scvI);
209 }
210
211 return residual;
212 }
213
217 // \{
218
222 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
223 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
224 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
225 {
226 assert(scvf.isFrontal() && !scvf.isLateral() && !scvf.boundary());
227 return this->curSol(freeFlowMassIndex)[fvGeometry.elementIndex()][pressureIdx];
228 }
229
236 Scalar cellPressure(const Element<freeFlowMassIndex>& element,
237 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
238 {
239 return this->curSol(freeFlowMassIndex)[scvf.insideScvIdx()][pressureIdx];
240 }
241
245 Scalar density(const Element<freeFlowMomentumIndex>& element,
246 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
247 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
248 const bool considerPreviousTimeStep = false) const
249 {
250 assert(!(considerPreviousTimeStep && !isTransient_));
251 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
252 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
253 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
254
255 const auto rho = [&](const auto& elemVolVars)
256 {
257 if (scvf.boundary())
258 return elemVolVars[insideMassScv].density();
259 else
260 {
261 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
262 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
263 // TODO distance weighting
264 return 0.5*(elemVolVars[insideMassScv].density() + elemVolVars[outsideMassScv].density());
265 }
266 };
267
268 return considerPreviousTimeStep ? rho(momentumCouplingContext_()[0].prevElemVolVars)
269 : rho(momentumCouplingContext_()[0].curElemVolVars);
270 }
271
272 auto insideAndOutsideDensity(const Element<freeFlowMomentumIndex>& element,
273 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
274 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
275 const bool considerPreviousTimeStep = false) const
276 {
277 assert(!(considerPreviousTimeStep && !isTransient_));
278 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
279 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
280 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
281
282 const auto result = [&](const auto& elemVolVars)
283 {
284 if (scvf.boundary())
285 return std::make_pair(elemVolVars[insideMassScv].density(), elemVolVars[insideMassScv].density());
286 else
287 {
288 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
289 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
290 return std::make_pair(elemVolVars[insideMassScv].density(), elemVolVars[outsideMassScv].density());
291 }
292 };
293
294 return considerPreviousTimeStep ? result(momentumCouplingContext_()[0].prevElemVolVars)
295 : result(momentumCouplingContext_()[0].curElemVolVars);
296 }
297
301 Scalar density(const Element<freeFlowMomentumIndex>& element,
302 const SubControlVolume<freeFlowMomentumIndex>& scv,
303 const bool considerPreviousTimeStep = false) const
304 {
305 assert(!(considerPreviousTimeStep && !isTransient_));
306 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
307 const auto& massScv = (*scvs(momentumCouplingContext_()[0].fvGeometry).begin());
308
309 return considerPreviousTimeStep ? momentumCouplingContext_()[0].prevElemVolVars[massScv].density()
310 : momentumCouplingContext_()[0].curElemVolVars[massScv].density();
311 }
312
316 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
317 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
318 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
319 {
320 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
321
322 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
323 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
324
325 if (scvf.boundary())
326 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
327
328 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
329 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
330
331 const auto mu = [&](const auto& elemVolVars)
332 {
333 // TODO distance weighting
334 return 0.5*(elemVolVars[insideMassScv].viscosity() + elemVolVars[outsideMassScv].viscosity());
335 };
336
337 return mu(momentumCouplingContext_()[0].curElemVolVars);
338 }
339
343 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
344 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
345 {
346 // TODO: rethink this! Maybe we only need scvJ.dofIndex()
347 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, scvf.insideScvIdx()/*eIdx*/);
348
349 // the TPFA scvf index corresponds to the staggered scv index (might need mapping)
350 const auto localMomentumScvIdx = massScvfToMomentumScvIdx_(scvf, massAndEnergyCouplingContext_()[0].fvGeometry);
351 const auto& scvJ = massAndEnergyCouplingContext_()[0].fvGeometry.scv(localMomentumScvIdx);
352
353 // create a unit normal vector oriented in positive coordinate direction
354 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition velocity;
355 velocity[scvJ.dofAxis()] = 1.0;
356
357 // create the actual velocity vector
358 velocity *= this->curSol(freeFlowMomentumIndex)[scvJ.dofIndex()];
359
360 return velocity;
361 }
362
367 template<std::size_t j>
368 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
369 const Element<freeFlowMomentumIndex>& elementI,
370 const SubControlVolume<freeFlowMomentumIndex>& scvI,
371 Dune::index_constant<j> domainJ) const
372 { return emptyStencil_; }
373
388 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
389 const Element<freeFlowMassIndex>& elementI,
390 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
391 {
392 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
393 return massAndEnergyToMomentumStencils_[eIdx];
394 }
395
405 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
406 const Element<freeFlowMomentumIndex>& elementI,
407 const SubControlVolume<freeFlowMomentumIndex>& scvI,
408 Dune::index_constant<freeFlowMassIndex> domainJ) const
409 {
410 return momentumToMassAndEnergyStencils_[scvI.index()];
411 }
412
413 // \}
414
418 // \{
419
438 template<std::size_t i, std::size_t j, class LocalAssemblerI>
439 void updateCouplingContext(Dune::index_constant<i> domainI,
440 const LocalAssemblerI& localAssemblerI,
441 Dune::index_constant<j> domainJ,
442 std::size_t dofIdxGlobalJ,
443 const PrimaryVariables<j>& priVarsJ,
444 int pvIdxJ)
445 {
446 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
447
448 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
449 {
450 bindCouplingContext_(domainI, localAssemblerI.element());
451
452 const auto& problem = this->problem(domainJ);
453 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
454 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
455 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
456 const auto scvIdxJ = dofIdxGlobalJ;
457 const auto& scv = fvGeometry.scv(scvIdxJ);
458
459 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
460 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
461 else
462 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
463 }
464 }
465
466 // \}
467
475 {
476 // use coloring of the fc staggered discretization for both domains
477 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
478 }
479
486 template<std::size_t i, class AssembleElementFunc>
487 void assembleMultithreaded(Dune::index_constant<i> domainI, AssembleElementFunc&& assembleElement) const
488 {
489 if (elementSets_.empty())
490 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
491
492 // make this element loop run in parallel
493 // for this we have to color the elements so that we don't get
494 // race conditions when writing into the global matrix
495 // each color can be assembled using multiple threads
496 const auto& grid = this->problem(freeFlowMomentumIndex).gridGeometry().gridView().grid();
497 for (const auto& elements : elementSets_)
498 {
499 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
500 {
501 const auto element = grid.entity(elements[eIdx]);
502 assembleElement(element);
503 });
504 }
505 }
506
507private:
508 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
509 const Element<freeFlowMomentumIndex>& elementI) const
510 {
511 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
512 bindCouplingContext_(domainI, elementI, eIdx);
513 }
514
515 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
516 const Element<freeFlowMomentumIndex>& elementI,
517 const std::size_t eIdx) const
518 {
519 if (momentumCouplingContext_().empty())
520 {
521 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
522 fvGeometry.bind(elementI);
523
524 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
525 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
526
527 auto prevElemVolVars = isTransient_ ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
528 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
529
530 if (isTransient_)
531 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[freeFlowMassIndex]);
532
533 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
534 }
535 else if (eIdx != momentumCouplingContext_()[0].eIdx)
536 {
537 momentumCouplingContext_()[0].eIdx = eIdx;
538 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
539 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->curSol(freeFlowMassIndex));
540
541 if (isTransient_)
542 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[freeFlowMassIndex]);
543 }
544 }
545
546 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
547 const Element<freeFlowMassIndex>& elementI) const
548 {
549 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
550 bindCouplingContext_(domainI, elementI, eIdx);
551 }
552
553 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
554 const Element<freeFlowMassIndex>& elementI,
555 const std::size_t eIdx) const
556 {
557 if (massAndEnergyCouplingContext_().empty())
558 {
559 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
560 auto fvGeometry = localView(gridGeometry);
561 fvGeometry.bindElement(elementI);
562 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
563 }
564 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
565 {
566 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
567 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
568 }
569 }
570
575 template<std::size_t i>
576 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
577 {
578 if (std::get<i>(gridVariables_))
579 return *std::get<i>(gridVariables_);
580 else
581 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
582 }
583
588 template<std::size_t i>
589 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
590 {
591 if (std::get<i>(gridVariables_))
592 return *std::get<i>(gridVariables_);
593 else
594 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
595 }
596
597
598 void computeCouplingStencils_()
599 {
600 // TODO higher order
601 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
602 auto momentumFvGeometry = localView(momentumGridGeometry);
603 massAndEnergyToMomentumStencils_.clear();
604 massAndEnergyToMomentumStencils_.resize(momentumGridGeometry.gridView().size(0));
605
606 momentumToMassAndEnergyStencils_.clear();
607 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.numScv());
608
609 for (const auto& element : elements(momentumGridGeometry.gridView()))
610 {
611 const auto eIdx = momentumGridGeometry.elementMapper().index(element);
612 momentumFvGeometry.bind(element);
613 for (const auto& scv : scvs(momentumFvGeometry))
614 {
615 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
616 momentumToMassAndEnergyStencils_[scv.index()].push_back(eIdx);
617
618 // extend the stencil for fluids with variable viscosity and density,
619 if constexpr (FluidSystem::isCompressible(0/*phaseIdx*/))
620 // if constexpr (FluidSystem::isCompressible(0/*phaseIdx*/) || !FluidSystem::viscosityIsConstant(0/*phaseIdx*/)) // TODO fix on master
621 {
622 for (const auto& scvf : scvfs(momentumFvGeometry, scv))
623 {
624 if (scvf.isLateral() && !scvf.boundary())
625 {
626 const auto& outsideScv = momentumFvGeometry.scv(scvf.outsideScvIdx());
627 momentumToMassAndEnergyStencils_[scv.index()].push_back(outsideScv.elementIndex());
628 }
629 }
630 }
631 }
632 }
633 }
634
635 std::size_t massScvfToMomentumScvIdx_(const SubControlVolumeFace<freeFlowMassIndex>& massScvf,
636 [[maybe_unused]] const FVElementGeometry<freeFlowMomentumIndex>& momentumFVGeometry) const
637 {
638 if constexpr (ConsistentlyOrientedGrid<typename GridView<freeFlowMomentumIndex>::Grid>{})
639 return massScvf.index();
640 else
641 {
642 static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true);
643 if (!makeConsistentlyOriented)
644 return massScvf.index();
645
646 for (const auto& momentumScv : scvs(momentumFVGeometry))
647 {
648 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition momentumUnitOuterNormal(0.0);
649 momentumUnitOuterNormal[momentumScv.dofAxis()] = momentumScv.directionSign();
650 if (Dune::FloatCmp::eq<typename GridView<freeFlowMomentumIndex>::ctype>(massScvf.unitOuterNormal()*momentumUnitOuterNormal, 1.0))
651 return momentumScv.index();
652 }
653 DUNE_THROW(Dune::InvalidStateException, "No Momentum SCV found");
654 }
655 }
656
657 CouplingStencilType emptyStencil_;
658 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
659 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
660
661 // the coupling context exists for each thread
662 // TODO this is a bad pattern, just like mutable caches
663 // we should really construct and pass the context and not store it globally
664 std::vector<MomentumCouplingContext>& momentumCouplingContext_() const
665 {
666 thread_local static std::vector<MomentumCouplingContext> c;
667 return c;
668 }
669
670 // the coupling context exists for each thread
671 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_() const
672 {
673 thread_local static std::vector<MassAndEnergyCouplingContext> c;
674 return c;
675 }
676
678 GridVariablesTuple gridVariables_;
679
680 const SolutionVector* prevSol_;
681 bool isTransient_;
682
683 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
684};
685
687template<class T>
689: public std::true_type {};
690
691} // end namespace Dumux
692
693#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)
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
decltype(auto) evalCouplingResidual(Dune::index_constant< freeFlowMomentumIndex > domainI, const LocalAssemblerI &localAssemblerI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
evaluates the element residual of a coupled element of domain i which depends on the variables at the...
Definition: couplingmanager_staggered.hh:186
decltype(auto) evalCouplingResidual(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
evaluates the element residual of a coupled element of domain i which depends on the variables at the...
Definition: multidomain/couplingmanager.hh:261
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
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_staggered.hh:60
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given frontal sub control volume face.
Definition: couplingmanager_staggered.hh:222
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_staggered.hh:63
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_staggered.hh:316
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_staggered.hh:439
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_staggered.hh:64
void assembleMultithreaded(Dune::index_constant< i > domainI, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_staggered.hh:487
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_staggered.hh:301
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_staggered.hh:126
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_staggered.hh:368
auto insideAndOutsideDensity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Definition: couplingmanager_staggered.hh:272
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, 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_staggered.hh:405
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_staggered.hh:139
static constexpr auto pressureIdx
Definition: couplingmanager_staggered.hh:118
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_staggered.hh:388
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_staggered.hh:245
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_staggered.hh:151
Scalar cellPressure(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the pressure at the center of a sub control volume corresponding to a given sub control volum...
Definition: couplingmanager_staggered.hh:236
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: couplingmanager_staggered.hh:343
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_staggered.hh:474
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.