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