version 3.8
couplingmanager1d3d_average.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_AVERAGE_HH
14#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_AVERAGE_HH
15
16#include <vector>
17
18#include <dune/common/timer.hh>
19#include <dune/geometry/quadraturerules.hh>
20
21#include <dumux/common/tag.hh>
24
26
31
32namespace Dumux {
33
34namespace Embedded1d3dCouplingMode {
35struct Average : public Utility::Tag<Average> {
36 static std::string name() { return "average"; }
37};
38
39inline constexpr Average average{};
40} // end namespace Embedded1d3dCouplingMode
41
42// forward declaration
43template<class MDTraits, class CouplingMode>
44class Embedded1d3dCouplingManager;
45
52template<class MDTraits>
53class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Average>
54: public EmbeddedCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Average>,
55 CircleAveragePointSourceTraits<MDTraits>>
56{
59 using Scalar = typename MDTraits::Scalar;
60 using SolutionVector = typename MDTraits::SolutionVector;
61 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
62
63 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
64 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
65
66 // the sub domain type aliases
67 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
68 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
69 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
70 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
71 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
72 template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex;
73
74 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
75
76 template<std::size_t id>
77 static constexpr bool isBox()
78 { return GridGeometry<id>::discMethod == DiscretizationMethods::box; }
79
80
81public:
82 enum {
83 bulkDim = GridView<bulkIdx>::dimension,
84 lowDimDim = GridView<lowDimIdx>::dimension,
85 dimWorld = GridView<bulkIdx>::dimensionworld
86 };
87
88 static constexpr Embedded1d3dCouplingMode::Average couplingMode{};
89
90 using ParentType::ParentType;
91
92 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
93 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
94 const SolutionVector& curSol)
95 {
96 ParentType::init(bulkProblem, lowDimProblem, curSol);
97 computeLowDimVolumeFractions();
98 }
99
106 template<std::size_t id, class JacobianPattern>
107 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
108 {
109 extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern);
110 }
111
119 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
120 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
121 const LocalAssemblerI& localAssemblerI,
122 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
123 JacobianMatrixDiagBlock& A,
124 GridVariables& gridVariables)
125 {
126 extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, A, gridVariables);
127 }
128
129 /* \brief Compute integration point point sources and associated data
130 *
131 * This method uses grid glue to intersect the given grids. Over each intersection
132 * we later need to integrate a source term. This method places point sources
133 * at each quadrature point and provides the point source with the necessary
134 * information to compute integrals (quadrature weight and integration element)
135 * \param order The order of the quadrature rule for integration of sources over an intersection
136 * \param verbose If the point source computation is verbose
137 */
138 void computePointSourceData(std::size_t order = 1, bool verbose = false)
139 {
140 // Initialize the bulk bounding box tree
141 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
142 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
143 const auto& bulkTree = bulkGridGeometry.boundingBoxTree();
144
145 // initialize the maps
146 // do some logging and profiling
147 Dune::Timer watch;
148 std::cout << "Initializing the point sources..." << std::endl;
149
150 // clear all internal members like pointsource vectors and stencils
151 // initializes the point source id counter
152 this->clear();
153 extendedSourceStencil_.clear();
154
155 // precompute the vertex indices for efficiency
156 this->precomputeVertexIndices(bulkIdx);
157 this->precomputeVertexIndices(lowDimIdx);
158
159 // intersect the bounding box trees
160 this->glueGrids();
161
162 // iterate over all intersection and add point sources
163 const auto& lowDimProblem = this->problem(lowDimIdx);
164 for (const auto& is : intersections(this->glue()))
165 {
166 // all inside elements are identical...
167 const auto& lowDimElement = is.targetEntity(0);
168 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(lowDimElement);
169
170 // get the intersection geometry
171 const auto intersectionGeometry = is.geometry();
172 // get the Gaussian quadrature rule for the local intersection
173 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order);
174
175 // apply the Gaussian quadrature rule and define point sources at each quadrature point
176 // note that the approximation is not optimal if
177 // (a) the one-dimensional elements are too large,
178 // (b) whenever a one-dimensional element is split between two or more elements,
179 // (c) when gradients of important quantities in the three-dimensional domain are large.
180
181 // iterate over all quadrature points
182 for (auto&& qp : quad)
183 {
184 // global position of the quadrature point
185 const auto globalPos = intersectionGeometry.global(qp.position());
186
187 const auto bulkElementIndices = intersectingEntities(globalPos, bulkTree);
188
189 // do not add a point source if the qp is outside of the 3d grid
190 // this is equivalent to having a source of zero for that qp
191 if (bulkElementIndices.empty())
192 continue;
193
195 // get circle average connectivity and interpolation data
197
198 static const auto numIp = getParam<int>("MixedDimension.NumCircleSegments");
199 const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx);
200 const auto normal = intersectionGeometry.corner(1)-intersectionGeometry.corner(0);
201 const auto circleAvgWeight = 2*M_PI*radius/numIp;
202 const auto integrationElement = intersectionGeometry.integrationElement(qp.position());
203 const auto qpweight = qp.weight();
204
205 const auto circlePoints = EmbeddedCoupling::circlePoints(globalPos, normal, radius, numIp);
206 std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(circlePoints.size());
207 std::vector<GridIndex<bulkIdx>> circleStencil; circleStencil.reserve(circlePoints.size());
208
209 // for box
210 std::vector<const std::vector<GridIndex<bulkIdx>>*> circleCornerIndices;
211 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
212 std::vector<ShapeValues> circleShapeValues;
213
214 // go over all points of the average operator
215 int insideCirclePoints = 0;
216 for (int k = 0; k < circlePoints.size(); ++k)
217 {
218 const auto circleBulkElementIndices = intersectingEntities(circlePoints[k], bulkTree);
219 if (circleBulkElementIndices.empty())
220 continue;
221
222 ++insideCirclePoints;
223 const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices.size();
224 for (const auto bulkElementIdx : circleBulkElementIndices)
225 {
226 circleStencil.push_back(bulkElementIdx);
227 circleIpWeight.push_back(localCircleAvgWeight);
228
229 // precompute interpolation data for box scheme for each cut bulk element
230 if constexpr (isBox<bulkIdx>())
231 {
232 const auto bulkElement = bulkGridGeometry.element(bulkElementIdx);
233 circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx)));
234
235 // evaluate shape functions at the integration point
236 const auto bulkGeometry = bulkElement.geometry();
237 ShapeValues shapeValues;
238 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePoints[k], shapeValues);
239 circleShapeValues.emplace_back(std::move(shapeValues));
240 }
241 }
242 }
243
244 // if the circle stencil is empty (that is the circle is entirely outside the domain)
245 // we do not add a (zero) point source
246 if (circleStencil.empty())
247 continue;
248
249 // export low dim circle stencil
250 if constexpr (isBox<bulkIdx>())
251 {
252 // we insert all vertices and make it unique later
253 for (const auto& vertices : circleCornerIndices)
254 {
255 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
256 vertices->begin(), vertices->end());
257
258 }
259 }
260 else
261 {
262 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
263 circleStencil.begin(), circleStencil.end());
264 }
265
266 // surface fraction that is inside the domain (only different than 1.0 on the boundary)
267 const auto surfaceFraction = Scalar(insideCirclePoints)/Scalar(circlePoints.size());
268
269 // loop over the bulk elements at the integration points (usually one except when it is on a face or edge or vertex)
270 for (auto bulkElementIdx : bulkElementIndices)
271 {
272 const auto id = this->idCounter_++;
273
274 this->pointSources(bulkIdx).emplace_back(globalPos, id, qpweight, integrationElement*surfaceFraction, bulkElementIdx);
275 this->pointSources(bulkIdx).back().setEmbeddings(bulkElementIndices.size());
276 this->pointSources(lowDimIdx).emplace_back(globalPos, id, qpweight, integrationElement*surfaceFraction, lowDimElementIdx);
277 this->pointSources(lowDimIdx).back().setEmbeddings(bulkElementIndices.size());
278
279 // pre compute additional data used for the evaluation of
280 // the actual solution dependent source term
281 PointSourceData psData;
282
283 if constexpr (isBox<lowDimIdx>())
284 {
285 ShapeValues shapeValues;
286 this->getShapeValues(lowDimIdx, lowDimGridGeometry, intersectionGeometry, globalPos, shapeValues);
287 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
288 }
289 else
290 {
291 psData.addLowDimInterpolation(lowDimElementIdx);
292 }
293
294 // add data needed to compute integral over the circle
295 if constexpr (isBox<bulkIdx>())
296 {
297 psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
298
299 const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
300 ShapeValues shapeValues;
301 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, globalPos, shapeValues);
302 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
303 }
304 else
305 {
306 psData.addCircleInterpolation(circleIpWeight, circleStencil);
307 psData.addBulkInterpolation(bulkElementIdx);
308 }
309
310 // publish point source data in the global vector
311 this->pointSourceData().emplace_back(std::move(psData));
312
313 // mean distance to outside element for source correction schemes
314 const auto outsideGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
315 this->averageDistanceToBulkCell().push_back(averageDistancePointGeometry(globalPos, outsideGeometry));
316
317 // export the bulk coupling stencil
318 if constexpr (isBox<lowDimIdx>())
319 {
320 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
321 this->vertexIndices(lowDimIdx, lowDimElementIdx).begin(),
322 this->vertexIndices(lowDimIdx, lowDimElementIdx).end());
323
324 }
325 else
326 {
327 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
328 }
329
330 // export bulk circle stencil
331 if constexpr (isBox<bulkIdx>())
332 {
333 // we insert all vertices and make it unique later
334 for (const auto& vertices : circleCornerIndices)
335 {
336 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
337 vertices->begin(), vertices->end());
338
339 }
340 }
341 else
342 {
343 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
344 circleStencil.begin(), circleStencil.end());
345 }
346 }
347 }
348 }
349
350 // make the circle stencil unique (for source derivatives)
351 for (auto&& stencil : extendedSourceStencil_.stencil())
352 {
353 std::sort(stencil.second.begin(), stencil.second.end());
354 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
355
356 // remove the vertices element (box)
357 if constexpr (isBox<bulkIdx>())
358 {
359 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
360 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
361 [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
362 stencil.second.end());
363 }
364 // remove the own element (cell-centered)
365 else
366 {
367 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
368 [&](auto i){ return i == stencil.first; }),
369 stencil.second.end());
370 }
371 }
372
373 // make stencils unique
374 using namespace Dune::Hybrid;
375 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
376 {
377 for (auto&& stencil : this->couplingStencils(domainIdx))
378 {
379 std::sort(stencil.second.begin(), stencil.second.end());
380 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
381 }
382 });
383
384 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
385 }
386
389 {
390 // resize the storage vector
391 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
392 // get references to the grid geometries
393 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
394 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
395
396 // compute the low dim volume fractions
397 for (const auto& is : intersections(this->glue()))
398 {
399 // all inside elements are identical...
400 const auto& inside = is.targetEntity(0);
401 const auto intersectionGeometry = is.geometry();
402 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
403
404 // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional
405 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
406 for (int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
407 {
408 const auto& outside = is.domainEntity(outsideIdx);
409 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
410 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
411 }
412 }
413 }
414
418 // \{
419
421 Scalar radius(std::size_t id) const
422 {
423 const auto& data = this->pointSourceData()[id];
424 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
425 }
426
428 // For one-dimensional low dim domain we assume radial tubes
429 Scalar lowDimVolume(const Element<bulkIdx>& element) const
430 {
431 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
432 return lowDimVolumeInBulkElement_[eIdx];
433 }
434
436 // For one-dimensional low dim domain we assume radial tubes
437 Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const
438 {
439 const auto totalVolume = element.geometry().volume();
440 return lowDimVolume(element) / totalVolume;
441 }
442
443 // \}
444
448 const typename ParentType::template CouplingStencils<bulkIdx>::mapped_type&
449 extendedSourceStencil(std::size_t eIdx) const
450 {
451 const auto& sourceStencils = extendedSourceStencil_.stencil();
452 if (auto stencil = sourceStencils.find(eIdx); stencil != sourceStencils.end())
453 return stencil->second;
454
455 return this->emptyStencil(bulkIdx);
456 }
457
458private:
461
463 std::vector<Scalar> lowDimVolumeInBulkElement_;
464};
465
467template<class MDTraits>
468struct CouplingManagerSupportsMultithreadedAssembly<Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Average>>
469: public std::true_type {};
470
471} // end namespace Dumux
472
473#endif
Point source traits for average-based coupling modes.
Helper function to compute points on a circle.
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d_average.hh:56
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_average.hh:107
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_average.hh:92
const ParentType::template CouplingStencils< bulkIdx >::mapped_type & extendedSourceStencil(std::size_t eIdx) const
Extended source stencil (for the bulk domain)
Definition: couplingmanager1d3d_average.hh:449
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d_average.hh:388
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_average.hh:138
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_average.hh:437
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_average.hh:120
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_average.hh:421
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_average.hh:429
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:24
A class managing an extended source stencil.
Definition: embedded/extendedsourcestencil.hh:36
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:71
Defines all properties used in Dumux.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
Helper functions for distance queries.
Extended source stencil helper class for coupling managers.
void circlePoints(std::vector< GlobalPosition > &points, const std::vector< Scalar > &sincos, const GlobalPosition &center, const GlobalPosition &normal, const Scalar radius)
Definition: circlepoints.hh:38
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
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:102
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:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Defines the index types used for grid and local indices.
constexpr Box box
Definition: method.hh:147
constexpr Average average
Definition: couplingmanager1d3d_average.hh:39
Definition: adapt.hh:17
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78
Definition: couplingmanager1d3d_average.hh:35
static std::string name()
Definition: couplingmanager1d3d_average.hh:36
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26
Helper class to create (named and comparable) tagged types Tags any given type. The tagged type is eq...
Definition: tag.hh:30
Helper class to create (named and comparable) tagged types.