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