version 3.9
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_STAGGERED_HH
13#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_STAGGERED_HH
14
15#include <memory>
16#include <tuple>
17#include <vector>
18#include <deque>
19
20#include <dune/common/exceptions.hh>
21#include <dune/common/indices.hh>
22#include <dune/common/float_cmp.hh>
23
26
30
34
37
38namespace Dumux {
39
45template<class Traits>
47: public CouplingManager<Traits>
48{
50public:
51 static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index();
52 static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index();
53
54 // this can be used if the coupling manager is used inside a meta-coupling manager (e.g. multi-binary)
55 // to manager the solution vector storage outside this class
57private:
58 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
59 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
60 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
61 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
62 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
63 template<std::size_t id> using ElementSeed = typename GridView<id>::Grid::template Codim<0>::EntitySeed;
64 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
65 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
66 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
67 template<std::size_t id> using GridVariables = typename Traits::template SubDomain<id>::GridVariables;
68 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
69 template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache;
70 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
71 template<std::size_t id> using VolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::VolumeVariables>;
72
73 using Scalar = typename Traits::Scalar;
74 using SolutionVector = typename Traits::SolutionVector;
75
76 using CouplingStencilType = std::vector<std::size_t>;
77
78 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
79
80 using FluidSystem = typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
81
82 using VelocityVector = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
83 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
84
85 struct MomentumCouplingContext
86 {
87 FVElementGeometry<freeFlowMassIndex> fvGeometry;
88 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
89 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
90 std::size_t eIdx;
91 };
92
93 struct MassAndEnergyCouplingContext
94 {
95 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f, const std::size_t i)
96 : fvGeometry(std::move(f))
97 , eIdx(i)
98 {}
99
100 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
101 std::size_t eIdx;
102 };
103
104public:
105
106 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
107
111 // \{
112
114 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
115 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
116 GridVariablesTuple&& gridVariables,
117 const SolutionVector& curSol)
118 {
119 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
120 gridVariables_ = gridVariables;
121 this->updateSolution(curSol);
122
123 computeCouplingStencils_();
124 }
125
127 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
128 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
129 GridVariablesTuple&& gridVariables,
130 const SolutionVector& curSol,
131 const SolutionVector& prevSol)
132 {
133 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
134 prevSol_ = &prevSol;
135 isTransient_ = true;
136 }
137
139 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
140 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
141 GridVariablesTuple&& gridVariables,
143 {
144 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
145 gridVariables_ = gridVariables;
146 this->attachSolution(curSol);
147
148 computeCouplingStencils_();
149 }
150
151 // \}
152
154
172 template<std::size_t j, class LocalAssemblerI>
173 decltype(auto) evalCouplingResidual(Dune::index_constant<freeFlowMomentumIndex> domainI,
174 const LocalAssemblerI& localAssemblerI,
175 const SubControlVolume<freeFlowMomentumIndex>& scvI,
176 Dune::index_constant<j> domainJ,
177 std::size_t dofIdxGlobalJ) const
178 {
179 const auto& problem = localAssemblerI.problem();
180 const auto& element = localAssemblerI.element();
181 const auto& fvGeometry = localAssemblerI.fvGeometry();
182 const auto& curElemVolVars = localAssemblerI.curElemVolVars();
183 const auto& prevElemVolVars = localAssemblerI.prevElemVolVars();
184 typename LocalAssemblerI::ElementResidualVector residual(localAssemblerI.element().subEntities(1));
185 const auto& localResidual = localAssemblerI.localResidual();
186
187 localResidual.evalSource(residual, problem, element, fvGeometry, curElemVolVars, scvI);
188
189 for (const auto& scvf : scvfs(fvGeometry, scvI))
190 localResidual.evalFlux(residual, problem, element, fvGeometry, curElemVolVars, localAssemblerI.elemBcTypes(), localAssemblerI.elemFluxVarsCache(), scvf);
191
192 if (!localAssemblerI.assembler().isStationaryProblem())
193 {
194 assert(isTransient_);
195 localResidual.evalStorage(residual, problem, element, fvGeometry, prevElemVolVars, curElemVolVars, scvI);
196 }
197
198 return residual;
199 }
200
204 // \{
205
209 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
210 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
211 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
212 {
213 assert(scvf.isFrontal() && !scvf.isLateral() && !scvf.boundary());
214 return this->curSol(freeFlowMassIndex)[fvGeometry.elementIndex()][pressureIdx];
215 }
216
223 Scalar cellPressure(const Element<freeFlowMassIndex>& element,
224 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
225 {
226 return this->curSol(freeFlowMassIndex)[scvf.insideScvIdx()][pressureIdx];
227 }
228
232 Scalar density(const Element<freeFlowMomentumIndex>& element,
233 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
234 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
235 const bool considerPreviousTimeStep = false) const
236 {
237 assert(!(considerPreviousTimeStep && !isTransient_));
238 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
239 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
240 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
241
242 const auto rho = [&](const auto& elemVolVars)
243 {
244 if (scvf.boundary())
245 return elemVolVars[insideMassScv].density();
246 else
247 {
248 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
249 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
250 // TODO distance weighting
251 return 0.5*(elemVolVars[insideMassScv].density() + elemVolVars[outsideMassScv].density());
252 }
253 };
254
255 return considerPreviousTimeStep ? rho(momentumCouplingContext_()[0].prevElemVolVars)
256 : rho(momentumCouplingContext_()[0].curElemVolVars);
257 }
258
259 auto insideAndOutsideDensity(const Element<freeFlowMomentumIndex>& element,
260 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
261 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
262 const bool considerPreviousTimeStep = false) const
263 {
264 assert(!(considerPreviousTimeStep && !isTransient_));
265 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
266 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
267 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
268
269 const auto result = [&](const auto& elemVolVars)
270 {
271 if (scvf.boundary())
272 return std::make_pair(elemVolVars[insideMassScv].density(), elemVolVars[insideMassScv].density());
273 else
274 {
275 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
276 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
277 return std::make_pair(elemVolVars[insideMassScv].density(), elemVolVars[outsideMassScv].density());
278 }
279 };
280
281 return considerPreviousTimeStep ? result(momentumCouplingContext_()[0].prevElemVolVars)
282 : result(momentumCouplingContext_()[0].curElemVolVars);
283 }
284
288 Scalar density(const Element<freeFlowMomentumIndex>& element,
289 const SubControlVolume<freeFlowMomentumIndex>& scv,
290 const bool considerPreviousTimeStep = false) const
291 {
292 assert(!(considerPreviousTimeStep && !isTransient_));
293 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
294 const auto& massScv = (*scvs(momentumCouplingContext_()[0].fvGeometry).begin());
295
296 return considerPreviousTimeStep ? momentumCouplingContext_()[0].prevElemVolVars[massScv].density()
297 : momentumCouplingContext_()[0].curElemVolVars[massScv].density();
298 }
299
303 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
304 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
305 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
306 {
307 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
308
309 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
310 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
311
312 if (scvf.boundary())
313 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
314
315 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
316 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
317
318 const auto mu = [&](const auto& elemVolVars)
319 {
320 // TODO distance weighting
321 return 0.5*(elemVolVars[insideMassScv].viscosity() + elemVolVars[outsideMassScv].viscosity());
322 };
323
324 return mu(momentumCouplingContext_()[0].curElemVolVars);
325 }
326
330 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
331 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
332 const SubControlVolume<freeFlowMomentumIndex>& scv) const
333 {
334 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
335 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(scv.elementIndex());
336 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
337 }
338
342 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
343 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
344 {
345 // TODO: rethink this! Maybe we only need scvJ.dofIndex()
346 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, scvf.insideScvIdx()/*eIdx*/);
347
348 // the TPFA scvf index corresponds to the staggered scv index (might need mapping)
349 const auto localMomentumScvIdx = massScvfToMomentumScvIdx_(scvf, massAndEnergyCouplingContext_()[0].fvGeometry);
350 const auto& scvJ = massAndEnergyCouplingContext_()[0].fvGeometry.scv(localMomentumScvIdx);
351
352 // create a unit normal vector oriented in positive coordinate direction
353 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition velocity;
354 velocity[scvJ.dofAxis()] = 1.0;
355
356 // create the actual velocity vector
357 velocity *= this->curSol(freeFlowMomentumIndex)[scvJ.dofIndex()];
358
359 return velocity;
360 }
361
371 template<std::size_t j>
372 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
373 const Element<freeFlowMomentumIndex>& elementI,
374 const SubControlVolume<freeFlowMomentumIndex>& scvI,
375 Dune::index_constant<j> domainJ) const
376 { return emptyStencil_; }
377
392 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
393 const Element<freeFlowMassIndex>& elementI,
394 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
395 {
396 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
397 return massAndEnergyToMomentumStencils_[eIdx];
398 }
399
409 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
410 const Element<freeFlowMomentumIndex>& elementI,
411 const SubControlVolume<freeFlowMomentumIndex>& scvI,
412 Dune::index_constant<freeFlowMassIndex> domainJ) const
413 {
414 return momentumToMassAndEnergyStencils_[scvI.index()];
415 }
416
417 // \}
418
422 // \{
423
425 template<std::size_t i, std::size_t j, class LocalAssemblerI>
426 void updateCouplingContext(Dune::index_constant<i> domainI,
427 const LocalAssemblerI& localAssemblerI,
428 Dune::index_constant<j> domainJ,
429 std::size_t dofIdxGlobalJ,
430 const PrimaryVariables<j>& priVarsJ,
431 int pvIdxJ)
432 {
433 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
434
435 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
436 {
437 bindCouplingContext_(domainI, localAssemblerI.element());
438
439 const auto& problem = this->problem(domainJ);
440 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
441 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
442 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
443 const auto scvIdxJ = dofIdxGlobalJ;
444 const auto& scv = fvGeometry.scv(scvIdxJ);
445
446 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
447 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
448 else
449 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
450 }
451 }
452
453 // \}
454
459 {
460 // use coloring of the fc staggered discretization for both domains
461 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
462 }
463
470 template<std::size_t i, class AssembleElementFunc>
471 void assembleMultithreaded(Dune::index_constant<i> domainI, AssembleElementFunc&& assembleElement) const
472 {
473 if (elementSets_.empty())
474 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
475
476 // make this element loop run in parallel
477 // for this we have to color the elements so that we don't get
478 // race conditions when writing into the global matrix
479 // each color can be assembled using multiple threads
480 const auto& grid = this->problem(freeFlowMomentumIndex).gridGeometry().gridView().grid();
481 for (const auto& elements : elementSets_)
482 {
483 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
484 {
485 const auto element = grid.entity(elements[eIdx]);
486 assembleElement(element);
487 });
488 }
489 }
490
491private:
492 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
493 const Element<freeFlowMomentumIndex>& elementI) const
494 {
495 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
496 bindCouplingContext_(domainI, elementI, eIdx);
497 }
498
499 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
500 const Element<freeFlowMomentumIndex>& elementI,
501 const std::size_t eIdx) const
502 {
503 if (momentumCouplingContext_().empty())
504 {
505 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
506 fvGeometry.bind(elementI);
507
508 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
509 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
510
511 auto prevElemVolVars = isTransient_ ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
512 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
513
514 if (isTransient_)
515 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[freeFlowMassIndex]);
516
517 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
518 }
519 else if (eIdx != momentumCouplingContext_()[0].eIdx)
520 {
521 momentumCouplingContext_()[0].eIdx = eIdx;
522 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
523 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->curSol(freeFlowMassIndex));
524
525 if (isTransient_)
526 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[freeFlowMassIndex]);
527 }
528 }
529
530 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
531 const Element<freeFlowMassIndex>& elementI) const
532 {
533 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
534 bindCouplingContext_(domainI, elementI, eIdx);
535 }
536
537 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
538 const Element<freeFlowMassIndex>& elementI,
539 const std::size_t eIdx) const
540 {
541 if (massAndEnergyCouplingContext_().empty())
542 {
543 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
544 auto fvGeometry = localView(gridGeometry);
545 fvGeometry.bindElement(elementI);
546 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
547 }
548 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
549 {
550 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
551 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
552 }
553 }
554
559 template<std::size_t i>
560 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
561 {
562 if (std::get<i>(gridVariables_))
563 return *std::get<i>(gridVariables_);
564 else
565 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
566 }
567
572 template<std::size_t i>
573 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
574 {
575 if (std::get<i>(gridVariables_))
576 return *std::get<i>(gridVariables_);
577 else
578 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
579 }
580
581
582 void computeCouplingStencils_()
583 {
584 // TODO higher order
585 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
586 auto momentumFvGeometry = localView(momentumGridGeometry);
587 massAndEnergyToMomentumStencils_.clear();
588 massAndEnergyToMomentumStencils_.resize(momentumGridGeometry.gridView().size(0));
589
590 momentumToMassAndEnergyStencils_.clear();
591 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.numScv());
592
593 for (const auto& element : elements(momentumGridGeometry.gridView()))
594 {
595 const auto eIdx = momentumGridGeometry.elementMapper().index(element);
596 momentumFvGeometry.bind(element);
597 for (const auto& scv : scvs(momentumFvGeometry))
598 {
599 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
600 momentumToMassAndEnergyStencils_[scv.index()].push_back(eIdx);
601
602 // extend the stencil for fluids with variable viscosity and density,
603 if constexpr (FluidSystem::isCompressible(0/*phaseIdx*/))
604 // if constexpr (FluidSystem::isCompressible(0/*phaseIdx*/) || !FluidSystem::viscosityIsConstant(0/*phaseIdx*/)) // TODO fix on master
605 {
606 for (const auto& scvf : scvfs(momentumFvGeometry, scv))
607 {
608 if (scvf.isLateral() && !scvf.boundary())
609 {
610 const auto& outsideScv = momentumFvGeometry.scv(scvf.outsideScvIdx());
611 momentumToMassAndEnergyStencils_[scv.index()].push_back(outsideScv.elementIndex());
612 }
613 }
614 }
615 }
616 }
617 }
618
619 std::size_t massScvfToMomentumScvIdx_(const SubControlVolumeFace<freeFlowMassIndex>& massScvf,
620 [[maybe_unused]] const FVElementGeometry<freeFlowMomentumIndex>& momentumFVGeometry) const
621 {
622 if constexpr (ConsistentlyOrientedGrid<typename GridView<freeFlowMomentumIndex>::Grid>{})
623 return massScvf.index();
624 else
625 {
626 static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true);
627 if (!makeConsistentlyOriented)
628 return massScvf.index();
629
630 for (const auto& momentumScv : scvs(momentumFVGeometry))
631 {
632 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition momentumUnitOuterNormal(0.0);
633 momentumUnitOuterNormal[momentumScv.dofAxis()] = momentumScv.directionSign();
634 if (Dune::FloatCmp::eq<typename GridView<freeFlowMomentumIndex>::ctype>(massScvf.unitOuterNormal()*momentumUnitOuterNormal, 1.0))
635 return momentumScv.index();
636 }
637 DUNE_THROW(Dune::InvalidStateException, "No Momentum SCV found");
638 }
639 }
640
641 CouplingStencilType emptyStencil_;
642 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
643 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
644
645 // the coupling context exists for each thread
646 // TODO this is a bad pattern, just like mutable caches
647 // we should really construct and pass the context and not store it globally
648 std::vector<MomentumCouplingContext>& momentumCouplingContext_() const
649 {
650 thread_local static std::vector<MomentumCouplingContext> c;
651 return c;
652 }
653
654 // the coupling context exists for each thread
655 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_() const
656 {
657 thread_local static std::vector<MassAndEnergyCouplingContext> c;
658 return c;
659 }
660
662 GridVariablesTuple gridVariables_;
663
664 const SolutionVector* prevSol_;
665 bool isTransient_;
666
667 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
668};
669
671template<class T>
673: public std::false_type {};
674
675} // end namespace Dumux
676
677#endif
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:322
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:287
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:309
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:64
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:338
void updateSolution(const SolutionVector &curSol)
Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want t...
Definition: multidomain/couplingmanager.hh:219
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:71
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_staggered.hh:48
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:209
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_staggered.hh:51
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:303
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:426
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:173
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_staggered.hh:52
void assembleMultithreaded(Dune::index_constant< i > domainI, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_staggered.hh:471
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:288
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:114
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 scv's residual depends ...
Definition: couplingmanager_staggered.hh:372
auto insideAndOutsideDensity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Definition: couplingmanager_staggered.hh:259
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:409
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:127
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv) const
Returns the pressure at a given sub control volume.
Definition: couplingmanager_staggered.hh:330
static constexpr auto pressureIdx
Definition: couplingmanager_staggered.hh:106
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:392
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:232
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:139
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:223
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:342
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_staggered.hh:458
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Type traits.
Element solution classes and factory functions.
free functions for the evaluation of primary variables inside elements.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
The interface of the coupling manager for multi domain problems.
A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
Definition: adapt.hh:17
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:239
Parallel for loop (multithreading)
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78