version 3.11-dev
couplingmanager1d3d_surface.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_SURFACE_HH
14#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_SURFACE_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
29
30namespace Dumux {
31
32namespace Embedded1d3dCouplingMode {
33struct Surface : public Utility::Tag<Surface> {
34 static std::string name() { return "surface"; }
35};
36
37inline constexpr Surface surface{};
38} // end namespace Embedded1d3dCouplingMode
39
40// forward declaration
41template<class MDTraits, class CouplingMode>
42class Embedded1d3dCouplingManager;
43
51template<class MDTraits>
52class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface>
53: public Embedded1d3dCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface>,
54 CircleAveragePointSourceTraits<MDTraits>>
55{
58 using Scalar = typename MDTraits::Scalar;
59 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
60
61 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
62 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
63
64 // the sub domain type aliases
65 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
66 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
67 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
68 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
69 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
70 template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex;
71
72 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
73
74 template<std::size_t id>
75 static constexpr bool isBox()
76 { return GridGeometry<id>::discMethod == DiscretizationMethods::box; }
77
78 enum {
79 bulkDim = GridView<bulkIdx>::dimension,
80 lowDimDim = GridView<lowDimIdx>::dimension,
81 dimWorld = GridView<bulkIdx>::dimensionworld
82 };
83public:
84 static constexpr Embedded1d3dCouplingMode::Surface couplingMode{};
85
86 using ParentType::ParentType;
87
94 template<std::size_t id, class JacobianPattern>
95 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
96 {
97 extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern);
98 }
99
107 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
108 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
109 const LocalAssemblerI& localAssemblerI,
110 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
111 JacobianMatrixDiagBlock& A,
112 GridVariables& gridVariables)
113 {
114 extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, A, gridVariables);
115 }
116
117 /* \brief Compute integration point point sources and associated data
118 *
119 * This method uses grid glue to intersect the given grids. Over each intersection
120 * we later need to integrate a source term. This method places point sources
121 * at each quadrature point and provides the point source with the necessary
122 * information to compute integrals (quadrature weight and integration element)
123 * \param order The order of the quadrature rule for integration of sources over an intersection
124 * \param verbose If the point source computation is verbose
125 */
126 void computePointSourceData(std::size_t order = 1, bool verbose = false)
127 {
128 // if we use the circle average as the 3D values or a point evaluation
129 static const bool useCircleAverage = getParam<bool>("MixedDimension.UseCircleAverage", true);
130
131 // Initialize the bulk bounding box tree
132 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
133 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
134 const auto& bulkTree = bulkGridGeometry.boundingBoxTree();
135
136 // initialize the maps
137 // do some logging and profiling
138 Dune::Timer watch;
139 std::cout << "Initializing the point sources..." << std::endl;
140
141 // clear all internal members like pointsource vectors and stencils
142 // initializes the point source id counter
143 this->clear();
144 extendedSourceStencil_.stencil().clear();
145
146 // precompute the vertex indices for efficiency
147 this->precomputeVertexIndices(bulkIdx);
148 this->precomputeVertexIndices(lowDimIdx);
149
150 // intersect the bounding box trees
151 this->glueGrids();
152
153 // iterate over all intersection and add point sources
154 const auto& lowDimProblem = this->problem(lowDimIdx);
155 for (const auto& is : intersections(this->glue()))
156 {
157 // all inside elements are identical...
158 const auto& lowDimElement = is.targetEntity(0);
159 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(lowDimElement);
160
161 // get the intersection geometry
162 const auto intersectionGeometry = is.geometry();
163 // get the Gaussian quadrature rule for the local intersection
164 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order);
165
166 // apply the Gaussian quadrature rule and define point sources at each quadrature point
167 // note that the approximation is not optimal if
168 // (a) the one-dimensional elements are too large,
169 // (b) whenever a one-dimensional element is split between two or more elements,
170 // (c) when gradients of important quantities in the three-dimensional domain are large.
171
172 // iterate over all quadrature points
173 for (auto&& qp : quad)
174 {
175 // global position of the quadrature point
176 const auto globalPos = intersectionGeometry.global(qp.position());
177
178 const auto bulkElementIndices = intersectingEntities(globalPos, bulkTree);
179
180 // do not add a point source if the qp is outside of the 3d grid
181 // this is equivalent to having a source of zero for that qp
182 if (bulkElementIndices.empty())
183 continue;
184
186 // get points on the cylinder surface at the integration point
188
189 static const auto numIp = getParam<int>("MixedDimension.NumCircleSegments");
190 const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx);
191 const auto normal = intersectionGeometry.corner(1)-intersectionGeometry.corner(0);
192 const auto integrationElement = intersectionGeometry.integrationElement(qp.position())*2*M_PI*radius/Scalar(numIp);
193 const auto qpweight = qp.weight()/(2*M_PI*radius);
194 const auto circleAvgWeight = 2*M_PI*radius/numIp;
195
196 const auto circlePoints = EmbeddedCoupling::circlePoints(globalPos, normal, radius, numIp);
197 std::vector<std::vector<std::size_t>> circleBulkElementIndices(circlePoints.size());
198 std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(circlePoints.size());
199 std::vector<GridIndex<bulkIdx>> circleStencil; circleStencil.reserve(circlePoints.size());
200 // for box
201 std::vector<const std::vector<GridIndex<bulkIdx>>*> circleCornerIndices;
202 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
203 std::vector<ShapeValues> circleShapeValues;
204
205 // go over all points of the average operator
206 for (int k = 0; k < circlePoints.size(); ++k)
207 {
208 circleBulkElementIndices[k] = intersectingEntities(circlePoints[k], bulkTree);
209 if (circleBulkElementIndices[k].empty())
210 continue;
211
212 const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices[k].size();
213 for (const auto bulkElementIdx : circleBulkElementIndices[k])
214 {
215 circleStencil.push_back(bulkElementIdx);
216 circleIpWeight.push_back(localCircleAvgWeight);
217
218 // precompute interpolation data for box scheme for each cut bulk element
219 if constexpr (isBox<bulkIdx>())
220 {
221 const auto bulkElement = bulkGridGeometry.element(bulkElementIdx);
222 circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx)));
223
224 // evaluate shape functions at the integration point
225 const auto bulkGeometry = bulkElement.geometry();
226 ShapeValues shapeValues;
227 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePoints[k], shapeValues);
228 circleShapeValues.emplace_back(std::move(shapeValues));
229 }
230 }
231 }
232
233 // export low dim circle stencil
234 if constexpr (isBox<bulkIdx>())
235 {
236 // we insert all vertices and make it unique later
237 for (const auto& vertices : circleCornerIndices)
238 {
239 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
240 vertices->begin(), vertices->end());
241
242 }
243 }
244 else
245 {
246 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
247 circleStencil.begin(), circleStencil.end());
248 }
249
250 for (int k = 0; k < circlePoints.size(); ++k)
251 {
252 const auto& circlePos = circlePoints[k];
253
254 // find bulk elements intersection with the circle elements
255 if (circleBulkElementIndices[k].empty())
256 continue;
257
258 // loop over the bulk elements at the integration points (usually one except when it is on a face or edge or vertex)
259 // and add a point source at every point on the circle
260 for (const auto bulkElementIdx : circleBulkElementIndices[k])
261 {
262 const auto id = this->idCounter_++;
263
264 this->pointSources(bulkIdx).emplace_back(circlePos, id, qpweight, integrationElement, bulkElementIdx);
265 this->pointSources(bulkIdx).back().setEmbeddings(circleBulkElementIndices[k].size());
266 this->pointSources(lowDimIdx).emplace_back(globalPos, id, qpweight, integrationElement, lowDimElementIdx);
267 this->pointSources(lowDimIdx).back().setEmbeddings(circleBulkElementIndices[k].size());
268
269 // pre compute additional data used for the evaluation of
270 // the actual solution dependent source term
271 PointSourceData psData;
272
273 if constexpr (isBox<lowDimIdx>())
274 {
275 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
276 ShapeValues shapeValues;
277 this->getShapeValues(lowDimIdx, lowDimGridGeometry, lowDimElement.geometry(), globalPos, shapeValues);
278 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
279 }
280 else
281 {
282 psData.addLowDimInterpolation(lowDimElementIdx);
283 }
284
285 // add data needed to compute integral over the circle
286 if constexpr (isBox<bulkIdx>())
287 {
288 if (useCircleAverage)
289 psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
290
291 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
292 const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
293 ShapeValues shapeValues;
294 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePos, shapeValues);
295 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
296 }
297 else
298 {
299 if (useCircleAverage)
300 psData.addCircleInterpolation(circleIpWeight, circleStencil);
301
302 psData.addBulkInterpolation(bulkElementIdx);
303 }
304
305 // publish point source data in the global vector
306 this->pointSourceData().emplace_back(std::move(psData));
307
308 // export the bulk coupling stencil
309 // we insert all vertices / elements and make it unique later
310 if constexpr (isBox<lowDimIdx>())
311 {
312 const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx);
313 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
314 vertices.begin(), vertices.end());
315
316 }
317 else
318 {
319 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
320 }
321
322 // export bulk circle stencil (only needed for circle average)
323 if (useCircleAverage)
324 {
325 if constexpr (isBox<bulkIdx>())
326 {
327 // we insert all vertices and make it unique later
328 for (const auto& vertices : circleCornerIndices)
329 {
330 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
331 vertices->begin(), vertices->end());
332
333 }
334 }
335 else
336 {
337 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
338 circleStencil.begin(), circleStencil.end());
339 }
340 }
341 }
342 }
343 }
344 }
345
346 // make the circle stencil unique (for source derivatives)
347 for (auto&& stencil : extendedSourceStencil_.stencil())
348 {
349 std::sort(stencil.second.begin(), stencil.second.end());
350 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
351
352 // remove the vertices element (box)
353 if constexpr (isBox<bulkIdx>())
354 {
355 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
356 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
357 [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
358 stencil.second.end());
359 }
360 // remove the own element (cell-centered)
361 else
362 {
363 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
364 [&](auto i){ return i == stencil.first; }),
365 stencil.second.end());
366 }
367 }
368
369 // make stencils unique
370 using namespace Dune::Hybrid;
371 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
372 {
373 for (auto&& stencil : this->couplingStencils(domainIdx))
374 {
375 std::sort(stencil.second.begin(), stencil.second.end());
376 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
377 }
378 });
379
380 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
381 }
382
383
387 const typename ParentType::template CouplingStencils<bulkIdx>::mapped_type&
388 extendedSourceStencil(std::size_t eIdx) const
389 {
390 const auto& sourceStencils = extendedSourceStencil_.stencil();
391 if (auto stencil = sourceStencils.find(eIdx); stencil != sourceStencils.end())
392 return stencil->second;
393
394 return this->emptyStencil(bulkIdx);
395 }
396
397private:
400
402 std::vector<Scalar> lowDimVolumeInBulkElement_;
403};
404
406template<class MDTraits>
407struct CouplingManagerSupportsMultithreadedAssembly<Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface>>
408: public std::true_type {};
409
410} // end namespace Dumux
411
412#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_surface.hh:55
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_surface.hh:126
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_surface.hh:108
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_surface.hh:95
const ParentType::template CouplingStencils< bulkIdx >::mapped_type & extendedSourceStencil(std::size_t eIdx) const
Extended source stencil (for the bulk domain)
Definition: couplingmanager1d3d_surface.hh:388
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.
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
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 Surface surface
Definition: couplingmanager1d3d_surface.hh:37
Definition: adapt.hh:17
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78
Definition: couplingmanager1d3d_surface.hh:33
static std::string name()
Definition: couplingmanager1d3d_surface.hh:34
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.