version 3.10.0
Loading...
Searching...
No Matches
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-FileCopyrightText: 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{
49 using ParentType = CouplingManager<Traits>;
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 template<std::size_t id>
77 using SubSolutionVector
78 = std::decay_t<decltype(std::declval<SolutionVector>()[Dune::index_constant<id>()])>;
79
80 template<std::size_t id>
81 using ConstSubSolutionVectorPtr = const SubSolutionVector<id>*;
82
83 using PrevSolutionVectorStorage = typename Traits::template Tuple<ConstSubSolutionVectorPtr>;
84
85 using CouplingStencilType = std::vector<std::size_t>;
86
87 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
88
89 using FluidSystem = typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
90
91 using VelocityVector = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
92 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
93
94 struct MomentumCouplingContext
95 {
96 FVElementGeometry<freeFlowMassIndex> fvGeometry;
97 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
98 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
99 std::size_t eIdx;
100 };
101
102 struct MassAndEnergyCouplingContext
103 {
104 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f, const std::size_t i)
105 : fvGeometry(std::move(f))
106 , eIdx(i)
107 {}
108
109 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
110 std::size_t eIdx;
111 };
112
113public:
115 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
116
120 // \{
121
122 //! use as regular coupling manager
123 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
124 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
125 GridVariablesTuple&& gridVariables,
126 const SolutionVector& curSol)
127 {
128 this->momentumCouplingContext_().clear();
129 this->massAndEnergyCouplingContext_().clear();
130
131 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
132 gridVariables_ = gridVariables;
133 this->updateSolution(curSol);
134
135 computeCouplingStencils_();
136 }
137
138 //! use as regular coupling manager in a transient setting
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
147 Dune::Hybrid::forEach(std::make_index_sequence<Traits::numSubDomains>{}, [&](auto i)
148 { std::get<i>(prevSolutions_) = &prevSol[i]; });
149 }
150
151 //! use as binary coupling manager in multi model context
152 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
153 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
154 GridVariablesTuple&& gridVariables,
156 {
157 this->momentumCouplingContext_().clear();
158 this->massAndEnergyCouplingContext_().clear();
159
160 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
161 gridVariables_ = gridVariables;
162 this->attachSolution(curSol);
163
164 computeCouplingStencils_();
165 }
166
167 //! use as binary coupling manager in multi model context and for transient problems
168 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
169 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
170 GridVariablesTuple&& gridVariables,
172 const PrevSolutionVectorStorage& prevSol)
173 {
174 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
175 prevSolutions_ = prevSol;
176 }
177
178 // \}
179
181
199 template<std::size_t j, class LocalAssemblerI>
200 decltype(auto) evalCouplingResidual(Dune::index_constant<freeFlowMomentumIndex> domainI,
201 const LocalAssemblerI& localAssemblerI,
202 const SubControlVolume<freeFlowMomentumIndex>& scvI,
203 Dune::index_constant<j> domainJ,
204 std::size_t dofIdxGlobalJ) const
205 {
206 const auto& problem = localAssemblerI.problem();
207 const auto& element = localAssemblerI.element();
208 const auto& fvGeometry = localAssemblerI.fvGeometry();
209 const auto& curElemVolVars = localAssemblerI.curElemVolVars();
210 const auto& prevElemVolVars = localAssemblerI.prevElemVolVars();
211 typename LocalAssemblerI::ElementResidualVector residual(localAssemblerI.element().subEntities(1));
212 const auto& localResidual = localAssemblerI.localResidual();
213
214 localResidual.evalSource(residual, problem, element, fvGeometry, curElemVolVars, scvI);
215
216 for (const auto& scvf : scvfs(fvGeometry, scvI))
217 localResidual.evalFlux(residual, problem, element, fvGeometry, curElemVolVars, localAssemblerI.elemBcTypes(), localAssemblerI.elemFluxVarsCache(), scvf);
218
219 if (!localAssemblerI.assembler().isStationaryProblem())
220 {
221 assert(isTransient_());
222 localResidual.evalStorage(residual, problem, element, fvGeometry, prevElemVolVars, curElemVolVars, scvI);
223 }
224
225 return residual;
226 }
227
231 // \{
232
235 */
236 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
237 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
238 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
239 {
240 assert(scvf.isFrontal() && !scvf.isLateral() && !scvf.boundary());
241 return this->curSol(freeFlowMassIndex)[fvGeometry.elementIndex()][pressureIdx];
242 }
243
249 */
250 Scalar cellPressure(const Element<freeFlowMassIndex>& element,
251 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
252 {
253 return this->curSol(freeFlowMassIndex)[scvf.insideScvIdx()][pressureIdx];
254 }
255
258 */
259 Scalar density(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 rho = [&](const auto& elemVolVars)
270 {
271 if (scvf.boundary())
272 return 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 // TODO distance weighting
278 return 0.5*(elemVolVars[insideMassScv].density() + elemVolVars[outsideMassScv].density());
279 }
280 };
281
282 return considerPreviousTimeStep ? rho(momentumCouplingContext_()[0].prevElemVolVars)
283 : rho(momentumCouplingContext_()[0].curElemVolVars);
284 }
286 auto insideAndOutsideDensity(const Element<freeFlowMomentumIndex>& element,
287 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
288 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
289 const bool considerPreviousTimeStep = false) const
290 {
291 assert(!(considerPreviousTimeStep && !isTransient_()));
292 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
293 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
294 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
295
296 const auto result = [&](const auto& elemVolVars)
297 {
298 if (scvf.boundary())
299 return std::make_pair(elemVolVars[insideMassScv].density(), elemVolVars[insideMassScv].density());
300 else
301 {
302 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
303 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
304 return std::make_pair(elemVolVars[insideMassScv].density(), elemVolVars[outsideMassScv].density());
305 }
306 };
307
308 return considerPreviousTimeStep ? result(momentumCouplingContext_()[0].prevElemVolVars)
309 : result(momentumCouplingContext_()[0].curElemVolVars);
310 }
311
314 */
315 Scalar density(const Element<freeFlowMomentumIndex>& element,
316 const SubControlVolume<freeFlowMomentumIndex>& scv,
317 const bool considerPreviousTimeStep = false) const
318 {
319 assert(!(considerPreviousTimeStep && !isTransient_()));
320 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
321 const auto& massScv = (*scvs(momentumCouplingContext_()[0].fvGeometry).begin());
322
323 return considerPreviousTimeStep ? momentumCouplingContext_()[0].prevElemVolVars[massScv].density()
324 : momentumCouplingContext_()[0].curElemVolVars[massScv].density();
325 }
326
329 */
330 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
331 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
332 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
333 {
334 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
335
336 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
337 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
338
339 if (scvf.boundary())
340 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
341
342 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
343 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
344
345 const auto mu = [&](const auto& elemVolVars)
346 {
347 // TODO distance weighting
348 return 0.5*(elemVolVars[insideMassScv].viscosity() + elemVolVars[outsideMassScv].viscosity());
349 };
350
351 return mu(momentumCouplingContext_()[0].curElemVolVars);
352 }
353
356 */
357 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
358 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
359 const SubControlVolume<freeFlowMomentumIndex>& scv) const
360 {
361 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
362 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(scv.elementIndex());
363 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
364 }
365
368 */
369 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
370 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
371 {
372 // TODO: rethink this! Maybe we only need scvJ.dofIndex()
373 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, scvf.insideScvIdx()/*eIdx*/);
374
375 // the TPFA scvf index corresponds to the staggered scv index (might need mapping)
376 const auto localMomentumScvIdx = massScvfToMomentumScvIdx_(scvf, massAndEnergyCouplingContext_()[0].fvGeometry);
377 const auto& scvJ = massAndEnergyCouplingContext_()[0].fvGeometry.scv(localMomentumScvIdx);
378
379 // create a unit normal vector oriented in positive coordinate direction
380 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition velocity;
381 velocity[scvJ.dofAxis()] = 1.0;
382
383 // create the actual velocity vector
384 velocity *= this->curSol(freeFlowMomentumIndex)[scvJ.dofIndex()];
385
386 return velocity;
387 }
388
398 template<std::size_t j>
399 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
400 const Element<freeFlowMomentumIndex>& elementI,
401 const SubControlVolume<freeFlowMomentumIndex>& scvI,
402 Dune::index_constant<j> domainJ) const
403 { return emptyStencil_; }
404
418 */
419 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
420 const Element<freeFlowMassIndex>& elementI,
421 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
422 {
423 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
424 return massAndEnergyToMomentumStencils_[eIdx];
425 }
426
435 */
436 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
437 const Element<freeFlowMomentumIndex>& elementI,
438 const SubControlVolume<freeFlowMomentumIndex>& scvI,
439 Dune::index_constant<freeFlowMassIndex> domainJ) const
440 {
441 return momentumToMassAndEnergyStencils_[scvI.index()];
442 }
443
444 // \}
445
449 // \{
450
452 template<std::size_t i, std::size_t j, class LocalAssemblerI>
453 void updateCouplingContext(Dune::index_constant<i> domainI,
454 const LocalAssemblerI& localAssemblerI,
455 Dune::index_constant<j> domainJ,
456 std::size_t dofIdxGlobalJ,
457 const PrimaryVariables<j>& priVarsJ,
458 int pvIdxJ)
459 {
460 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
461
462 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
463 {
464 bindCouplingContext_(domainI, localAssemblerI.element());
465
466 const auto& problem = this->problem(domainJ);
467 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
468 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
469 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
470 const auto scvIdxJ = dofIdxGlobalJ;
471 const auto& scv = fvGeometry.scv(scvIdxJ);
472
473 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
474 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
475 else
476 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
477 }
478 }
479
480 // \}
481
484 */
486 {
487 // use coloring of the fc staggered discretization for both domains
488 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
489 }
490
497 template<std::size_t i, class AssembleElementFunc>
498 void assembleMultithreaded(Dune::index_constant<i> domainI, AssembleElementFunc&& assembleElement) const
499 {
500 if (elementSets_.empty())
501 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
502
503 // make this element loop run in parallel
504 // for this we have to color the elements so that we don't get
505 // race conditions when writing into the global matrix
506 // each color can be assembled using multiple threads
507 const auto& grid = this->problem(freeFlowMomentumIndex).gridGeometry().gridView().grid();
508 for (const auto& elements : elementSets_)
509 {
510 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
511 {
512 const auto element = grid.entity(elements[eIdx]);
513 assembleElement(element);
514 });
515 }
516 }
517
518private:
519 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
520 const Element<freeFlowMomentumIndex>& elementI) const
521 {
522 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
523 bindCouplingContext_(domainI, elementI, eIdx);
524 }
525
526 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
527 const Element<freeFlowMomentumIndex>& elementI,
528 const std::size_t eIdx) const
529 {
530 if (momentumCouplingContext_().empty())
531 {
532 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
533 fvGeometry.bind(elementI);
534
535 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
536 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
537
538 auto prevElemVolVars = isTransient_() ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
539 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
540
541 if (isTransient_())
542 prevElemVolVars.bindElement(elementI, fvGeometry, prevSol_(freeFlowMassIndex));
543
544 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
545 }
546 else if (eIdx != momentumCouplingContext_()[0].eIdx)
547 {
548 momentumCouplingContext_()[0].eIdx = eIdx;
549 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
550 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->curSol(freeFlowMassIndex));
551
552 if (isTransient_())
553 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, prevSol_(freeFlowMassIndex));
554 }
555 }
556
557 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
558 const Element<freeFlowMassIndex>& elementI) const
559 {
560 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
561 bindCouplingContext_(domainI, elementI, eIdx);
562 }
563
564 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
565 const Element<freeFlowMassIndex>& elementI,
566 const std::size_t eIdx) const
567 {
568 if (massAndEnergyCouplingContext_().empty())
569 {
570 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
571 auto fvGeometry = localView(gridGeometry);
572 fvGeometry.bindElement(elementI);
573 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
574 }
575 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
576 {
577 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
578 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
579 }
580 }
581
586 template<std::size_t i>
587 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
588 {
589 if (std::get<i>(gridVariables_))
590 return *std::get<i>(gridVariables_);
591 else
592 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
593 }
594
599 template<std::size_t i>
600 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
601 {
602 if (std::get<i>(gridVariables_))
603 return *std::get<i>(gridVariables_);
604 else
605 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
606 }
607
608
609 void computeCouplingStencils_()
610 {
611 // TODO higher order
612 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
613 auto momentumFvGeometry = localView(momentumGridGeometry);
614 massAndEnergyToMomentumStencils_.clear();
615 massAndEnergyToMomentumStencils_.resize(momentumGridGeometry.gridView().size(0));
616
617 momentumToMassAndEnergyStencils_.clear();
618 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.numScv());
619
620 for (const auto& element : elements(momentumGridGeometry.gridView()))
621 {
622 const auto eIdx = momentumGridGeometry.elementMapper().index(element);
623 momentumFvGeometry.bind(element);
624 for (const auto& scv : scvs(momentumFvGeometry))
625 {
626 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
627 momentumToMassAndEnergyStencils_[scv.index()].push_back(eIdx);
628
629 // extend the stencil for fluids with variable viscosity and density,
630 if constexpr (FluidSystem::isCompressible(0/*phaseIdx*/))
631 // if constexpr (FluidSystem::isCompressible(0/*phaseIdx*/) || !FluidSystem::viscosityIsConstant(0/*phaseIdx*/)) // TODO fix on master
632 {
633 for (const auto& scvf : scvfs(momentumFvGeometry, scv))
634 {
635 if (scvf.isLateral() && !scvf.boundary())
636 {
637 const auto& outsideScv = momentumFvGeometry.scv(scvf.outsideScvIdx());
638 momentumToMassAndEnergyStencils_[scv.index()].push_back(outsideScv.elementIndex());
639 }
640 }
641 }
642 }
643 }
644 }
645
646 std::size_t massScvfToMomentumScvIdx_(const SubControlVolumeFace<freeFlowMassIndex>& massScvf,
647 [[maybe_unused]] const FVElementGeometry<freeFlowMomentumIndex>& momentumFVGeometry) const
648 {
649 if constexpr (ConsistentlyOrientedGrid<typename GridView<freeFlowMomentumIndex>::Grid>{})
650 return massScvf.index();
651 else
652 {
653 static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true);
654 if (!makeConsistentlyOriented)
655 return massScvf.index();
656
657 for (const auto& momentumScv : scvs(momentumFVGeometry))
658 {
659 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition momentumUnitOuterNormal(0.0);
660 momentumUnitOuterNormal[momentumScv.dofAxis()] = momentumScv.directionSign();
661 if (Dune::FloatCmp::eq<typename GridView<freeFlowMomentumIndex>::ctype>(massScvf.unitOuterNormal()*momentumUnitOuterNormal, 1.0))
662 return momentumScv.index();
663 }
664 DUNE_THROW(Dune::InvalidStateException, "No Momentum SCV found");
665 }
666 }
667
668 CouplingStencilType emptyStencil_;
669 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
670 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
671
672 std::vector<MomentumCouplingContext>& momentumCouplingContext_() const
673 { return momentumCouplingContextImpl_; }
674
675 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_() const
676 { return massAndEnergyCouplingContextImpl_; }
677
678 mutable std::vector<MassAndEnergyCouplingContext> massAndEnergyCouplingContextImpl_;
679 mutable std::vector<MomentumCouplingContext> momentumCouplingContextImpl_;
680
682 GridVariablesTuple gridVariables_;
683
684 bool isTransient_() const
685 { return std::get<0>(prevSolutions_) != nullptr; }
686
687 template<std::size_t i>
688 const SubSolutionVector<i>& prevSol_(Dune::index_constant<i>) const
689 { return *std::get<i>(prevSolutions_); }
690
691 PrevSolutionVectorStorage prevSolutions_;
692
693 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
694};
695
696// multi-threading is not supported because we have only one coupling context instance and a mutable cache
697template<class T>
700
701} // end namespace Dumux
702
703#endif
The interface of the coupling manager for multi domain problems.
Definition multidomain/couplingmanager.hh:37
void attachSolution(const SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition multidomain/couplingmanager.hh:310
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition multidomain/couplingmanager.hh:275
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition multidomain/couplingmanager.hh:297
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition multidomain/couplingmanager.hh:326
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:207
CouplingManager()
Default constructor.
Definition multidomain/couplingmanager.hh:70
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition multidomain/couplingmanager.hh:59
The interface of the coupling manager for free flow systems.
Definition couplingmanager_staggered.hh:48
typename ParentType::SolutionVectorStorage SolutionVectorStorage
Definition couplingmanager_staggered.hh:56
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:235
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:329
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:452
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:199
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:497
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:122
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:398
auto insideAndOutsideDensity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Definition couplingmanager_staggered.hh:285
static constexpr auto pressureIdx
Definition couplingmanager_staggered.hh:114
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:258
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:249
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:368
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition couplingmanager_staggered.hh:484
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Type traits.
Helper type to determine whether a grid is guaranteed to be oriented consistently....
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
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:139
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