3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_SURFACE_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_SURFACE_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
41
42namespace Dumux {
43
44namespace Embedded1d3dCouplingMode {
45struct Surface : public Utility::Tag<Surface> {
46 static std::string name() { return "surface"; }
47};
48
49inline constexpr Surface surface{};
50} // end namespace Embedded1d3dCouplingMode
51
52// forward declaration
53template<class MDTraits, class CouplingMode>
54class Embedded1d3dCouplingManager;
55
63template<class MDTraits>
64class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface>
65: public EmbeddedCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface>,
66 CircleAveragePointSourceTraits<MDTraits>>
67{
70 using Scalar = typename MDTraits::Scalar;
71 using SolutionVector = typename MDTraits::SolutionVector;
72 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
73
74 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
75 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
76
77 // the sub domain type aliases
78 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
79 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
80 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
81 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
82 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
83 template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex;
84
85 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
86
87 template<std::size_t id>
88 static constexpr bool isBox()
89 { return GridGeometry<id>::discMethod == DiscretizationMethods::box; }
90
91 enum {
92 bulkDim = GridView<bulkIdx>::dimension,
93 lowDimDim = GridView<lowDimIdx>::dimension,
94 dimWorld = GridView<bulkIdx>::dimensionworld
95 };
96public:
97 static constexpr Embedded1d3dCouplingMode::Surface couplingMode{};
98
99 using ParentType::ParentType;
100
101 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
102 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
103 const SolutionVector& curSol)
104 {
105 ParentType::init(bulkProblem, lowDimProblem, curSol);
106 computeLowDimVolumeFractions();
107 }
108
115 template<std::size_t id, class JacobianPattern>
116 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
117 {
118 extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern);
119 }
120
128 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
129 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
130 const LocalAssemblerI& localAssemblerI,
131 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
132 JacobianMatrixDiagBlock& A,
133 GridVariables& gridVariables)
134 {
135 extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, A, gridVariables);
136 }
137
138 /* \brief Compute integration point point sources and associated data
139 *
140 * This method uses grid glue to intersect the given grids. Over each intersection
141 * we later need to integrate a source term. This method places point sources
142 * at each quadrature point and provides the point source with the necessary
143 * information to compute integrals (quadrature weight and integration element)
144 * \param order The order of the quadrature rule for integration of sources over an intersection
145 * \param verbose If the point source computation is verbose
146 */
147 void computePointSourceData(std::size_t order = 1, bool verbose = false)
148 {
149 // if we use the circle average as the 3D values or a point evaluation
150 static const bool useCircleAverage = getParam<bool>("MixedDimension.UseCircleAverage", true);
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_.stencil().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 points on the cylinder surface at the integration point
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 integrationElement = intersectionGeometry.integrationElement(qp.position())*2*M_PI*radius/Scalar(numIp);
214 const auto qpweight = qp.weight()/(2*M_PI*radius);
215 const auto circleAvgWeight = 2*M_PI*radius/numIp;
216
217 const auto circlePoints = EmbeddedCoupling::circlePoints(globalPos, normal, radius, numIp);
218 std::vector<std::vector<std::size_t>> circleBulkElementIndices(circlePoints.size());
219 std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(circlePoints.size());
220 std::vector<GridIndex<bulkIdx>> circleStencil; circleStencil.reserve(circlePoints.size());
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 for (int k = 0; k < circlePoints.size(); ++k)
228 {
229 circleBulkElementIndices[k] = intersectingEntities(circlePoints[k], bulkTree);
230 if (circleBulkElementIndices[k].empty())
231 continue;
232
233 const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices[k].size();
234 for (const auto bulkElementIdx : circleBulkElementIndices[k])
235 {
236 circleStencil.push_back(bulkElementIdx);
237 circleIpWeight.push_back(localCircleAvgWeight);
238
239 // precompute interpolation data for box scheme for each cut bulk element
240 if constexpr (isBox<bulkIdx>())
241 {
242 const auto bulkElement = bulkGridGeometry.element(bulkElementIdx);
243 circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx)));
244
245 // evaluate shape functions at the integration point
246 const auto bulkGeometry = bulkElement.geometry();
247 ShapeValues shapeValues;
248 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePoints[k], shapeValues);
249 circleShapeValues.emplace_back(std::move(shapeValues));
250 }
251 }
252 }
253
254 // export low dim circle stencil
255 if constexpr (isBox<bulkIdx>())
256 {
257 // we insert all vertices and make it unique later
258 for (const auto& vertices : circleCornerIndices)
259 {
260 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
261 vertices->begin(), vertices->end());
262
263 }
264 }
265 else
266 {
267 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
268 circleStencil.begin(), circleStencil.end());
269 }
270
271 for (int k = 0; k < circlePoints.size(); ++k)
272 {
273 const auto& circlePos = circlePoints[k];
274
275 // find bulk elements intersection with the circle elements
276 if (circleBulkElementIndices[k].empty())
277 continue;
278
279 // loop over the bulk elements at the integration points (usually one except when it is on a face or edge or vertex)
280 // and add a point source at every point on the circle
281 for (const auto bulkElementIdx : circleBulkElementIndices[k])
282 {
283 const auto id = this->idCounter_++;
284
285 this->pointSources(bulkIdx).emplace_back(circlePos, id, qpweight, integrationElement, bulkElementIdx);
286 this->pointSources(bulkIdx).back().setEmbeddings(circleBulkElementIndices[k].size());
287 this->pointSources(lowDimIdx).emplace_back(globalPos, id, qpweight, integrationElement, lowDimElementIdx);
288 this->pointSources(lowDimIdx).back().setEmbeddings(circleBulkElementIndices[k].size());
289
290 // pre compute additional data used for the evaluation of
291 // the actual solution dependent source term
292 PointSourceData psData;
293
294 if constexpr (isBox<lowDimIdx>())
295 {
296 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
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 if (useCircleAverage)
310 psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
311
312 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
313 const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
314 ShapeValues shapeValues;
315 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePos, shapeValues);
316 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
317 }
318 else
319 {
320 if (useCircleAverage)
321 psData.addCircleInterpolation(circleIpWeight, circleStencil);
322
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 // export the bulk coupling stencil
330 // we insert all vertices / elements and make it unique later
331 if constexpr (isBox<lowDimIdx>())
332 {
333 const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx);
334 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
335 vertices.begin(), vertices.end());
336
337 }
338 else
339 {
340 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
341 }
342
343 // export bulk circle stencil (only needed for circle average)
344 if (useCircleAverage)
345 {
346 if constexpr (isBox<bulkIdx>())
347 {
348 // we insert all vertices and make it unique later
349 for (const auto& vertices : circleCornerIndices)
350 {
351 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
352 vertices->begin(), vertices->end());
353
354 }
355 }
356 else
357 {
358 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
359 circleStencil.begin(), circleStencil.end());
360 }
361 }
362 }
363 }
364 }
365 }
366
367 // make the circle stencil unique (for source derivatives)
368 for (auto&& stencil : extendedSourceStencil_.stencil())
369 {
370 std::sort(stencil.second.begin(), stencil.second.end());
371 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
372
373 // remove the vertices element (box)
374 if constexpr (isBox<bulkIdx>())
375 {
376 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
377 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
378 [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
379 stencil.second.end());
380 }
381 // remove the own element (cell-centered)
382 else
383 {
384 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
385 [&](auto i){ return i == stencil.first; }),
386 stencil.second.end());
387 }
388 }
389
390 // make stencils unique
391 using namespace Dune::Hybrid;
392 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
393 {
394 for (auto&& stencil : this->couplingStencils(domainIdx))
395 {
396 std::sort(stencil.second.begin(), stencil.second.end());
397 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
398 }
399 });
400
401 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
402 }
403
406 {
407 // resize the storage vector
408 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
409 // get references to the grid geometries
410 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
411 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
412
413 // compute the low dim volume fractions
414 for (const auto& is : intersections(this->glue()))
415 {
416 // all inside elements are identical...
417 const auto& inside = is.targetEntity(0);
418 const auto intersectionGeometry = is.geometry();
419 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
420
421 // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional
422 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
423 for (int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
424 {
425 const auto& outside = is.domainEntity(outsideIdx);
426 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
427 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
428 }
429 }
430 }
431
435 // \{
436
438 Scalar radius(std::size_t id) const
439 {
440 const auto& data = this->pointSourceData()[id];
441 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
442 }
443
445 // For one-dimensional low dim domain we assume radial tubes
446 Scalar lowDimVolume(const Element<bulkIdx>& element) const
447 {
448 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
449 return lowDimVolumeInBulkElement_[eIdx];
450 }
451
453 // For one-dimensional low dim domain we assume radial tubes
454 Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const
455 {
456 const auto totalVolume = element.geometry().volume();
457 return lowDimVolume(element) / totalVolume;
458 }
459
460 // \}
461
462private:
465
467 std::vector<Scalar> lowDimVolumeInBulkElement_;
468};
469
470} // end namespace Dumux
471
472#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.
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 Box box
Definition: method.hh:139
constexpr Surface surface
Definition: couplingmanager1d3d_surface.hh:49
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:102
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_surface.hh:45
static std::string name()
Definition: couplingmanager1d3d_surface.hh:46
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d_surface.hh:67
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_surface.hh:147
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:129
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_surface.hh:446
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:116
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d_surface.hh:405
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_surface.hh:438
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_surface.hh:454
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_surface.hh:101
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.