Loading [MathJax]/extensions/tex2jax.js
3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
couplingmanagerbase.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGERBASE_HH
27#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGERBASE_HH
28
29#include <iostream>
30#include <fstream>
31#include <string>
32#include <utility>
33#include <unordered_map>
34
35#include <dune/common/timer.hh>
36#include <dune/geometry/quadraturerules.hh>
37
38#include <dumux/io/format.hh>
51
52namespace Dumux {
53
55template<class MDTraits>
57{
58private:
59 template<std::size_t i> using SubDomainTypeTag = typename MDTraits::template SubDomain<i>::TypeTag;
60 template<std::size_t i> using GridGeometry = GetPropType<SubDomainTypeTag<i>, Properties::GridGeometry>;
61 template<std::size_t i> using NumEqVector = Dumux::NumEqVector<GetPropType<SubDomainTypeTag<i>, Properties::PrimaryVariables>>;
62public:
64 template<std::size_t i>
66
68 template<std::size_t i>
70
73};
74
80template<class MDTraits, class Implementation, class PSTraits = DefaultPointSourceTraits<MDTraits>>
82: public CouplingManager<MDTraits>
83{
85 using Scalar = typename MDTraits::Scalar;
86 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
87 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
88 using SolutionVector = typename MDTraits::SolutionVector;
89 using PointSourceData = typename PSTraits::PointSourceData;
90
91 // the sub domain type tags
92 template<std::size_t id> using PointSource = typename PSTraits::template PointSource<id>;
93 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
94 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
95 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
96 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
97 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
98 template<std::size_t id> using ElementMapper = typename GridGeometry<id>::ElementMapper;
99 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
100 template<std::size_t id> using ElementSeed = typename GridView<id>::Grid::template Codim<0>::EntitySeed;
101 template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex;
102 template<std::size_t id> using CouplingStencil = std::vector<GridIndex<id>>;
103
104 static constexpr int bulkDim = GridView<bulkIdx>::dimension;
105 static constexpr int lowDimDim = GridView<lowDimIdx>::dimension;
106 static constexpr int dimWorld = GridView<bulkIdx>::dimensionworld;
107
108 template<std::size_t id>
109 static constexpr bool isBox()
110 { return GridGeometry<id>::discMethod == DiscretizationMethods::box; }
111
112 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
113 using GlueType = MultiDomainGlue<GridView<bulkIdx>, GridView<lowDimIdx>, ElementMapper<bulkIdx>, ElementMapper<lowDimIdx>>;
114
115public:
117 using MultiDomainTraits = MDTraits;
119 using PointSourceTraits = PSTraits;
121 template<std::size_t id> using CouplingStencils = std::unordered_map<GridIndex<id>, CouplingStencil<id>>;
122
126 void updateAfterGridAdaption(std::shared_ptr<const GridGeometry<bulkIdx>> bulkGridGeometry,
127 std::shared_ptr<const GridGeometry<lowDimIdx>> lowDimGridGeometry)
128 {
129 glue_ = std::make_shared<GlueType>();
130 }
131
135 EmbeddedCouplingManagerBase(std::shared_ptr<const GridGeometry<bulkIdx>> bulkGridGeometry,
136 std::shared_ptr<const GridGeometry<lowDimIdx>> lowDimGridGeometry)
137 {
138 updateAfterGridAdaption(bulkGridGeometry, lowDimGridGeometry);
139 }
140
144 // \{
145
146 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
147 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
148 const SolutionVector& curSol)
149 {
150 this->updateSolution(curSol);
151 this->setSubProblems(std::make_tuple(bulkProblem, lowDimProblem));
152
153 integrationOrder_ = getParam<int>("MixedDimension.IntegrationOrder", 1);
154 asImp_().computePointSourceData(integrationOrder_);
155 }
156
157 // \}
158
162 // \{
163
178 template<std::size_t i, std::size_t j>
179 const CouplingStencil<j>& couplingStencil(Dune::index_constant<i> domainI,
180 const Element<i>& element,
181 Dune::index_constant<j> domainJ) const
182 {
183 static_assert(i != j, "A domain cannot be coupled to itself!");
184
185 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
186 if (couplingStencils(domainI).count(eIdx))
187 return couplingStencils(domainI).at(eIdx);
188 else
189 return emptyStencil(domainI);
190 }
191
204 template<std::size_t i, std::size_t j, class LocalAssemblerI>
205 decltype(auto) evalCouplingResidual(Dune::index_constant<i> domainI,
206 const LocalAssemblerI& localAssemblerI,
207 Dune::index_constant<j> domainJ,
208 std::size_t dofIdxGlobalJ)
209 {
210 static_assert(i != j, "A domain cannot be coupled to itself!");
211
212 typename LocalAssemblerI::LocalResidual::ElementResidualVector residual;
213
214 const auto& element = localAssemblerI.element();
215 const auto& fvGeometry = localAssemblerI.fvGeometry();
216 const auto& curElemVolVars = localAssemblerI.curElemVolVars();
217
218 residual.resize(fvGeometry.numScv());
219 for (const auto& scv : scvs(fvGeometry))
220 {
221 auto couplingSource = this->problem(domainI).scvPointSources(element, fvGeometry, curElemVolVars, scv);
222 couplingSource += this->problem(domainI).source(element, fvGeometry, curElemVolVars, scv);
223 couplingSource *= -GridGeometry<i>::Extrusion::volume(fvGeometry, scv)*curElemVolVars[scv].extrusionFactor();
224 residual[scv.indexInElement()] = couplingSource;
225 }
226 return residual;
227 }
228
229 // \}
230
231 /* \brief Compute integration point point sources and associated data
232 *
233 * This method uses grid glue to intersect the given grids. Over each intersection
234 * we later need to integrate a source term. This method places point sources
235 * at each quadrature point and provides the point source with the necessary
236 * information to compute integrals (quadrature weight and integration element)
237 * \param order The order of the quadrature rule for integration of sources over an intersection
238 * \param verbose If the point source computation is verbose
239 */
240 void computePointSourceData(std::size_t order = 1, bool verbose = false)
241 {
242 // initialize the maps
243 // do some logging and profiling
244 Dune::Timer watch;
245 std::cout << "Initializing the point sources..." << std::endl;
246
247 // clear all internal members like pointsource vectors and stencils
248 // initializes the point source id counter
249 clear();
250
251 // precompute the vertex indices for efficiency for the box method
252 this->precomputeVertexIndices(bulkIdx);
253 this->precomputeVertexIndices(lowDimIdx);
254
255 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
256 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
257
258 // intersect the bounding box trees
259 glueGrids();
260
261 pointSourceData_.reserve(this->glue().size());
262 averageDistanceToBulkCell_.reserve(this->glue().size());
263 for (const auto& is : intersections(this->glue()))
264 {
265 // all inside elements are identical...
266 const auto& inside = is.targetEntity(0);
267 // get the intersection geometry for integrating over it
268 const auto intersectionGeometry = is.geometry();
269
270 // get the Gaussian quadrature rule for the local intersection
271 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order);
272 const std::size_t lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
273
274 // iterate over all quadrature points
275 for (auto&& qp : quad)
276 {
277 // compute the coupling stencils
278 for (std::size_t outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
279 {
280 const auto& outside = is.domainEntity(outsideIdx);
281 const std::size_t bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
282
283 // each quadrature point will be a point source for the sub problem
284 const auto globalPos = intersectionGeometry.global(qp.position());
285 const auto id = idCounter_++;
286 const auto qpweight = qp.weight();
287 const auto ie = intersectionGeometry.integrationElement(qp.position());
288 pointSources(bulkIdx).emplace_back(globalPos, id, qpweight, ie, bulkElementIdx);
289 pointSources(bulkIdx).back().setEmbeddings(is.numDomainNeighbors());
290 pointSources(lowDimIdx).emplace_back(globalPos, id, qpweight, ie, lowDimElementIdx);
291 pointSources(lowDimIdx).back().setEmbeddings(is.numDomainNeighbors());
292
293 // pre compute additional data used for the evaluation of
294 // the actual solution dependent source term
295 PointSourceData psData;
296
297 if constexpr (isBox<lowDimIdx>())
298 {
299 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
300 const auto lowDimGeometry = this->problem(lowDimIdx).gridGeometry().element(lowDimElementIdx).geometry();
301 ShapeValues shapeValues;
302 this->getShapeValues(lowDimIdx, this->problem(lowDimIdx).gridGeometry(), lowDimGeometry, globalPos, shapeValues);
303 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
304 }
305 else
306 {
307 psData.addLowDimInterpolation(lowDimElementIdx);
308 }
309
310 // add data needed to compute integral over the circle
311 if constexpr (isBox<bulkIdx>())
312 {
313 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
314 const auto bulkGeometry = this->problem(bulkIdx).gridGeometry().element(bulkElementIdx).geometry();
315 ShapeValues shapeValues;
316 this->getShapeValues(bulkIdx, this->problem(bulkIdx).gridGeometry(), bulkGeometry, globalPos, shapeValues);
317 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
318 }
319 else
320 {
321 psData.addBulkInterpolation(bulkElementIdx);
322 }
323
324 // publish point source data in the global vector
325 this->pointSourceData().emplace_back(std::move(psData));
326
327 // compute average distance to bulk cell
328 averageDistanceToBulkCell_.push_back(averageDistancePointGeometry(globalPos, outside.geometry()));
329
330 // export the lowdim coupling stencil
331 // we insert all vertices / elements and make it unique later
332 if (isBox<bulkIdx>())
333 {
334 const auto& vertices = this->vertexIndices(bulkIdx, bulkElementIdx);
335 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
336 vertices.begin(), vertices.end());
337 }
338 else
339 {
340 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
341 }
342
343 // export the bulk coupling stencil
344 // we insert all vertices / elements and make it unique later
345 if (isBox<lowDimIdx>())
346 {
347 const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx);
348 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
349 vertices.begin(), vertices.end());
350
351 }
352 else
353 {
354 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
355 }
356 }
357 }
358 }
359
360 // make stencils unique
361 using namespace Dune::Hybrid;
362 forEach(integralRange(std::integral_constant<std::size_t, MDTraits::numSubDomains>{}), [&](const auto domainIdx)
363 {
364 for (auto&& stencil : this->couplingStencils(domainIdx))
365 {
366 std::sort(stencil.second.begin(), stencil.second.end());
367 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
368 }
369 });
370
371 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
372 }
373
377 // \{
378
380 const PointSourceData& pointSourceData(std::size_t id) const
381 { return pointSourceData_[id]; }
382
384 template<std::size_t id>
385 const GridView<id>& gridView(Dune::index_constant<id> domainIdx) const
386 { return this->problem(domainIdx).gridGeometry().gridView(); }
387
389 PrimaryVariables<bulkIdx> bulkPriVars(std::size_t id) const
390 { return pointSourceData_[id].interpolateBulk(this->curSol(bulkIdx)); }
391
393 PrimaryVariables<lowDimIdx> lowDimPriVars(std::size_t id) const
394 { return pointSourceData_[id].interpolateLowDim(this->curSol(lowDimIdx)); }
395
397 Scalar averageDistance(std::size_t id) const
398 { return averageDistanceToBulkCell_[id]; }
399
401 const std::vector<PointSource<bulkIdx>>& bulkPointSources() const
402 { return std::get<bulkIdx>(pointSources_); }
403
405 const std::vector<PointSource<lowDimIdx>>& lowDimPointSources() const
406 { return std::get<lowDimIdx>(pointSources_); }
407
409 template<std::size_t i>
410 const std::vector<PointSource<i>>& pointSources(Dune::index_constant<i> dom) const
411 { return std::get<i>(pointSources_); }
412
414 template<std::size_t i>
415 const CouplingStencils<i>& couplingStencils(Dune::index_constant<i> dom) const
416 { return std::get<i>(couplingStencils_); }
417
419 const std::vector<PointSourceData>& pointSourceData() const
420 { return pointSourceData_; }
421
423 template<std::size_t i>
424 const CouplingStencil<i>& emptyStencil(Dune::index_constant<i> dom) const
425 { return std::get<i>(emptyStencil_); }
426
432 template<std::size_t i>
433 auto& curSol(Dune::index_constant<i> domainIdx)
434 { return ParentType::curSol(domainIdx); }
435
441 template<std::size_t i>
442 const auto& curSol(Dune::index_constant<i> domainIdx) const
443 { return ParentType::curSol(domainIdx); }
444
449 {
450 Dune::Timer timer;
451
452 // start with elements connected within the subdomains
453 const auto& bulkGG = asImp_().problem(bulkIdx).gridGeometry();
454 const auto& lowDimGG = asImp_().problem(lowDimIdx).gridGeometry();
455
456 auto connectedElementsBulk = Detail::computeConnectedElements(bulkGG);
457 auto connectedElementsLowDim = Detail::computeConnectedElements(lowDimGG);
458
459 // to the elements from the own domain, we have to add the extended source stencil (non-local couplings)
460 // we consider this only in the bulk domain
461 for (const auto& element : elements(bulkGG.gridView()))
462 {
463 const auto eIdx = bulkGG.elementMapper().index(element);
464 const auto& extendedSourceStencil = asImp_().extendedSourceStencil(eIdx);
465 for (auto dofIdx : extendedSourceStencil)
466 {
467 const auto& elems = connectedElementsBulk[dofIdx];
468 if (std::find(elems.begin(), elems.end(), eIdx) == elems.end())
469 connectedElementsBulk[dofIdx].push_back(eIdx);
470 }
471 }
472
473 // then we also need a similar container that tells us
474 // which elements _from the other domain_ modify our dofs.
475 // That, we can obtain via the coupling stencils
476 std::array<std::vector<std::vector<std::size_t>>, 2> connectedElementsCoupling;
477
478 using namespace Dune::Hybrid;
479 forEach(integralRange(Dune::index_constant<MDTraits::numSubDomains>{}), [&](const auto i)
480 {
481 connectedElementsCoupling[i()].resize(asImp_().problem(i).gridGeometry().numDofs());
482 forEach(integralRange(Dune::index_constant<MDTraits::numSubDomains>{}), [&](const auto j)
483 {
484 if constexpr (i != j)
485 {
486 for (const auto& element : elements(asImp_().problem(j).gridGeometry().gridView()))
487 {
488 const auto eIdx = asImp_().problem(j).gridGeometry().elementMapper().index(element);
489 for (auto dofIdx : asImp_().couplingStencil(j, element, i))
490 {
491 const auto& elems = connectedElementsCoupling[i()][dofIdx];
492 if (std::find(elems.begin(), elems.end(), eIdx) == elems.end())
493 connectedElementsCoupling[i()][dofIdx].push_back(eIdx);
494 }
495 }
496 }
497 });
498 });
499
500 forEach(integralRange(Dune::index_constant<MDTraits::numSubDomains>{}), [&](const auto i)
501 {
502 const auto& gg = asImp_().problem(i).gridGeometry();
503
504 std::vector<int> colors(gg.gridView().size(0), -1);
505
506 // pre-reserve some memory for helper arrays to avoid reallocation
507 std::vector<int> neighborColors; neighborColors.reserve(200);
508 std::vector<bool> colorUsed; colorUsed.reserve(200);
509
510 auto& elementSets = std::get<i>(elementSets_);
511 const auto& connectedElements = std::get<i>(std::tie(connectedElementsBulk, connectedElementsLowDim));
512
513 for (const auto& element : elements(gg.gridView()))
514 {
515 // compute neighbor colors based on discretization-dependent stencil
516 neighborColors.clear();
517 Detail::addNeighborColors(gg, element, colors, connectedElements, neighborColors);
518
519 if constexpr (i == bulkIdx)
520 {
521 // make sure we also add all elements in the extended source stencil
522 // Detail::addNeighborColors only adds direct neighbors
523 // add everyone that also modifies the extended source stencil dofs
524 const auto eIdx = bulkGG.elementMapper().index(element);
525 const auto& extendedSourceStencil = asImp_().extendedSourceStencil(eIdx);
526 for (auto dofIdx : extendedSourceStencil)
527 for (const auto nIdx : connectedElementsBulk[dofIdx])
528 neighborColors.push_back(colors[nIdx]);
529 }
530
531 // add neighbor colors based on coupling stencil
532 forEach(integralRange(Dune::index_constant<MDTraits::numSubDomains>{}), [&](const auto j)
533 {
534 if constexpr (i != j)
535 {
536 // add all elements that manipulate the same coupled dofs
537 for (auto dofIdx : asImp_().couplingStencil(i, element, j))
538 for (auto eIdx : connectedElementsCoupling[j][dofIdx])
539 neighborColors.push_back(colors[eIdx]);
540 }
541 });
542
543 // find smallest color (positive integer) not in neighborColors
544 const auto color = Detail::smallestAvailableColor(neighborColors, colorUsed);
545
546 // assign color to element
547 colors[gg.elementMapper().index(element)] = color;
548
549 // add element to the set of elements with the same color
550 if (color < elementSets.size())
551 elementSets[color].push_back(element.seed());
552 else
553 elementSets.push_back(std::vector<ElementSeed<i>>{ element.seed() });
554 }
555 });
556
557 std::cout << Fmt::format("Colored in {} seconds:\n", timer.elapsed());
558 forEach(integralRange(Dune::index_constant<MDTraits::numSubDomains>{}), [&](const auto i)
559 {
560 std::cout << Fmt::format(
561 "-- {} elements in subdomain {} with {} colors\n",
562 asImp_().problem(i).gridGeometry().gridView().size(0),
563 i(), std::get<i>(elementSets_).size()
564 );
565 });
566 }
567
574 template<std::size_t i, class AssembleElementFunc>
575 void assembleMultithreaded(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
576 {
577 if (std::get<i>(elementSets_).empty())
578 DUNE_THROW(Dune::InvalidStateException, "Call computeColorsForAssembly before assembling in parallel!");
579
580 // make this element loop run in parallel
581 // for this we have to color the elements so that we don't get
582 // race conditions when writing into the global matrix or modifying grid variable caches
583 // each color can be assembled using multiple threads
584 const auto& grid = this->problem(domainId).gridGeometry().gridView().grid();
585 for (const auto& elements : std::get<i>(elementSets_))
586 {
587 Dumux::parallelFor(elements.size(), [&](const std::size_t n)
588 {
589 const auto element = grid.entity(elements[n]);
590 assembleElement(element);
591 });
592 }
593 }
594
600 const CouplingStencil<bulkIdx>& extendedSourceStencil(std::size_t eIdx) const
601 { return asImp_().emptyStencil(bulkIdx); }
602
603protected:
604 using ParentType::curSol;
605
607 template<std::size_t id>
608 void precomputeVertexIndices(Dune::index_constant<id> domainIdx)
609 {
610 // fill helper structure for box discretization
611 if constexpr (isBox<domainIdx>())
612 {
613 this->vertexIndices(domainIdx).resize(gridView(domainIdx).size(0));
614 for (const auto& element : elements(gridView(domainIdx)))
615 {
616 constexpr int dim = GridView<domainIdx>::dimension;
617 const auto eIdx = this->problem(domainIdx).gridGeometry().elementMapper().index(element);
618 this->vertexIndices(domainIdx, eIdx).resize(element.subEntities(dim));
619 for (int i = 0; i < element.subEntities(dim); ++i)
620 this->vertexIndices(domainIdx, eIdx)[i] = this->problem(domainIdx).gridGeometry().vertexMapper().subIndex(element, i, dim);
621 }
622 }
623 }
624
626 template<std::size_t i, class FVGG, class Geometry, class ShapeValues>
627 void getShapeValues(Dune::index_constant<i> domainI, const FVGG& gridGeometry, const Geometry& geo, const GlobalPosition& globalPos, ShapeValues& shapeValues)
628 {
629 if constexpr (FVGG::discMethod == DiscretizationMethods::box)
630 {
631 const auto ipLocal = geo.local(globalPos);
632 const auto& localBasis = this->problem(domainI).gridGeometry().feCache().get(geo.type()).localBasis();
633 localBasis.evaluateFunction(ipLocal, shapeValues);
634 }
635 else
636 DUNE_THROW(Dune::InvalidStateException, "Shape values requested for other discretization than box!");
637 }
638
640 void clear()
641 {
642 pointSources(bulkIdx).clear();
643 pointSources(lowDimIdx).clear();
644 couplingStencils(bulkIdx).clear();
645 couplingStencils(lowDimIdx).clear();
646 vertexIndices(bulkIdx).clear();
647 vertexIndices(lowDimIdx).clear();
648 pointSourceData_.clear();
649 averageDistanceToBulkCell_.clear();
650
651 idCounter_ = 0;
652 }
653
656 {
657 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
658 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
659
660 // intersect the bounding box trees
661 glue_->build(bulkGridGeometry.boundingBoxTree(), lowDimGridGeometry.boundingBoxTree());
662 }
663
665 std::vector<PointSourceData>& pointSourceData()
666 { return pointSourceData_; }
667
669 std::vector<Scalar>& averageDistanceToBulkCell()
670 { return averageDistanceToBulkCell_; }
671
673 template<std::size_t i>
674 std::vector<PointSource<i>>& pointSources(Dune::index_constant<i> dom)
675 { return std::get<i>(pointSources_); }
676
678 template<std::size_t i>
679 CouplingStencils<i>& couplingStencils(Dune::index_constant<i> dom)
680 { return std::get<i>(couplingStencils_); }
681
683 template<std::size_t i>
684 std::vector<GridIndex<i>>& vertexIndices(Dune::index_constant<i> dom, GridIndex<i> eIdx)
685 { return std::get<i>(vertexIndices_)[eIdx]; }
686
688 template<std::size_t i>
689 std::vector<std::vector<GridIndex<i>>>& vertexIndices(Dune::index_constant<i> dom)
690 { return std::get<i>(vertexIndices_); }
691
692 const GlueType& glue() const
693 { return *glue_; }
694
696 Implementation &asImp_()
697 { return *static_cast<Implementation *>(this); }
698
700 const Implementation &asImp_() const
701 { return *static_cast<const Implementation *>(this); }
702
704 std::size_t idCounter_ = 0;
705
706private:
707
709 std::tuple<std::vector<PointSource<bulkIdx>>, std::vector<PointSource<lowDimIdx>>> pointSources_;
710 std::vector<PointSourceData> pointSourceData_;
711 std::vector<Scalar> averageDistanceToBulkCell_;
712
714 std::tuple<std::vector<std::vector<GridIndex<bulkIdx>>>,
715 std::vector<std::vector<GridIndex<lowDimIdx>>>> vertexIndices_;
716 std::tuple<CouplingStencils<bulkIdx>, CouplingStencils<lowDimIdx>> couplingStencils_;
717 std::tuple<CouplingStencil<bulkIdx>, CouplingStencil<lowDimIdx>> emptyStencil_;
718
720 std::shared_ptr<GlueType> glue_;
721
723 int integrationOrder_ = 1;
724
726 std::tuple<
727 std::deque<std::vector<ElementSeed<bulkIdx>>>,
728 std::deque<std::vector<ElementSeed<lowDimIdx>>>
729 > elementSets_;
730};
731
732} // end namespace Dumux
733
734#endif
Coloring schemes for shared-memory-parallel assembly.
Parallel for loop (multithreading)
A helper to deduce a vector with the same size as numbers of equations.
The available discretization methods in Dumux.
Formatting based on the fmt-library which implements std::format of C++20.
Algorithms that finds which geometric entities intersect.
Helper functions for distance queries.
An integration point source class, i.e. sources located at a single point in space associated with a ...
Data associated with a point source.
A class glueing two grids of potentially different dimension geometrically. Intersections are compute...
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
static Geometry::ctype averageDistancePointGeometry(const typename Geometry::GlobalCoordinate &p, const Geometry &geometry, std::size_t integrationOrder=2)
Compute the average distance from a point to a geometry by integration.
Definition: distance.hh:41
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:172
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:38
A vector of primary variables.
Definition: common/properties.hh:49
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
Definition: common/properties.hh:100
A class representing the intersection entities two geometric entity sets.
Definition: intersectionentityset.hh:55
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:299
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:321
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:350
void updateSolution(const SolutionVector &curSol)
Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want t...
Definition: multidomain/couplingmanager.hh:231
the default point source traits
Definition: couplingmanagerbase.hh:57
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:83
std::vector< GridIndex< i > > & vertexIndices(Dune::index_constant< i > dom, GridIndex< i > eIdx)
Return a reference to the vertex indices.
Definition: couplingmanagerbase.hh:684
Scalar averageDistance(std::size_t id) const
return the average distance to the coupled bulk cell center
Definition: couplingmanagerbase.hh:397
Implementation & asImp_()
Returns the implementation of the problem (i.e. static polymorphism)
Definition: couplingmanagerbase.hh:696
const std::vector< PointSource< i > > & pointSources(Dune::index_constant< i > dom) const
Return the point source if domain i.
Definition: couplingmanagerbase.hh:410
const CouplingStencil< bulkIdx > & extendedSourceStencil(std::size_t eIdx) const
Extended source stencil (for the bulk domain)
Definition: couplingmanagerbase.hh:600
PSTraits PointSourceTraits
export the point source traits
Definition: couplingmanagerbase.hh:119
EmbeddedCouplingManagerBase(std::shared_ptr< const GridGeometry< bulkIdx > > bulkGridGeometry, std::shared_ptr< const GridGeometry< lowDimIdx > > lowDimGridGeometry)
Constructor.
Definition: couplingmanagerbase.hh:135
const CouplingStencils< i > & couplingStencils(Dune::index_constant< i > dom) const
Return reference to bulk coupling stencil member of domain i.
Definition: couplingmanagerbase.hh:415
void updateAfterGridAdaption(std::shared_ptr< const GridGeometry< bulkIdx > > bulkGridGeometry, std::shared_ptr< const GridGeometry< lowDimIdx > > lowDimGridGeometry)
call this after grid adaption
Definition: couplingmanagerbase.hh:126
const GridView< id > & gridView(Dune::index_constant< id > domainIdx) const
Return a reference to the bulk problem.
Definition: couplingmanagerbase.hh:385
const Implementation & asImp_() const
Returns the implementation of the problem (i.e. static polymorphism)
Definition: couplingmanagerbase.hh:700
std::vector< Scalar > & averageDistanceToBulkCell()
Return reference to average distances to bulk cell.
Definition: couplingmanagerbase.hh:669
const auto & curSol(Dune::index_constant< i > domainIdx) const
the solution vector of the subproblem
Definition: couplingmanagerbase.hh:442
const CouplingStencil< i > & emptyStencil(Dune::index_constant< i > dom) const
Return a reference to an empty stencil.
Definition: couplingmanagerbase.hh:424
const std::vector< PointSource< bulkIdx > > & bulkPointSources() const
Return reference to bulk point sources.
Definition: couplingmanagerbase.hh:401
auto & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: couplingmanagerbase.hh:433
const PointSourceData & pointSourceData(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanagerbase.hh:380
std::vector< std::vector< GridIndex< i > > > & vertexIndices(Dune::index_constant< i > dom)
Return a reference to the vertex indices container.
Definition: couplingmanagerbase.hh:689
void glueGrids()
compute the intersections between the two grids
Definition: couplingmanagerbase.hh:655
MDTraits MultiDomainTraits
export traits
Definition: couplingmanagerbase.hh:117
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanagerbase.hh:240
std::vector< PointSourceData > & pointSourceData()
Return reference to point source data vector member.
Definition: couplingmanagerbase.hh:665
const CouplingStencil< j > & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &element, Dune::index_constant< j > domainJ) const
Methods to be accessed by the assembly.
Definition: couplingmanagerbase.hh:179
PrimaryVariables< lowDimIdx > lowDimPriVars(std::size_t id) const
Return data for a low dim point source with the identifier id.
Definition: couplingmanagerbase.hh:393
void getShapeValues(Dune::index_constant< i > domainI, const FVGG &gridGeometry, const Geometry &geo, const GlobalPosition &globalPos, ShapeValues &shapeValues)
compute the shape function for a given point and geometry
Definition: couplingmanagerbase.hh:627
const std::vector< PointSource< lowDimIdx > > & lowDimPointSources() const
Return reference to low dim point sources.
Definition: couplingmanagerbase.hh:405
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanagerbase.hh:575
void clear()
Clear all internal data members.
Definition: couplingmanagerbase.hh:640
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanagerbase.hh:448
const GlueType & glue() const
Definition: couplingmanagerbase.hh:692
decltype(auto) evalCouplingResidual(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ)
evaluates the element residual of a coupled element of domain i which depends on the variables at the...
Definition: couplingmanagerbase.hh:205
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Methods to be accessed by main.
Definition: couplingmanagerbase.hh:146
CouplingStencils< i > & couplingStencils(Dune::index_constant< i > dom)
Return reference to bulk coupling stencil member of domain i.
Definition: couplingmanagerbase.hh:679
std::vector< PointSource< i > > & pointSources(Dune::index_constant< i > dom)
Return the point source if domain i.
Definition: couplingmanagerbase.hh:674
std::unordered_map< GridIndex< id >, CouplingStencil< id > > CouplingStencils
export stencil types
Definition: couplingmanagerbase.hh:121
PrimaryVariables< bulkIdx > bulkPriVars(std::size_t id) const
Return data for a bulk point source with the identifier id.
Definition: couplingmanagerbase.hh:389
void precomputeVertexIndices(Dune::index_constant< id > domainIdx)
computes the vertex indices per element for the box method
Definition: couplingmanagerbase.hh:608
std::size_t idCounter_
id generator for point sources
Definition: couplingmanagerbase.hh:704
const std::vector< PointSourceData > & pointSourceData() const
Return reference to point source data vector member.
Definition: couplingmanagerbase.hh:419
An integration point source class with an identifier to attach data and a quadrature weight and integ...
Definition: integrationpointsource.hh:45
A helper class calculating a DOF-index to point source map.
Definition: integrationpointsource.hh:112
A point source data class used for integration in multidimensional models.
Definition: pointsourcedata.hh:43
A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
Declares all properties used in Dumux.
The interface of the coupling manager for multi domain problems.