version 3.11-dev
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-FileCopyrightText: 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 Embedded1d3dCouplingManagerBase<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
98 template<std::size_t id, class JacobianPattern>
99 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
100 {
101 extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern);
102 }
103
111 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
112 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
113 const LocalAssemblerI& localAssemblerI,
114 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
115 JacobianMatrixDiagBlock& A,
116 GridVariables& gridVariables)
117 {
118 extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, A, gridVariables);
119 }
120
121 /* \brief Compute integration point point sources and associated data
122 *
123 * This method uses grid glue to intersect the given grids. Over each intersection
124 * we later need to integrate a source term. This method places point sources
125 * at each quadrature point and provides the point source with the necessary
126 * information to compute integrals (quadrature weight and integration element)
127 * \param order The order of the quadrature rule for integration of sources over an intersection
128 * \param verbose If the point source computation is verbose
129 */
130 void computePointSourceData(std::size_t order = 1, bool verbose = false)
131 {
132 // Initialize the bulk bounding box tree
133 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
134 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
135 const auto& bulkTree = bulkGridGeometry.boundingBoxTree();
136
137 // initialize the maps
138 // do some logging and profiling
139 Dune::Timer watch;
140 std::cout << "Initializing the point sources..." << std::endl;
141
142 // clear all internal members like pointsource vectors and stencils
143 // initializes the point source id counter
144 this->clear();
145 extendedSourceStencil_.clear();
146
147 // precompute the vertex indices for efficiency
148 this->precomputeVertexIndices(bulkIdx);
149 this->precomputeVertexIndices(lowDimIdx);
150
151 // intersect the bounding box trees
152 this->glueGrids();
153
154 // iterate over all intersection and add point sources
155 const auto& lowDimProblem = this->problem(lowDimIdx);
156 for (const auto& is : intersections(this->glue()))
157 {
158 // all inside elements are identical...
159 const auto& lowDimElement = is.targetEntity(0);
160 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(lowDimElement);
161
162 // get the intersection geometry
163 const auto intersectionGeometry = is.geometry();
164 // get the Gaussian quadrature rule for the local intersection
165 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order);
166
167 // apply the Gaussian quadrature rule and define point sources at each quadrature point
168 // note that the approximation is not optimal if
169 // (a) the one-dimensional elements are too large,
170 // (b) whenever a one-dimensional element is split between two or more elements,
171 // (c) when gradients of important quantities in the three-dimensional domain are large.
172
173 // iterate over all quadrature points
174 for (auto&& qp : quad)
175 {
176 // global position of the quadrature point
177 const auto globalPos = intersectionGeometry.global(qp.position());
178
179 const auto bulkElementIndices = intersectingEntities(globalPos, bulkTree);
180
181 // do not add a point source if the qp is outside of the 3d grid
182 // this is equivalent to having a source of zero for that qp
183 if (bulkElementIndices.empty())
184 continue;
185
187 // get circle average connectivity and interpolation data
189
190 static const auto numIp = getParam<int>("MixedDimension.NumCircleSegments");
191 const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx);
192 const auto normal = intersectionGeometry.corner(1)-intersectionGeometry.corner(0);
193 const auto circleAvgWeight = 2*M_PI*radius/numIp;
194 const auto integrationElement = intersectionGeometry.integrationElement(qp.position());
195 const auto qpweight = qp.weight();
196
197 const auto circlePoints = EmbeddedCoupling::circlePoints(globalPos, normal, radius, numIp);
198 std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(circlePoints.size());
199 std::vector<GridIndex<bulkIdx>> circleStencil; circleStencil.reserve(circlePoints.size());
200
201 // for box
202 std::vector<const std::vector<GridIndex<bulkIdx>>*> circleCornerIndices;
203 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
204 std::vector<ShapeValues> circleShapeValues;
205
206 // go over all points of the average operator
207 int insideCirclePoints = 0;
208 for (int k = 0; k < circlePoints.size(); ++k)
209 {
210 const auto circleBulkElementIndices = intersectingEntities(circlePoints[k], bulkTree);
211 if (circleBulkElementIndices.empty())
212 continue;
213
214 ++insideCirclePoints;
215 const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices.size();
216 for (const auto bulkElementIdx : circleBulkElementIndices)
217 {
218 circleStencil.push_back(bulkElementIdx);
219 circleIpWeight.push_back(localCircleAvgWeight);
220
221 // precompute interpolation data for box scheme for each cut bulk element
222 if constexpr (isBox<bulkIdx>())
223 {
224 const auto bulkElement = bulkGridGeometry.element(bulkElementIdx);
225 circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx)));
226
227 // evaluate shape functions at the integration point
228 const auto bulkGeometry = bulkElement.geometry();
229 ShapeValues shapeValues;
230 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePoints[k], shapeValues);
231 circleShapeValues.emplace_back(std::move(shapeValues));
232 }
233 }
234 }
235
236 // if the circle stencil is empty (that is the circle is entirely outside the domain)
237 // we do not add a (zero) point source
238 if (circleStencil.empty())
239 continue;
240
241 // export low dim circle stencil
242 if constexpr (isBox<bulkIdx>())
243 {
244 // we insert all vertices and make it unique later
245 for (const auto& vertices : circleCornerIndices)
246 {
247 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
248 vertices->begin(), vertices->end());
249
250 }
251 }
252 else
253 {
254 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
255 circleStencil.begin(), circleStencil.end());
256 }
257
258 // surface fraction that is inside the domain (only different than 1.0 on the boundary)
259 const auto surfaceFraction = Scalar(insideCirclePoints)/Scalar(circlePoints.size());
260
261 // loop over the bulk elements at the integration points (usually one except when it is on a face or edge or vertex)
262 for (auto bulkElementIdx : bulkElementIndices)
263 {
264 const auto id = this->idCounter_++;
265
266 this->pointSources(bulkIdx).emplace_back(globalPos, id, qpweight, integrationElement*surfaceFraction, bulkElementIdx);
267 this->pointSources(bulkIdx).back().setEmbeddings(bulkElementIndices.size());
268 this->pointSources(lowDimIdx).emplace_back(globalPos, id, qpweight, integrationElement*surfaceFraction, lowDimElementIdx);
269 this->pointSources(lowDimIdx).back().setEmbeddings(bulkElementIndices.size());
270
271 // pre compute additional data used for the evaluation of
272 // the actual solution dependent source term
273 PointSourceData psData;
274
275 if constexpr (isBox<lowDimIdx>())
276 {
277 ShapeValues shapeValues;
278 this->getShapeValues(lowDimIdx, lowDimGridGeometry, lowDimElement.geometry(), globalPos, shapeValues);
279 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
280 }
281 else
282 {
283 psData.addLowDimInterpolation(lowDimElementIdx);
284 }
285
286 // add data needed to compute integral over the circle
287 if constexpr (isBox<bulkIdx>())
288 {
289 psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
290
291 const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
292 ShapeValues shapeValues;
293 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, globalPos, shapeValues);
294 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
295 }
296 else
297 {
298 psData.addCircleInterpolation(circleIpWeight, circleStencil);
299 psData.addBulkInterpolation(bulkElementIdx);
300 }
301
302 // publish point source data in the global vector
303 this->pointSourceData().emplace_back(std::move(psData));
304
305 // mean distance to outside element for source correction schemes
306 const auto outsideGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
307 this->averageDistanceToBulkCell().push_back(averageDistancePointGeometry(globalPos, outsideGeometry));
308
309 // export the bulk coupling stencil
310 if constexpr (isBox<lowDimIdx>())
311 {
312 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
313 this->vertexIndices(lowDimIdx, lowDimElementIdx).begin(),
314 this->vertexIndices(lowDimIdx, lowDimElementIdx).end());
315
316 }
317 else
318 {
319 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
320 }
321
322 // export bulk circle stencil
323 if constexpr (isBox<bulkIdx>())
324 {
325 // we insert all vertices and make it unique later
326 for (const auto& vertices : circleCornerIndices)
327 {
328 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
329 vertices->begin(), vertices->end());
330
331 }
332 }
333 else
334 {
335 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
336 circleStencil.begin(), circleStencil.end());
337 }
338 }
339 }
340 }
341
342 // make the circle stencil unique (for source derivatives)
343 for (auto&& stencil : extendedSourceStencil_.stencil())
344 {
345 std::sort(stencil.second.begin(), stencil.second.end());
346 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
347
348 // remove the vertices element (box)
349 if constexpr (isBox<bulkIdx>())
350 {
351 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
352 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
353 [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
354 stencil.second.end());
355 }
356 // remove the own element (cell-centered)
357 else
358 {
359 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
360 [&](auto i){ return i == stencil.first; }),
361 stencil.second.end());
362 }
363 }
364
365 // make stencils unique
366 using namespace Dune::Hybrid;
367 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
368 {
369 for (auto&& stencil : this->couplingStencils(domainIdx))
370 {
371 std::sort(stencil.second.begin(), stencil.second.end());
372 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
373 }
374 });
375
376 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
377 }
378
382 const typename ParentType::template CouplingStencils<bulkIdx>::mapped_type&
383 extendedSourceStencil(std::size_t eIdx) const
384 {
385 const auto& sourceStencils = extendedSourceStencil_.stencil();
386 if (auto stencil = sourceStencils.find(eIdx); stencil != sourceStencils.end())
387 return stencil->second;
388
389 return this->emptyStencil(bulkIdx);
390 }
391
392private:
395
396};
397
399template<class MDTraits>
400struct CouplingManagerSupportsMultithreadedAssembly<Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Average>>
401: public std::true_type {};
402
403} // end namespace Dumux
404
405#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:99
const ParentType::template CouplingStencils< bulkIdx >::mapped_type & extendedSourceStencil(std::size_t eIdx) const
Extended source stencil (for the bulk domain)
Definition: couplingmanager1d3d_average.hh:383
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_average.hh:130
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:112
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3dbase.hh:35
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
Defines all properties used in Dumux.
Coupling manager for low-dimensional domains embedded in the bulk domain.
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.