3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
couplingmanager1d3d.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_COUPLINGMANAGER_1D3D_HH
27#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_HH
28
29#include <vector>
30
31#include <dune/common/timer.hh>
32#include <dune/geometry/quadraturerules.hh>
33
40
41namespace Dumux {
42
49
51template<class MDTraits>
53{
54private:
55 template<std::size_t i> using SubDomainTypeTag = typename MDTraits::template SubDomain<i>::TypeTag;
56 template<std::size_t i> using GridGeometry = GetPropType<SubDomainTypeTag<i>, Properties::GridGeometry>;
57 template<std::size_t i> using NumEqVector = GetPropType<SubDomainTypeTag<i>, Properties::NumEqVector>;
58public:
60 template<std::size_t i>
62
64 template<std::size_t i>
66
69};
70
76template<class MDTraits, EmbeddedCouplingMode mode>
78
87template<class MDTraits>
89: public EmbeddedCouplingManagerBase<MDTraits, EmbeddedCouplingManager1d3d<MDTraits, EmbeddedCouplingMode::line>>
90{
93
94 using Scalar = typename MDTraits::Scalar;
95 using SolutionVector = typename MDTraits::SolutionVector;
96
97 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
98 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
99
100 // the sub domain type aliases
101 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
102 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
103 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
104 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
105 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
106
107public:
108 static constexpr EmbeddedCouplingMode couplingMode = EmbeddedCouplingMode::line;
109
111
112 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
113 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
114 const SolutionVector& curSol)
115 {
116 ParentType::init(bulkProblem, lowDimProblem, curSol);
117 computeLowDimVolumeFractions();
118 }
119
122 {
123 // resize the storage vector
124 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
125
126 // compute the low dim volume fractions
127 for (const auto& is : intersections(this->glue()))
128 {
129 // all inside elements are identical...
130 const auto& inside = is.targetEntity(0);
131 const auto intersectionGeometry = is.geometry();
132 const std::size_t lowDimElementIdx = this->problem(lowDimIdx).gridGeometry().elementMapper().index(inside);
133
134 // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional
135 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
136 for (std::size_t outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
137 {
138 const auto& outside = is.domainEntity(outsideIdx);
139 const std::size_t bulkElementIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(outside);
140 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
141 }
142 }
143 }
144
148 // \{
149
151 Scalar radius(std::size_t id) const
152 {
153 const auto& data = this->pointSourceData()[id];
154 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
155 }
156
158 // For one-dimensional low dim domain we assume radial tubes
159 Scalar lowDimVolume(const Element<bulkIdx>& element) const
160 {
161 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
162 return lowDimVolumeInBulkElement_[eIdx];
163 }
164
166 // For one-dimensional low dim domain we assume radial tubes
167 Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const
168 {
169 const auto totalVolume = element.geometry().volume();
170 return lowDimVolume(element) / totalVolume;
171 }
172
173 // \}
174
175private:
177 std::vector<Scalar> lowDimVolumeInBulkElement_;
178};
179
186template<class MDTraits>
188: public EmbeddedCouplingManagerBase<MDTraits, EmbeddedCouplingManager1d3d<MDTraits, EmbeddedCouplingMode::average>,
189 CircleAveragePointSourceTraits<MDTraits>>
190{
193 using Scalar = typename MDTraits::Scalar;
194 using SolutionVector = typename MDTraits::SolutionVector;
195 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
196
197 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
198 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
199
200 // the sub domain type aliases
201 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
202 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
203 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
204 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
205 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
206
207 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
208
209 template<std::size_t id>
210 static constexpr bool isBox()
211 { return GridGeometry<id>::discMethod == DiscretizationMethod::box; }
212
213
214public:
215 enum {
216 bulkDim = GridView<bulkIdx>::dimension,
217 lowDimDim = GridView<lowDimIdx>::dimension,
218 dimWorld = GridView<bulkIdx>::dimensionworld
219 };
220
222
224
225 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
226 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
227 const SolutionVector& curSol)
228 {
229 ParentType::init(bulkProblem, lowDimProblem, curSol);
230 computeLowDimVolumeFractions();
231 }
232
239 template<std::size_t id, class JacobianPattern>
240 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
241 {
242 extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern);
243 }
244
252 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
253 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
254 const LocalAssemblerI& localAssemblerI,
255 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
256 JacobianMatrixDiagBlock& A,
257 GridVariables& gridVariables)
258 {
259 extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, this->curSol(), A, gridVariables);
260 }
261
262 /* \brief Compute integration point point sources and associated data
263 *
264 * This method uses grid glue to intersect the given grids. Over each intersection
265 * we later need to integrate a source term. This method places point sources
266 * at each quadrature point and provides the point source with the necessary
267 * information to compute integrals (quadrature weight and integration element)
268 * \param order The order of the quadrature rule for integration of sources over an intersection
269 * \param verbose If the point source computation is verbose
270 */
271 void computePointSourceData(std::size_t order = 1, bool verbose = false)
272 {
273 // Initialize the bulk bounding box tree
274 const auto& bulkTree = this->problem(bulkIdx).gridGeometry().boundingBoxTree();
275
276 // initilize the maps
277 // do some logging and profiling
278 Dune::Timer watch;
279 std::cout << "Initializing the point sources..." << std::endl;
280
281 // clear all internal members like pointsource vectors and stencils
282 // initializes the point source id counter
283 this->clear();
284 extendedSourceStencil_.stencil().clear();
285
286 // precompute the vertex indices for efficiency
287 this->preComputeVertexIndices(bulkIdx);
288 this->preComputeVertexIndices(lowDimIdx);
289
290 // iterate over all lowdim elements
291 const auto& lowDimProblem = this->problem(lowDimIdx);
292 for (const auto& lowDimElement : elements(this->gridView(lowDimIdx)))
293 {
294 // get the Gaussian quadrature rule for the low dim element
295 const auto lowDimGeometry = lowDimElement.geometry();
296 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(lowDimGeometry.type(), order);
297
298 const auto lowDimElementIdx = lowDimProblem.gridGeometry().elementMapper().index(lowDimElement);
299
300 // apply the Gaussian quadrature rule and define point sources at each quadrature point
301 // note that the approximation is not optimal if
302 // (a) the one-dimensional elements are too large,
303 // (b) whenever a one-dimensional element is split between two or more elements,
304 // (c) when gradients of important quantities in the three-dimensional domain are large.
305
306 // iterate over all quadrature points
307 for (auto&& qp : quad)
308 {
309 // global position of the quadrature point
310 const auto globalPos = lowDimGeometry.global(qp.position());
311
312 const auto bulkElementIndices = intersectingEntities(globalPos, bulkTree);
313
314 // do not add a point source if the qp is outside of the 3d grid
315 // this is equivalent to having a source of zero for that qp
316 if (bulkElementIndices.empty())
317 continue;
318
320 // get circle average connectivity and interpolation data
322
323 static const auto numIp = getParam<int>("MixedDimension.NumCircleSegments");
324 const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx);
325 const auto normal = lowDimGeometry.corner(1)-lowDimGeometry.corner(0);
326 const auto weight = 2*M_PI*radius/numIp;
327
328 const auto circlePoints = EmbeddedCoupling::circlePoints(globalPos, normal, radius, numIp);
329 std::vector<Scalar> circleIpWeight(circlePoints.size());
330 std::vector<std::size_t> circleStencil(circlePoints.size());
331 // for box
332 std::unordered_map<std::size_t, std::vector<std::size_t> > circleCornerIndices;
333 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
334 std::unordered_map<std::size_t, ShapeValues> circleShapeValues;
335
336 for (int k = 0; k < circlePoints.size(); ++k)
337 {
338 const auto circleBulkElementIndices = intersectingEntities(circlePoints[k], bulkTree);
339 if (circleBulkElementIndices.empty())
340 continue;
341
342 const auto bulkElementIdx = circleBulkElementIndices[0];
343 circleStencil[k] = bulkElementIdx;
344 circleIpWeight[k] = weight;
345
346 if (isBox<bulkIdx>())
347 {
348 if (!static_cast<bool>(circleCornerIndices.count(bulkElementIdx)))
349 {
350 const auto bulkElement = this->problem(bulkIdx).gridGeometry().element(bulkElementIdx);
351 circleCornerIndices[bulkElementIdx] = this->vertexIndices(bulkIdx, bulkElementIdx);
352
353 // evaluate shape functions at the integration point
354 const auto bulkGeometry = bulkElement.geometry();
355 this->getShapeValues(bulkIdx, this->problem(bulkIdx).gridGeometry(), bulkGeometry, circlePoints[k], circleShapeValues[bulkElementIdx]);
356 }
357 }
358 }
359
360 // export low dim circle stencil
361 if (isBox<bulkIdx>())
362 {
363 // we insert all vertices and make it unique later
364 for (const auto& vertices : circleCornerIndices)
365 {
366 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
367 vertices.second.begin(), vertices.second.end());
368
369 }
370 }
371 else
372 {
373 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
374 circleStencil.begin(), circleStencil.end());
375 }
376
377 // loop over the bulk elements at the integration points (usually one except when it is on a face or edge or vertex)
378 for (auto bulkElementIdx : bulkElementIndices)
379 {
380 const auto id = this->idCounter_++;
381 const auto ie = lowDimGeometry.integrationElement(qp.position());
382 const auto qpweight = qp.weight();
383
384 this->pointSources(bulkIdx).emplace_back(globalPos, id, qpweight, ie, std::vector<std::size_t>({bulkElementIdx}));
385 this->pointSources(bulkIdx).back().setEmbeddings(bulkElementIndices.size());
386 this->pointSources(lowDimIdx).emplace_back(globalPos, id, qpweight, ie, std::vector<std::size_t>({lowDimElementIdx}));
387 this->pointSources(lowDimIdx).back().setEmbeddings(bulkElementIndices.size());
388
389 // pre compute additional data used for the evaluation of
390 // the actual solution dependent source term
391 PointSourceData psData;
392
393 if (isBox<lowDimIdx>())
394 {
395 ShapeValues shapeValues;
396 this->getShapeValues(lowDimIdx, this->problem(lowDimIdx).gridGeometry(), lowDimGeometry, globalPos, shapeValues);
397 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
398 }
399 else
400 {
401 psData.addLowDimInterpolation(lowDimElementIdx);
402 }
403
404 // add data needed to compute integral over the circle
405 if (isBox<bulkIdx>())
406 {
407 psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
408
409 const auto bulkGeometry = this->problem(bulkIdx).gridGeometry().element(bulkElementIdx).geometry();
410 ShapeValues shapeValues;
411 this->getShapeValues(bulkIdx, this->problem(bulkIdx).gridGeometry(), bulkGeometry, globalPos, shapeValues);
412 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
413 }
414 else
415 {
416 psData.addCircleInterpolation(circleIpWeight, circleStencil);
417 psData.addBulkInterpolation(bulkElementIdx);
418 }
419
420 // publish point source data in the global vector
421 this->pointSourceData().emplace_back(std::move(psData));
422
423 // export the bulk coupling stencil
424 if (isBox<lowDimIdx>())
425 {
426 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
427 this->vertexIndices(lowDimIdx, lowDimElementIdx).begin(),
428 this->vertexIndices(lowDimIdx, lowDimElementIdx).end());
429
430 }
431 else
432 {
433 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
434 }
435
436 // export bulk circle stencil
437 if (isBox<bulkIdx>())
438 {
439 // we insert all vertices and make it unique later
440 for (const auto& vertices : circleCornerIndices)
441 {
442 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
443 vertices.second.begin(), vertices.second.end());
444
445 }
446 }
447 else
448 {
449 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
450 circleStencil.begin(), circleStencil.end());
451 }
452 }
453 }
454 }
455
456 // make the circle stencil unique (for source derivatives)
457 for (auto&& stencil : extendedSourceStencil_.stencil())
458 {
459 std::sort(stencil.second.begin(), stencil.second.end());
460 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
461
462 // remove the vertices element (box)
463 if (isBox<bulkIdx>())
464 {
465 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
466 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
467 [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
468 stencil.second.end());
469 }
470 // remove the own element (cell-centered)
471 else
472 {
473 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
474 [&](auto i){ return i == stencil.first; }),
475 stencil.second.end());
476 }
477 }
478
479 // make stencils unique
480 using namespace Dune::Hybrid;
481 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
482 {
483 for (auto&& stencil : this->couplingStencils(domainIdx))
484 {
485 std::sort(stencil.second.begin(), stencil.second.end());
486 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
487 }
488 });
489
490 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
491 }
492
495 {
496 // resize the storage vector
497 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
498
499 // compute the low dim volume fractions
500 for (const auto& is : intersections(this->glue()))
501 {
502 // all inside elements are identical...
503 const auto& inside = is.targetEntity(0);
504 const auto intersectionGeometry = is.geometry();
505 const std::size_t lowDimElementIdx = this->problem(lowDimIdx).gridGeometry().elementMapper().index(inside);
506
507 // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional
508 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
509 for (std::size_t outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
510 {
511 const auto& outside = is.domainEntity(outsideIdx);
512 const std::size_t bulkElementIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(outside);
513 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
514 }
515 }
516 }
517
521 // \{
522
524 Scalar radius(std::size_t id) const
525 {
526 const auto& data = this->pointSourceData()[id];
527 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
528 }
529
531 // For one-dimensional low dim domain we assume radial tubes
532 Scalar lowDimVolume(const Element<bulkIdx>& element) const
533 {
534 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
535 return lowDimVolumeInBulkElement_[eIdx];
536 }
537
539 // For one-dimensional low dim domain we assume radial tubes
540 Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const
541 {
542 const auto totalVolume = element.geometry().volume();
543 return lowDimVolume(element) / totalVolume;
544 }
545
546 // \}
547
548private:
551
553 std::vector<Scalar> lowDimVolumeInBulkElement_;
554};
555
556
564template<class MDTraits>
566: public EmbeddedCouplingManagerBase<MDTraits, EmbeddedCouplingManager1d3d<MDTraits, EmbeddedCouplingMode::cylindersources>>
567{
570 using Scalar = typename MDTraits::Scalar;
571 using SolutionVector = typename MDTraits::SolutionVector;
572 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
573
574 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
575 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
576
577 // the sub domain type aliases
578 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
579 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
580 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
581 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
582 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
583
584 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
585
586 template<std::size_t id>
587 static constexpr bool isBox()
588 { return GridGeometry<id>::discMethod == DiscretizationMethod::box; }
589
590 enum {
591 bulkDim = GridView<bulkIdx>::dimension,
592 lowDimDim = GridView<lowDimIdx>::dimension,
593 dimWorld = GridView<bulkIdx>::dimensionworld
594 };
595public:
597
599
600 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
601 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
602 const SolutionVector& curSol)
603 {
604 ParentType::init(bulkProblem, lowDimProblem, curSol);
605 computeLowDimVolumeFractions();
606 }
607
608 /* \brief Compute integration point point sources and associated data
609 *
610 * This method uses grid glue to intersect the given grids. Over each intersection
611 * we later need to integrate a source term. This method places point sources
612 * at each quadrature point and provides the point source with the necessary
613 * information to compute integrals (quadrature weight and integration element)
614 * \param order The order of the quadrature rule for integration of sources over an intersection
615 * \param verbose If the point source computation is verbose
616 */
617 void computePointSourceData(std::size_t order = 1, bool verbose = false)
618 {
619 // Initialize the bulk bounding box tree
620 const auto& bulkTree = this->problem(bulkIdx).gridGeometry().boundingBoxTree();
621
622 // initilize the maps
623 // do some logging and profiling
624 Dune::Timer watch;
625 std::cout << "Initializing the point sources..." << std::endl;
626
627 // clear all internal members like pointsource vectors and stencils
628 // initializes the point source id counter
629 this->clear();
630
631 // precompute the vertex indices for efficiency
632 this->preComputeVertexIndices(bulkIdx);
633 this->preComputeVertexIndices(lowDimIdx);
634
635 // iterate over all lowdim elements
636 const auto& lowDimProblem = this->problem(lowDimIdx);
637 for (const auto& lowDimElement : elements(this->gridView(lowDimIdx)))
638 {
639 // get the Gaussian quadrature rule for the low dim element
640 const auto lowDimGeometry = lowDimElement.geometry();
641 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(lowDimGeometry.type(), order);
642
643 const auto lowDimElementIdx = lowDimProblem.gridGeometry().elementMapper().index(lowDimElement);
644
645 // apply the Gaussian quadrature rule and define point sources at each quadrature point
646 // note that the approximation is not optimal if
647 // (a) the one-dimensional elements are too large,
648 // (b) whenever a one-dimensional element is split between two or more elements,
649 // (c) when gradients of important quantities in the three-dimensional domain are large.
650
651 // iterate over all quadrature points
652 for (auto&& qp : quad)
653 {
654 // global position of the quadrature point
655 const auto globalPos = lowDimGeometry.global(qp.position());
656
657 const auto bulkElementIndices = intersectingEntities(globalPos, bulkTree);
658
659 // do not add a point source if the qp is outside of the 3d grid
660 // this is equivalent to having a source of zero for that qp
661 if (bulkElementIndices.empty())
662 continue;
663
665 // get points on the cylinder surface at the integration point
667
668 static const auto numIp = getParam<int>("MixedDimension.NumCircleSegments", 25);
669 const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx);
670 const auto normal = lowDimGeometry.corner(1)-lowDimGeometry.corner(0);
671 const auto integrationElement = lowDimGeometry.integrationElement(qp.position())*2*M_PI*radius/Scalar(numIp);
672 const auto weight = qp.weight()/(2*M_PI*radius);
673
674 const auto circlePoints = EmbeddedCoupling::circlePoints(globalPos, normal, radius, numIp);
675 for (int k = 0; k < circlePoints.size(); ++k)
676 {
677 const auto& circlePos = circlePoints[k];
678 const auto circleBulkElementIndices = intersectingEntities(circlePos, bulkTree);
679 if (circleBulkElementIndices.empty())
680 continue;
681
682 // loop over the bulk elements at the integration points (usually one except when it is on a face or edge or vertex)
683 // and add a point source at every point on the circle
684 for (const auto bulkElementIdx : circleBulkElementIndices)
685 {
686 const auto id = this->idCounter_++;
687 this->pointSources(bulkIdx).emplace_back(circlePos, id, weight, integrationElement, std::vector<std::size_t>({bulkElementIdx}));
688 this->pointSources(bulkIdx).back().setEmbeddings(circleBulkElementIndices.size());
689 this->pointSources(lowDimIdx).emplace_back(globalPos, id, weight, integrationElement, std::vector<std::size_t>({lowDimElementIdx}));
690 this->pointSources(lowDimIdx).back().setEmbeddings(circleBulkElementIndices.size());
691
692 // pre compute additional data used for the evaluation of
693 // the actual solution dependent source term
694 PointSourceData psData;
695
696 if (isBox<lowDimIdx>())
697 {
698 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
699 ShapeValues shapeValues;
700 this->getShapeValues(lowDimIdx, this->problem(lowDimIdx).gridGeometry(), lowDimGeometry, globalPos, shapeValues);
701 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
702 }
703 else
704 {
705 psData.addLowDimInterpolation(lowDimElementIdx);
706 }
707
708 // add data needed to compute integral over the circle
709 if (isBox<bulkIdx>())
710 {
711 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
712 const auto bulkGeometry = this->problem(bulkIdx).gridGeometry().element(bulkElementIdx).geometry();
713 ShapeValues shapeValues;
714 this->getShapeValues(bulkIdx, this->problem(bulkIdx).gridGeometry(), bulkGeometry, circlePos, shapeValues);
715 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
716 }
717 else
718 {
719 psData.addBulkInterpolation(bulkElementIdx);
720 }
721
722 // publish point source data in the global vector
723 this->pointSourceData().emplace_back(std::move(psData));
724
725 // export the lowdim coupling stencil
726 // we insert all vertices / elements and make it unique later
727 if (isBox<bulkIdx>())
728 {
729 const auto& vertices = this->vertexIndices(bulkIdx, bulkElementIdx);
730 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
731 vertices.begin(), vertices.end());
732 }
733 else
734 {
735 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
736 }
737
738 // export the bulk coupling stencil
739 // we insert all vertices / elements and make it unique later
740 if (isBox<lowDimIdx>())
741 {
742 const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx);
743 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
744 vertices.begin(), vertices.end());
745
746 }
747 else
748 {
749 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
750 }
751
752 }
753 }
754 }
755 }
756
757 // make stencils unique
758 using namespace Dune::Hybrid;
759 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
760 {
761 for (auto&& stencil : this->couplingStencils(domainIdx))
762 {
763 std::sort(stencil.second.begin(), stencil.second.end());
764 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
765 }
766 });
767
768 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
769 }
770
773 {
774 // resize the storage vector
775 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
776
777 // compute the low dim volume fractions
778 for (const auto& is : intersections(this->glue()))
779 {
780 // all inside elements are identical...
781 const auto& inside = is.targetEntity(0);
782 const auto intersectionGeometry = is.geometry();
783 const std::size_t lowDimElementIdx = this->problem(lowDimIdx).gridGeometry().elementMapper().index(inside);
784
785 // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional
786 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
787 for (std::size_t outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
788 {
789 const auto& outside = is.domainEntity(outsideIdx);
790 const std::size_t bulkElementIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(outside);
791 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
792 }
793 }
794 }
795
799 // \{
800
802 Scalar radius(std::size_t id) const
803 {
804 const auto& data = this->pointSourceData()[id];
805 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
806 }
807
809 // For one-dimensional low dim domain we assume radial tubes
810 Scalar lowDimVolume(const Element<bulkIdx>& element) const
811 {
812 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
813 return lowDimVolumeInBulkElement_[eIdx];
814 }
815
817 // For one-dimensional low dim domain we assume radial tubes
818 Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const
819 {
820 const auto totalVolume = element.geometry().volume();
821 return lowDimVolume(element) / totalVolume;
822 }
823
824 // \}
825
826private:
827
829 std::vector<Scalar> lowDimVolumeInBulkElement_;
830};
831
832
840template<class MDTraits>
842: public EmbeddedCouplingManagerBase<MDTraits, EmbeddedCouplingManager1d3d<MDTraits, EmbeddedCouplingMode::kernel>>
843{
846 using Scalar = typename MDTraits::Scalar;
847 using SolutionVector = typename MDTraits::SolutionVector;
848 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
849
850 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
851 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
852
853 // the sub domain type aliases
854 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
855 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
856 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
857 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
858 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
859
860 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
861
862 template<std::size_t id>
863 static constexpr bool isBox()
864 { return GridGeometry<id>::discMethod == DiscretizationMethod::box; }
865
866 enum {
867 bulkDim = GridView<bulkIdx>::dimension,
868 lowDimDim = GridView<lowDimIdx>::dimension,
869 dimWorld = GridView<bulkIdx>::dimensionworld
870 };
871public:
872 static constexpr EmbeddedCouplingMode couplingMode = EmbeddedCouplingMode::kernel;
873
875
876 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
877 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
878 const SolutionVector& curSol)
879 {
880 ParentType::init(bulkProblem, lowDimProblem, curSol);
881 computeLowDimVolumeFractions();
882 }
883
890 template<std::size_t id, class JacobianPattern>
891 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
892 {
893 extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern);
894 }
895
903 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
904 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
905 const LocalAssemblerI& localAssemblerI,
906 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
907 JacobianMatrixDiagBlock& A,
908 GridVariables& gridVariables)
909 {
910 extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, this->curSol(), A, gridVariables);
911 }
912
913 /* \brief Compute integration point point sources and associated data
914 *
915 * This method uses grid glue to intersect the given grids. Over each intersection
916 * we later need to integrate a source term. This method places point sources
917 * at each quadrature point and provides the point source with the necessary
918 * information to compute integrals (quadrature weight and integration element)
919 * \param order The order of the quadrature rule for integration of sources over an intersection
920 * \param verbose If the point source computation is verbose
921 */
922 void computePointSourceData(std::size_t order = 1, bool verbose = false)
923 {
924 // initilize the maps
925 // do some logging and profiling
926 Dune::Timer watch;
927 std::cout << "Initializing the point sources..." << std::endl;
928
929 // clear all internal members like pointsource vectors and stencils
930 // initializes the point source id counter
931 this->clear();
932 bulkSourceIds_.clear();
933 bulkSourceWeights_.clear();
934 extendedSourceStencil_.stencil().clear();
935
936 // precompute the vertex indices for efficiency for the box method
937 this->preComputeVertexIndices(bulkIdx);
938 this->preComputeVertexIndices(lowDimIdx);
939
940 const auto& bulkFvGridGeometry = this->problem(bulkIdx).gridGeometry();
941 const auto& lowDimFvGridGeometry = this->problem(lowDimIdx).gridGeometry();
942
943 bulkSourceIds_.resize(this->gridView(bulkIdx).size(0));
944 bulkSourceWeights_.resize(this->gridView(bulkIdx).size(0));
945
946 // intersect the bounding box trees
947 this->glueGrids();
948
949 this->pointSourceData().reserve(this->glue().size());
950 this->averageDistanceToBulkCell().reserve(this->glue().size());
951 const Scalar kernelWidth = getParam<Scalar>("MixedDimension.KernelWidth");
952 for (const auto& is : intersections(this->glue()))
953 {
954 // all inside elements are identical...
955 const auto& inside = is.targetEntity(0);
956 // get the intersection geometry for integrating over it
957 const auto intersectionGeometry = is.geometry();
958
959 // get the Gaussian quadrature rule for the local intersection
960 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order);
961 const std::size_t lowDimElementIdx = lowDimFvGridGeometry.elementMapper().index(inside);
962
963 // iterate over all quadrature points and place a source
964 // for 1d: make a new point source
965 // for 3d: make a new kernel volume source
966 for (auto&& qp : quad)
967 {
968 // compute the coupling stencils
969 for (std::size_t outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
970 {
971 const auto& outside = is.domainEntity(outsideIdx);
972 const std::size_t bulkElementIdx = bulkFvGridGeometry.elementMapper().index(outside);
973
974 // each quadrature point will be a point source for the sub problem
975 const auto globalPos = intersectionGeometry.global(qp.position());
976 const auto id = this->idCounter_++;
977 const auto qpweight = qp.weight();
978 const auto ie = intersectionGeometry.integrationElement(qp.position());
979 this->pointSources(lowDimIdx).emplace_back(globalPos, id, qpweight, ie, std::vector<std::size_t>({lowDimElementIdx}));
980 this->pointSources(lowDimIdx).back().setEmbeddings(is.numDomainNeighbors());
981 computeBulkSource(globalPos, kernelWidth, id, lowDimElementIdx, bulkElementIdx, qpweight*ie/is.numDomainNeighbors());
982
983 // pre compute additional data used for the evaluation of
984 // the actual solution dependent source term
985 PointSourceData psData;
986
987 if (isBox<lowDimIdx>())
988 {
989 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
990 const auto lowDimGeometry = this->problem(lowDimIdx).gridGeometry().element(lowDimElementIdx).geometry();
991 ShapeValues shapeValues;
992 this->getShapeValues(lowDimIdx, this->problem(lowDimIdx).gridGeometry(), lowDimGeometry, globalPos, shapeValues);
993 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
994 }
995 else
996 {
997 psData.addLowDimInterpolation(lowDimElementIdx);
998 }
999
1000 // add data needed to compute integral over the circle
1001 if (isBox<bulkIdx>())
1002 {
1003 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
1004 const auto bulkGeometry = this->problem(bulkIdx).gridGeometry().element(bulkElementIdx).geometry();
1005 ShapeValues shapeValues;
1006 this->getShapeValues(bulkIdx, this->problem(bulkIdx).gridGeometry(), bulkGeometry, globalPos, shapeValues);
1007 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
1008 }
1009 else
1010 {
1011 psData.addBulkInterpolation(bulkElementIdx);
1012 }
1013
1014 // publish point source data in the global vector
1015 this->pointSourceData().emplace_back(std::move(psData));
1016
1017 // compute average distance to bulk cell
1018 this->averageDistanceToBulkCell().push_back(this->computeDistance(outside.geometry(), globalPos));
1019
1020 // export the lowdim coupling stencil
1021 // we insert all vertices / elements and make it unique later
1022 if (isBox<bulkIdx>())
1023 {
1024 const auto& vertices = this->vertexIndices(bulkIdx, bulkElementIdx);
1025 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
1026 vertices.begin(), vertices.end());
1027 }
1028 else
1029 {
1030 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
1031 }
1032 }
1033 }
1034 }
1035
1036 // make extra stencils unique
1037 for (auto&& stencil : extendedSourceStencil_.stencil())
1038 {
1039 std::sort(stencil.second.begin(), stencil.second.end());
1040 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
1041
1042 // remove the vertices element (box)
1043 if (isBox<bulkIdx>())
1044 {
1045 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
1046 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
1047 [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
1048 stencil.second.end());
1049 }
1050 // remove the own element (cell-centered)
1051 else
1052 {
1053 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
1054 [&](auto i){ return i == stencil.first; }),
1055 stencil.second.end());
1056 }
1057 }
1058
1059 // make stencils unique
1060 using namespace Dune::Hybrid;
1061 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
1062 {
1063 for (auto&& stencil : this->couplingStencils(domainIdx))
1064 {
1065 std::sort(stencil.second.begin(), stencil.second.end());
1066 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
1067 }
1068 });
1069
1070 if (!this->pointSources(bulkIdx).empty())
1071 DUNE_THROW(Dune::InvalidStateException, "Kernel method shouldn't have point sources in the bulk domain but only volume sources!");
1072
1073 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
1074 }
1075
1078 {
1079 // resize the storage vector
1080 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
1081
1082 // compute the low dim volume fractions
1083 for (const auto& is : intersections(this->glue()))
1084 {
1085 // all inside elements are identical...
1086 const auto& inside = is.targetEntity(0);
1087 const auto intersectionGeometry = is.geometry();
1088 const std::size_t lowDimElementIdx = this->problem(lowDimIdx).gridGeometry().elementMapper().index(inside);
1089
1090 // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional
1091 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
1092 for (std::size_t outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
1093 {
1094 const auto& outside = is.domainEntity(outsideIdx);
1095 const std::size_t bulkElementIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(outside);
1096 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
1097 }
1098 }
1099 }
1100
1104 // \{
1105
1107 Scalar radius(std::size_t id) const
1108 {
1109 const auto& data = this->pointSourceData()[id];
1110 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
1111 }
1112
1114 // For one-dimensional low dim domain we assume radial tubes
1115 Scalar lowDimVolume(const Element<bulkIdx>& element) const
1116 {
1117 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
1118 return lowDimVolumeInBulkElement_[eIdx];
1119 }
1120
1122 // For one-dimensional low dim domain we assume radial tubes
1123 Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const
1124 {
1125 const auto totalVolume = element.geometry().volume();
1126 return lowDimVolume(element) / totalVolume;
1127 }
1128
1130 const std::vector<std::size_t> bulkSourceIds(std::size_t eIdx) const
1131 { return bulkSourceIds_[eIdx]; }
1132
1134 const std::vector<Scalar> bulkSourceWeights(std::size_t eIdx) const
1135 { return bulkSourceWeights_[eIdx]; }
1136
1137 // \}
1138
1139private:
1141 void computeBulkSource(const GlobalPosition& globalPos, const Scalar kernelWidth,
1142 std::size_t id, std::size_t lowDimElementIdx, std::size_t coupledBulkElementIdx,
1143 Scalar pointSourceWeight)
1144 {
1145 // make sure it is mass conservative
1146 // i.e. the point source in the 1d domain needs to have the exact same integral as the distributed
1147 // kernel source integral in the 3d domain. Correct the integration formula by balancing the error
1148 // by scaling the kernel with the checksum inverse
1149 Scalar checkSum = 0.0;
1150 std::vector<bool> mask(this->gridView(bulkIdx).size(0), false);
1151 for (const auto& element : elements(this->gridView(bulkIdx)))
1152 {
1153 Scalar weight = 0.0;
1154 const auto geometry = element.geometry();
1155 const auto& quad = Dune::QuadratureRules<Scalar, bulkDim>::rule(geometry.type(), 3);
1156 for (auto&& qp : quad)
1157 {
1158 const auto qpweight = qp.weight();
1159 const auto ie = geometry.integrationElement(qp.position());
1160 weight += evalKernel(globalPos, geometry.global(qp.position()), kernelWidth)*qpweight*ie;
1161 }
1162
1163 if (weight > 1e-13)
1164 {
1165 const auto bulkElementIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
1166 bulkSourceIds_[bulkElementIdx].push_back(id);
1167 bulkSourceWeights_[bulkElementIdx].push_back(weight*pointSourceWeight);
1168 mask[bulkElementIdx] = true;
1169
1170 // add lowDim dofs that the source is related to to the bulk stencil
1171 if (isBox<lowDimIdx>())
1172 {
1173 const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx);
1174 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
1175 vertices.begin(), vertices.end());
1176
1177 }
1178 else
1179 {
1180 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
1181 }
1182
1183 // tpfa
1184 extendedSourceStencil_.stencil()[bulkElementIdx].push_back(coupledBulkElementIdx);
1185
1186 // compute check sum -> should sum up to 1.0 to be mass conservative
1187 checkSum += weight;
1188 }
1189 }
1190
1191 for (std::size_t eIdx = 0; eIdx < bulkSourceWeights_.size(); ++eIdx)
1192 if (mask[eIdx]) bulkSourceWeights_[eIdx].back() /= checkSum;
1193
1194 // balance error of the quadrature rule -> TODO: what to do at boundaries
1195 // const auto diff = 1.0 - checkSum;
1196 // std::cout << "Integrated kernel with integration error of " << diff << std::endl;
1197 }
1198
1200 inline Scalar evalKernel(const GlobalPosition& origin,
1201 const GlobalPosition& pos,
1202 const Scalar width) const noexcept
1203 {
1204 const auto r = (pos-origin).two_norm();
1205 const auto r2 = r*r;
1206 const auto r3 = r2*r;
1207
1208 if (r > width)
1209 return 0.0;
1210
1211 const Scalar w2 = width*width;
1212 const Scalar w3 = w2*width;
1213 const Scalar k = 15.0/(4*M_PI*w3);
1214 const Scalar a = 2.0/w3;
1215 const Scalar b = 3.0/w2;
1216
1217 return k*(a*r3 - b*r2 + 1.0);
1218 }
1219
1221 EmbeddedCoupling::ExtendedSourceStencil<ThisType> extendedSourceStencil_;
1223 std::vector<Scalar> lowDimVolumeInBulkElement_;
1225 std::vector<std::vector<std::size_t>> bulkSourceIds_;
1227 std::vector<std::vector<Scalar>> bulkSourceWeights_;
1228};
1229
1230} // end namespace Dumux
1231
1232#endif
Helper function to compute points on a circle.
Extended source stencil helper class for coupling managers.
An integration point source class, i.e. sources located at a single point in space associated with a ...
Data associated with a point source.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
std::vector< std::size_t > intersectingEntities(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
Compute all intersections between entities and a point.
Definition: intersectingentities.hh:99
std::vector< GlobalPosition > circlePoints(const GlobalPosition &center, const GlobalPosition &normal, const typename GlobalPosition::value_type radius, const std::size_t numPoints=20)
returns a vector of points on a circle
Definition: circlepoints.hh:45
EmbeddedCouplingMode
The coupling mode.
Definition: couplingmanager1d3d.hh:48
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
A vector of size number equations that can be used for Neumann fluxes, sources, residuals,...
Definition: common/properties.hh:62
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:70
Definition: common/properties.hh:113
point source traits for the circle average coupling mode
Definition: couplingmanager1d3d.hh:53
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:77
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:90
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d.hh:159
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d.hh:151
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d.hh:167
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d.hh:121
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d.hh:112
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:190
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d.hh:540
void evalAdditionalDomainDerivatives(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const typename LocalAssemblerI::LocalResidual::ElementResidualVector &, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
evaluate additional derivatives of the element residual of a domain with respect to dofs in the same ...
Definition: couplingmanager1d3d.hh:253
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d.hh:494
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d.hh:271
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d.hh:524
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d.hh:532
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d.hh:225
void extendJacobianPattern(Dune::index_constant< id > domainI, JacobianPattern &pattern) const
extend the jacobian pattern of the diagonal block of domain i by those entries that are not already i...
Definition: couplingmanager1d3d.hh:240
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:567
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d.hh:810
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d.hh:818
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d.hh:772
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d.hh:802
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d.hh:617
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d.hh:600
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:843
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d.hh:922
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d.hh:1123
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d.hh:1077
void extendJacobianPattern(Dune::index_constant< id > domainI, JacobianPattern &pattern) const
extend the jacobian pattern of the diagonal block of domain i by those entries that are not already i...
Definition: couplingmanager1d3d.hh:891
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d.hh:1115
void evalAdditionalDomainDerivatives(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const typename LocalAssemblerI::LocalResidual::ElementResidualVector &, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
evaluate additional derivatives of the element residual of a domain with respect to dofs in the same ...
Definition: couplingmanager1d3d.hh:904
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d.hh:1107
const std::vector< Scalar > bulkSourceWeights(std::size_t eIdx) const
return all source ids for a bulk elements
Definition: couplingmanager1d3d.hh:1134
const std::vector< std::size_t > bulkSourceIds(std::size_t eIdx) const
return all source ids for a bulk elements
Definition: couplingmanager1d3d.hh:1130
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d.hh:876
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:77
A class managing an extended source stencil.
Definition: extendedsourcestencil.hh:47
An integration point source class with an identifier to attach data and a quadrature weight and integ...
Definition: integrationpointsource.hh:44
A helper class calculating a DOF-index to point source map.
Definition: integrationpointsource.hh:107
A point source data class used for integration in multidimension models.
Definition: pointsourcedata.hh:149
Declares all properties used in Dumux.