version 3.10.0
Loading...
Searching...
No Matches
couplingmanager_cvfe.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_CVFE_HH
13#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_CVFE_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#include <dune/geometry/referenceelements.hh>
24
27
31
34
37
38#include "typetraits.hh"
39
40namespace Dumux {
41
47template<class Traits>
49: public CouplingManager<Traits>
50{
51 using ParentType = CouplingManager<Traits>;
52public:
53 static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index();
54 static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index();
55
56 // this can be used if the coupling manager is used inside a meta-coupling manager (e.g. multi-binary)
57 // to manager the solution vector storage outside this class
59private:
60 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
61 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
62 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
63 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
64 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
65 template<std::size_t id> using ElementSeed = typename GridView<id>::Grid::template Codim<0>::EntitySeed;
66 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
67 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
68 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
69 template<std::size_t id> using GridVariables = typename Traits::template SubDomain<id>::GridVariables;
70 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
71 template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache;
72 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
73 template<std::size_t id> using VolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::VolumeVariables>;
74
75 using Scalar = typename Traits::Scalar;
76 using SolutionVector = typename Traits::SolutionVector;
77
78 using CouplingStencilType = std::vector<std::size_t>;
79
80 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
81
82 using FluidSystem = typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
83
84 using VelocityVector = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
85 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
86
87 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
88
89 struct MomentumCouplingContext
90 {
91 FVElementGeometry<freeFlowMassIndex> fvGeometry;
92 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
93 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
94 std::size_t eIdx;
95 };
96
97 struct MassAndEnergyCouplingContext
98 {
99 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f, const std::size_t i)
100 : fvGeometry(std::move(f))
101 , eIdx(i)
102 {}
103
104 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
105 std::size_t eIdx;
106 };
107
108 using MomentumDiscretizationMethod = typename GridGeometry<freeFlowMomentumIndex>::DiscretizationMethod;
109 using MassDiscretizationMethod = typename GridGeometry<freeFlowMassIndex>::DiscretizationMethod;
110
111public:
112
113 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
114
118 // \{
119
121 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
122 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
123 GridVariablesTuple&& gridVariables,
124 const SolutionVector& curSol)
125 {
126 this->momentumCouplingContext_().clear();
127 this->massAndEnergyCouplingContext_().clear();
128
129 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
130 gridVariables_ = gridVariables;
131 this->updateSolution(curSol);
132
133 computeCouplingStencils_();
134 }
135
137 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
138 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
139 GridVariablesTuple&& gridVariables,
140 const SolutionVector& curSol,
141 const SolutionVector& prevSol)
142 {
143 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
144 prevSol_ = &prevSol;
145 isTransient_ = true;
146 }
147
149 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
150 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
151 GridVariablesTuple&& gridVariables,
153 {
154 this->momentumCouplingContext_().clear();
155 this->massAndEnergyCouplingContext_().clear();
156
157 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
158 gridVariables_ = gridVariables;
159 this->attachSolution(curSol);
160
161 computeCouplingStencils_();
162 }
163
164 // \}
165
169 // \{
170
174 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
175 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
176 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
177 const bool considerPreviousTimeStep = false) const
178 {
179 assert(!(considerPreviousTimeStep && !this->isTransient_));
180 const auto& gg = this->problem(freeFlowMassIndex).gridGeometry();
181 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
182 : this->curSol(freeFlowMassIndex);
183 const auto elemSol = elementSolution(element, sol, gg);
184 return evalSolution(element, element.geometry(), gg, elemSol, scvf.ipGlobal())[pressureIdx];
185 }
186
190 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
191 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
192 const SubControlVolume<freeFlowMomentumIndex>& scv,
193 const bool considerPreviousTimeStep = false) const
194 {
195 assert(!(considerPreviousTimeStep && !this->isTransient_));
196 const auto& gg = this->problem(freeFlowMassIndex).gridGeometry();
197 const auto& sol = considerPreviousTimeStep ? (*prevSol_)[freeFlowMassIndex]
198 : this->curSol(freeFlowMassIndex);
199 const auto elemSol = elementSolution(element, sol, gg);
200 return evalSolution(element, element.geometry(), gg, elemSol, scv.dofPosition())[pressureIdx];
201 }
202
206 Scalar density(const Element<freeFlowMomentumIndex>& element,
207 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
208 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
209 const bool considerPreviousTimeStep = false) const
210 {
211 assert(!(considerPreviousTimeStep && !this->isTransient_));
212 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
213
214 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
215 {
216 const auto eIdx = fvGeometry.elementIndex();
217 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
218
219 const auto& volVars = considerPreviousTimeStep ?
220 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
221 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
222
223 return volVars.density();
224 }
225 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
226 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
227 {
228 // TODO: cache the shape values when Box method is used
229 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
230 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
231 std::vector<ShapeValue> shapeValues;
232 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
233 localBasis.evaluateFunction(ipLocal, shapeValues);
234
235 Scalar rho = 0.0;
236 for (const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
237 {
238 const auto& volVars = considerPreviousTimeStep ?
239 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
240 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
241 rho += volVars.density()*shapeValues[scv.indexInElement()][0];
242 }
243
244 return rho;
245 }
246 else
247 DUNE_THROW(Dune::NotImplemented,
248 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
249 );
250 }
251
255 Scalar density(const Element<freeFlowMomentumIndex>& element,
256 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
257 const SubControlVolume<freeFlowMomentumIndex>& scv,
258 const bool considerPreviousTimeStep = false) const
259 {
260 assert(!(considerPreviousTimeStep && !this->isTransient_));
261 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
262
263 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
264 {
265 const auto eIdx = scv.elementIndex();
266 const auto& scvI = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
267
268 const auto& volVars = considerPreviousTimeStep ?
269 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
270 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
271
272 return volVars.density();
273 }
274 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
275 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
276 {
277 // TODO: cache the shape values when Box method is used
278 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
279 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
280 std::vector<ShapeValue> shapeValues;
281 const auto ipLocal = element.geometry().local(scv.dofPosition());
282 localBasis.evaluateFunction(ipLocal, shapeValues);
283
284 Scalar rho = 0.0;
285 for (const auto& scvI : scvs(this->momentumCouplingContext_()[0].fvGeometry))
286 {
287 const auto& volVars = considerPreviousTimeStep ?
288 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
289 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
290 rho += volVars.density()*shapeValues[scvI.indexInElement()][0];
291 }
292 return rho;
293 }
294 else
295 DUNE_THROW(Dune::NotImplemented,
296 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
297 );
298 }
299
303 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
304 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
305 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
306 const bool considerPreviousTimeStep = false) const
307 {
308 assert(!(considerPreviousTimeStep && !this->isTransient_));
309 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
310
311 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
312 {
313 const auto eIdx = fvGeometry.elementIndex();
314 const auto& scv = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
315 const auto& volVars = considerPreviousTimeStep ?
316 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
317 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
318 return volVars.viscosity();
319 }
320 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
321 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
322 {
323 // TODO: cache the shape values when Box method is used
324 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
325 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
326 std::vector<ShapeValue> shapeValues;
327 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
328 localBasis.evaluateFunction(ipLocal, shapeValues);
329
330 Scalar mu = 0.0;
331 for (const auto& scv : scvs(this->momentumCouplingContext_()[0].fvGeometry))
332 {
333 const auto& volVars = considerPreviousTimeStep ?
334 this->momentumCouplingContext_()[0].prevElemVolVars[scv]
335 : this->momentumCouplingContext_()[0].curElemVolVars[scv];
336 mu += volVars.viscosity()*shapeValues[scv.indexInElement()][0];
337 }
338
339 return mu;
340 }
341 else
342 DUNE_THROW(Dune::NotImplemented,
343 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
344 );
345 }
346
350 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
351 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
352 const SubControlVolume<freeFlowMomentumIndex>& scv,
353 const bool considerPreviousTimeStep = false) const
354 {
355 assert(!(considerPreviousTimeStep && !this->isTransient_));
356 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
357
358 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
359 {
360 const auto eIdx = fvGeometry.elementIndex();
361 const auto& scvI = this->momentumCouplingContext_()[0].fvGeometry.scv(eIdx);
362 const auto& volVars = considerPreviousTimeStep ?
363 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
364 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
365 return volVars.viscosity();
366 }
367 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
368 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
369 {
370 // TODO: cache the shape values when Box method is used
371 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
372 const auto& localBasis = this->momentumCouplingContext_()[0].fvGeometry.feLocalBasis();
373 std::vector<ShapeValue> shapeValues;
374 const auto ipLocal = element.geometry().local(scv.dofPosition());
375 localBasis.evaluateFunction(ipLocal, shapeValues);
376
377 Scalar mu = 0.0;
378 for (const auto& scvI : scvs(this->momentumCouplingContext_()[0].fvGeometry))
379 {
380 const auto& volVars = considerPreviousTimeStep ?
381 this->momentumCouplingContext_()[0].prevElemVolVars[scvI]
382 : this->momentumCouplingContext_()[0].curElemVolVars[scvI];
383 mu += volVars.viscosity()*shapeValues[scvI.indexInElement()][0];
384 }
385
386 return mu;
387 }
388 else
389 DUNE_THROW(Dune::NotImplemented,
390 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
391 );
392 }
393
397 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
398 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
399 {
400 // TODO: optimize this function for tpfa where the scvf ip coincides with the dof location
401
402 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(element);
403 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, eIdx);
404
405 const auto& fvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
406 const auto& localBasis = fvGeometry.feLocalBasis();
407
408 std::vector<ShapeValue> shapeValues;
409 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
410 localBasis.evaluateFunction(ipLocal, shapeValues);
411
412 // interpolate velocity at scvf
413 VelocityVector velocity(0.0);
414 for (const auto& scv : scvs(fvGeometry))
415 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(freeFlowMomentumIndex)[scv.dofIndex()]);
416
417 return velocity;
418 }
419
423 VelocityVector elementVelocity(const FVElementGeometry<freeFlowMassIndex>& fvGeometry) const
424 {
425 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), fvGeometry.element());
426
427 const auto& momentumFvGeometry = this->massAndEnergyCouplingContext_()[0].fvGeometry;
428 const auto& localBasis = momentumFvGeometry.feLocalBasis();
429
430 // interpolate velocity at scvf
431 VelocityVector velocity(0.0);
432 std::vector<ShapeValue> shapeValues;
433 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
434
435 for (const auto& scv : scvs(momentumFvGeometry))
436 velocity.axpy(shapeValues[scv.localDofIndex()][0], this->curSol(freeFlowMomentumIndex)[scv.dofIndex()]);
437
438 return velocity;
439 }
440
445 template<std::size_t j>
446 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
447 const Element<freeFlowMomentumIndex>& elementI,
448 const SubControlVolume<freeFlowMomentumIndex>& scvI,
449 Dune::index_constant<j> domainJ) const
450 { return emptyStencil_; }
451
466 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
467 const Element<freeFlowMassIndex>& elementI,
468 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
469 {
470 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
471 return massAndEnergyToMomentumStencils_[eIdx];
472 }
473
482 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
483 const Element<freeFlowMomentumIndex>& elementI,
484 Dune::index_constant<freeFlowMassIndex> domainJ) const
485 {
486 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
487 return momentumToMassAndEnergyStencils_[eIdx];
488 }
489
490 // \}
491
495 // \{
496
516 template<std::size_t i, std::size_t j, class LocalAssemblerI>
517 void updateCouplingContext(Dune::index_constant<i> domainI,
518 const LocalAssemblerI& localAssemblerI,
519 Dune::index_constant<j> domainJ,
520 std::size_t dofIdxGlobalJ,
521 const PrimaryVariables<j>& priVarsJ,
522 int pvIdxJ)
523 {
524 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
525
526 if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::cctpfa)
527 {
528 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
529 {
530 bindCouplingContext_(domainI, localAssemblerI.element());
531
532 const auto& problem = this->problem(domainJ);
533 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
534 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
535 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
536 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
537
538 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
539 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
540 else
541 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
542 }
543 }
544 else if constexpr (MassDiscretizationMethod{} == DiscretizationMethods::box
545 || MassDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
546 {
547 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
548 {
549 bindCouplingContext_(domainI, localAssemblerI.element());
550
551 const auto& problem = this->problem(domainJ);
552 const auto& deflectedElement = problem.gridGeometry().element(this->momentumCouplingContext_()[0].eIdx);
553 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
554 const auto& fvGeometry = this->momentumCouplingContext_()[0].fvGeometry;
555
556 for (const auto& scv : scvs(fvGeometry))
557 {
558 if(scv.dofIndex() == dofIdxGlobalJ)
559 {
560 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
561 this->gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
562 else
563 this->momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
564 }
565 }
566 }
567 }
568 else
569 DUNE_THROW(Dune::NotImplemented,
570 "Context update for discretization scheme " << MassDiscretizationMethod{}
571 );
572 }
573
574 // \}
575
580 {
581 if constexpr (MomentumDiscretizationMethod{} == DiscretizationMethods::fcdiamond)
582 {
583 // use coloring of the mass discretization for both domains
584 // the diamond coloring is a subset (minimum amount of colors) of cctpfa/box coloring
585 elementSets_ = computeColoring(this->problem(freeFlowMassIndex).gridGeometry()).sets;
586 }
587 else
588 {
589 // use coloring of the momentum discretization for both domains
590 elementSets_ = computeColoring(this->problem(freeFlowMomentumIndex).gridGeometry()).sets;
591 }
592 }
593
600 template<std::size_t i, class AssembleElementFunc>
601 void assembleMultithreaded(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
602 {
603 if (elementSets_.empty())
604 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
605
606 // make this element loop run in parallel
607 // for this we have to color the elements so that we don't get
608 // race conditions when writing into the global matrix
609 // each color can be assembled using multiple threads
610 const auto& grid = this->problem(freeFlowMassIndex).gridGeometry().gridView().grid();
611 for (const auto& elements : elementSets_)
612 {
613 Dumux::parallelFor(elements.size(), [&](const std::size_t eIdx)
614 {
615 const auto element = grid.entity(elements[eIdx]);
616 assembleElement(element);
617 });
618 }
619 }
620
621private:
622 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
623 const Element<freeFlowMomentumIndex>& elementI) const
624 {
625 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
626 if (momentumCouplingContext_().empty())
627 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI));
628 else
629 bindCouplingContext_(domainI, elementI, momentumCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
630 }
631
632 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
633 const Element<freeFlowMomentumIndex>& elementI,
634 const std::size_t eIdx) const
635 {
636 if (momentumCouplingContext_().empty())
637 {
638 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
639 fvGeometry.bind(elementI);
640
641 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
642 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
643
644 auto prevElemVolVars = isTransient_ ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
645 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
646
647 if (isTransient_)
648 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[freeFlowMassIndex]);
649
650 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
651 }
652 else if (eIdx != momentumCouplingContext_()[0].eIdx)
653 {
654 momentumCouplingContext_()[0].eIdx = eIdx;
655 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
656 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->curSol(freeFlowMassIndex));
657
658 if (isTransient_)
659 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[freeFlowMassIndex]);
660 }
661 }
662
663 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
664 const Element<freeFlowMassIndex>& elementI) const
665 {
666 // The call to this->problem() is expensive because of std::weak_ptr (see base class). Here we try to avoid it if possible.
667 if (massAndEnergyCouplingContext_().empty())
668 bindCouplingContext_(domainI, elementI, this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI));
669 else
670 bindCouplingContext_(domainI, elementI, massAndEnergyCouplingContext_()[0].fvGeometry.gridGeometry().elementMapper().index(elementI));
671 }
672
673 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
674 const Element<freeFlowMassIndex>& elementI,
675 const std::size_t eIdx) const
676 {
677 if (massAndEnergyCouplingContext_().empty())
678 {
679 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
680 auto fvGeometry = localView(gridGeometry);
681 fvGeometry.bindElement(elementI);
682 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
683 }
684 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
685 {
686 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
687 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
688 }
689 }
690
695 template<std::size_t i>
696 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
697 {
698 if (std::get<i>(gridVariables_))
699 return *std::get<i>(gridVariables_);
700 else
701 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
702 }
703
708 template<std::size_t i>
709 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
710 {
711 if (std::get<i>(gridVariables_))
712 return *std::get<i>(gridVariables_);
713 else
714 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
715 }
716
717
718 void computeCouplingStencils_()
719 {
720 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
721 const auto& massGridGeometry = this->problem(freeFlowMassIndex).gridGeometry();
722 auto momentumFvGeometry = localView(momentumGridGeometry);
723 auto massFvGeometry = localView(massGridGeometry);
724
725 massAndEnergyToMomentumStencils_.clear();
726 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
727
728 momentumToMassAndEnergyStencils_.clear();
729 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
730
731 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
732
733 for (const auto& element : elements(momentumGridGeometry.gridView()))
734 {
735 momentumFvGeometry.bindElement(element);
736 massFvGeometry.bindElement(element);
737 const auto eIdx = momentumFvGeometry.elementIndex();
738
739 for (const auto& scv : scvs(momentumFvGeometry))
740 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
741
742 for (const auto& scv : scvs(massFvGeometry))
743 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
744 }
745 }
746
747 CouplingStencilType emptyStencil_;
748 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
749 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
750
751 std::vector<MomentumCouplingContext>& momentumCouplingContext_() const
752 { return momentumCouplingContextImpl_; }
753
754 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_() const
755 { return massAndEnergyCouplingContextImpl_; }
756
757 mutable std::vector<MassAndEnergyCouplingContext> massAndEnergyCouplingContextImpl_;
758 mutable std::vector<MomentumCouplingContext> momentumCouplingContextImpl_;
759
761 GridVariablesTuple gridVariables_;
762
763 const SolutionVector* prevSol_;
764 bool isTransient_;
765
766 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
767};
768
769namespace Detail {
770
771// declaration (specialize for different discretization types)
772template<class Traits, class DiscretizationMethod = typename Detail::MomentumDiscretizationMethod<Traits>::type>
774
775// multi-threading is not supported because we have only one coupling context instance and a mutable cache
776template<class Traits, class D>
778{ using type = std::false_type; };
779
780} // end namespace Detail
781
783template<class T>
787
788} // end namespace Dumux
789
790#endif
The interface of the coupling manager for free flow systems.
Definition couplingmanager_cvfe.hh:50
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition couplingmanager_cvfe.hh:397
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol)
Methods to be accessed by main.
Definition couplingmanager_cvfe.hh:121
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const typename ParentType::SolutionVectorStorage &curSol)
use as binary coupling manager in multi model context
Definition couplingmanager_cvfe.hh:149
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the effective viscosity at a given sub control volume.
Definition couplingmanager_cvfe.hh:350
static constexpr auto freeFlowMomentumIndex
Definition couplingmanager_cvfe.hh:53
static constexpr auto pressureIdx
Definition couplingmanager_cvfe.hh:113
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_cvfe.hh:466
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, Dune::index_constant< freeFlowMassIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition couplingmanager_cvfe.hh:482
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_cvfe.hh:206
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the effective viscosity at a given sub control volume face.
Definition couplingmanager_cvfe.hh:303
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition couplingmanager_cvfe.hh:423
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition couplingmanager_cvfe.hh:255
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the pressure at a given sub control volume.
Definition couplingmanager_cvfe.hh:190
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition couplingmanager_cvfe.hh:601
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the pressure at a given sub control volume face.
Definition couplingmanager_cvfe.hh:174
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_cvfe.hh:137
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition couplingmanager_cvfe.hh:579
static constexpr auto freeFlowMassIndex
Definition couplingmanager_cvfe.hh:54
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_cvfe.hh:446
typename ParentType::SolutionVectorStorage SolutionVectorStorage
Definition couplingmanager_cvfe.hh:58
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
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
PrimaryVariables evalSolution(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition evalsolution.hh:152
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 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_cvfe.hh:517
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.
Some useful type traits.
A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
Distance implementation details.
Definition cvfelocalresidual.hh:25
Definition method.hh:20
constexpr FCDiamond fcdiamond
Definition method.hh:152
constexpr CCTpfa cctpfa
Definition method.hh:145
constexpr Box box
Definition method.hh:147
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