ecfvstencil.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_ECFV_STENCIL_HH
29#define EWOMS_ECFV_STENCIL_HH
30
31#include <dune/common/fvector.hh>
32#include <dune/common/version.hh>
33
34#include <dune/geometry/type.hh>
35
36#include <dune/grid/common/mcmgmapper.hh>
37#include <dune/grid/common/intersectioniterator.hh>
38
39#include <opm/common/ErrorMacros.hpp>
40
41#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
42
43#include <opm/material/common/ConditionalStorage.hpp>
44
46
47#include <cassert>
48#include <cstddef>
49#include <stdexcept>
50#include <vector>
51
52namespace Opm {
53
64template <class Scalar,
65 class GridView,
66 bool needFaceIntegrationPos = true,
67 bool needFaceNormal = true>
69{
70 enum { dimWorld = GridView::dimensionworld };
71
72 using CoordScalar = typename GridView::ctype;
73 using Intersection = typename GridView::Intersection;
74 using Element = typename GridView::template Codim<0>::Entity;
75
76 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
77
78 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
79
80 using WorldVector = Dune::FieldVector<Scalar, dimWorld>;
81
82public:
83 using Entity = Element;
84 using Mapper = ElementMapper;
85
86 using LocalGeometry = typename Element::Geometry;
87
97 {
98 public:
99 // default construct an uninitialized object.
100 // this is only here because std::vector needs it...
102 {}
103
104 explicit SubControlVolume(const Element& element)
105 : element_(element)
106 { update(); }
107
108 void update(const Element& element)
109 { element_ = element; }
110
111 void update()
112 {}
113
117 decltype(auto) globalPos() const
118 { return element_.geometry().center(); }
119
123 decltype(auto) center() const
124 { return element_.geometry().center(); }
125
129 Scalar volume() const
130 { return element_.geometry().volume(); }
131
136 { return element_.geometry(); }
137
142 { return element_.geometryInFather(); }
143
144 private:
145 Element element_;
146 };
147
151 template <bool needIntegrationPos, bool needNormal>
153 {
154 public:
156 {}
157
158 EcfvSubControlVolumeFace(const Intersection& intersection, unsigned localNeighborIdx)
159 {
160 exteriorIdx_ = static_cast<unsigned short>(localNeighborIdx);
161
162 if constexpr (needNormal) {
163 (*normal_) = intersection.centerUnitOuterNormal();
164 }
165 const auto& geometry = intersection.geometry();
166 if constexpr (needIntegrationPos) {
167 (*integrationPos_) = geometry.center();
168 }
169 area_ = geometry.volume();
170 // TODO: facedir_ = intersection.boundaryId(); // This did not compile for some reason
171 // using indexInInside() as in ecltransmissibility.cc instead
172 // Note that indexInInside() returns -1 for NNC faces, that's the reason we
173 // cannot convert directly to a FaceDir::DirEnum (see FaceDir.hpp)
174 dirId_ = intersection.indexInInside(); // Legal values: -1, 0, 1, 2, 3, 4, 5
175 }
176
181 unsigned short interiorIndex() const
182 {
183 // The local index of the control volume in the interior
184 // of a face of the stencil in the element centered finite
185 // volume discretization is always the "central"
186 // element. In this implementation this element always has
187 // index 0....
188 return 0;
189 }
190
195 unsigned short exteriorIndex() const
196 { return exteriorIdx_; }
197
202 const GlobalPosition& integrationPos() const
203 { return *integrationPos_; }
204
209 const WorldVector& normal() const
210 { return *normal_; }
211
215 Scalar area() const
216 { return area_; }
217
224 int dirId() const
225 { return dirId_; }
226
230 FaceDir::DirEnum faceDirFromDirId() const
231 {
232 using Dir = FaceDir::DirEnum;
233 if (dirId_ == -1) {
234 OPM_THROW(std::runtime_error, "NNC faces does not have a face id");
235 }
236 switch(dirId_) {
237 case 0:
238 case 1:
239 return Dir::XPlus;
240 case 2:
241 case 3:
242 return Dir::YPlus;
243 case 4:
244 case 5:
245 return Dir::ZPlus;
246 default:
247 OPM_THROW(std::runtime_error,
248 "Unexpected face id" + std::to_string(dirId_));
249 }
250 }
251
252 private:
253 ConditionalStorage<needIntegrationPos, GlobalPosition> integrationPos_;
254 ConditionalStorage<needNormal, WorldVector> normal_;
255 Scalar area_;
256 int dirId_;
257
258 unsigned short exteriorIdx_;
259 };
260
262 using BoundaryFace = EcfvSubControlVolumeFace</*needFaceIntegrationPos=*/true, needFaceNormal>;
263
264 EcfvStencil(const GridView& gridView, const Mapper& mapper)
266 , elementMapper_(mapper)
267 {
268 // try to ensure that the mapper passed indeed maps elements
269 assert(int(gridView.size(/*codim=*/0)) == int(elementMapper_.size()));
270 }
271
272 void updateTopology(const Element& element)
273 {
274 // add the "center" element of the stencil
275 subControlVolumes_.clear();
276 subControlVolumes_.emplace_back(/*SubControlVolume(*/element/*)*/);
277 elements_.clear();
278 elements_.emplace_back(element);
279
280 interiorFaces_.clear();
281 boundaryFaces_.clear();
282
283 for (const auto& intersection : intersections(gridView_, element)) {
284 // if the current intersection has a neighbor, add a
285 // degree of freedom and an internal face, else add a
286 // boundary face
287 if (intersection.neighbor()) {
288 elements_.emplace_back(intersection.outside());
289 subControlVolumes_.emplace_back(/*SubControlVolume(*/elements_.back()/*)*/);
290 interiorFaces_.emplace_back(/*SubControlVolumeFace(*/intersection, subControlVolumes_.size() - 1/*)*/);
291 }
292 else {
293 boundaryFaces_.emplace_back(/*SubControlVolumeFace(*/intersection, - 10000/*)*/);
294 }
295 }
296 }
297
298 void updatePrimaryTopology(const Element& element)
299 {
300 // add the "center" element of the stencil
301 subControlVolumes_.clear();
302 subControlVolumes_.emplace_back(/*SubControlVolume(*/element/*)*/);
303 elements_.clear();
304 elements_.emplace_back(element);
305 }
306
307 void update(const Element& element)
308 {
310 }
311
313 {
314 assert(false); // not yet implemented
315 }
316
320 const Element& element() const
321 { return element(0); }
322
327 const GridView& gridView() const
328 { return *gridView_; }
329
334 std::size_t numDof() const
335 { return subControlVolumes_.size(); }
336
346 std::size_t numPrimaryDof() const
347 { return 1; }
348
353 unsigned globalSpaceIndex(unsigned dofIdx) const
354 {
355 assert(dofIdx < numDof());
356
357 return static_cast<unsigned>(elementMapper_.index(element(dofIdx)));
358 }
359
363 Dune::PartitionType partitionType(unsigned dofIdx) const
364 { return elements_[dofIdx]->partitionType(); }
365
373 const Element& element(unsigned dofIdx) const
374 {
375 assert(dofIdx < numDof());
376
377 return elements_[dofIdx];
378 }
379
384 const Entity& entity(unsigned dofIdx) const
385 { return element( dofIdx ); }
386
391 const SubControlVolume& subControlVolume(unsigned dofIdx) const
392 { return subControlVolumes_[dofIdx]; }
393
397 std::size_t numInteriorFaces() const
398 { return interiorFaces_.size(); }
399
404 const SubControlVolumeFace& interiorFace(unsigned faceIdx) const
405 { return interiorFaces_[faceIdx]; }
406
410 std::size_t numBoundaryFaces() const
411 { return boundaryFaces_.size(); }
412
417 const BoundaryFace& boundaryFace(unsigned bfIdx) const
418 { return boundaryFaces_[bfIdx]; }
419
420protected:
421 const GridView& gridView_;
422 const ElementMapper& elementMapper_;
423
424 std::vector<Element> elements_;
425 std::vector<SubControlVolume> subControlVolumes_;
426 std::vector<SubControlVolumeFace> interiorFaces_;
427 std::vector<BoundaryFace> boundaryFaces_;
428};
429
430} // namespace Opm
431
432#endif
Represents a face of a sub-control volume.
Definition: ecfvstencil.hh:153
const GlobalPosition & integrationPos() const
Returns the global position of the face's integration point.
Definition: ecfvstencil.hh:202
EcfvSubControlVolumeFace()
Definition: ecfvstencil.hh:155
unsigned short interiorIndex() const
Returns the local index of the degree of freedom to the face's interior.
Definition: ecfvstencil.hh:181
unsigned short exteriorIndex() const
Returns the local index of the degree of freedom to the face's outside.
Definition: ecfvstencil.hh:195
const WorldVector & normal() const
Returns the outer unit normal at the face's integration point.
Definition: ecfvstencil.hh:209
FaceDir::DirEnum faceDirFromDirId() const
Returns the direction of the face.
Definition: ecfvstencil.hh:230
Scalar area() const
Returns the area [m^2] of the face.
Definition: ecfvstencil.hh:215
int dirId() const
Returns the direction id of the face w.r.t the cell.
Definition: ecfvstencil.hh:224
EcfvSubControlVolumeFace(const Intersection &intersection, unsigned localNeighborIdx)
Definition: ecfvstencil.hh:158
Represents a sub-control volume.
Definition: ecfvstencil.hh:97
const LocalGeometry geometry() const
The geometry of the sub-control volume.
Definition: ecfvstencil.hh:135
Scalar volume() const
The volume [m^3] occupied by the sub-control volume.
Definition: ecfvstencil.hh:129
void update(const Element &element)
Definition: ecfvstencil.hh:108
SubControlVolume(const Element &element)
Definition: ecfvstencil.hh:104
const LocalGeometry localGeometry() const
Geometry of the sub-control volume relative to parent.
Definition: ecfvstencil.hh:141
SubControlVolume()
Definition: ecfvstencil.hh:101
decltype(auto) center() const
The center of the sub-control volume.
Definition: ecfvstencil.hh:123
decltype(auto) globalPos() const
The global position associated with the sub-control volume.
Definition: ecfvstencil.hh:117
void update()
Definition: ecfvstencil.hh:111
Represents the stencil (finite volume geometry) of a single element in the ECFV discretization.
Definition: ecfvstencil.hh:69
const SubControlVolumeFace & interiorFace(unsigned faceIdx) const
Returns the face object belonging to a given face index in the interior of the domain.
Definition: ecfvstencil.hh:404
std::vector< BoundaryFace > boundaryFaces_
Definition: ecfvstencil.hh:427
const BoundaryFace & boundaryFace(unsigned bfIdx) const
Returns the boundary face object belonging to a given boundary face index.
Definition: ecfvstencil.hh:417
const Element & element(unsigned dofIdx) const
Return the element given the index of a degree of freedom.
Definition: ecfvstencil.hh:373
std::vector< SubControlVolumeFace > interiorFaces_
Definition: ecfvstencil.hh:426
EcfvStencil(const GridView &gridView, const Mapper &mapper)
Definition: ecfvstencil.hh:264
const Element & element() const
Return the element to which the stencil refers.
Definition: ecfvstencil.hh:320
const Entity & entity(unsigned dofIdx) const
Return the entity given the index of a degree of freedom.
Definition: ecfvstencil.hh:384
std::vector< SubControlVolume > subControlVolumes_
Definition: ecfvstencil.hh:425
void updateTopology(const Element &element)
Definition: ecfvstencil.hh:272
std::size_t numDof() const
Returns the number of degrees of freedom which the current element interacts with.
Definition: ecfvstencil.hh:334
unsigned globalSpaceIndex(unsigned dofIdx) const
Return the global space index given the index of a degree of freedom.
Definition: ecfvstencil.hh:353
Element Entity
Definition: ecfvstencil.hh:83
typename Element::Geometry LocalGeometry
Definition: ecfvstencil.hh:86
void update(const Element &element)
Definition: ecfvstencil.hh:307
ElementMapper Mapper
Definition: ecfvstencil.hh:84
std::size_t numPrimaryDof() const
Returns the number of degrees of freedom which are contained by within the current element.
Definition: ecfvstencil.hh:346
const GridView & gridView() const
Return the grid view of the element to which the stencil refers.
Definition: ecfvstencil.hh:327
const ElementMapper & elementMapper_
Definition: ecfvstencil.hh:422
void updateCenterGradients()
Definition: ecfvstencil.hh:312
Dune::PartitionType partitionType(unsigned dofIdx) const
Return partition type of a given degree of freedom.
Definition: ecfvstencil.hh:363
void updatePrimaryTopology(const Element &element)
Definition: ecfvstencil.hh:298
std::size_t numInteriorFaces() const
Returns the number of interior faces of the stencil.
Definition: ecfvstencil.hh:397
std::size_t numBoundaryFaces() const
Returns the number of boundary faces of the stencil.
Definition: ecfvstencil.hh:410
std::vector< Element > elements_
Definition: ecfvstencil.hh:424
const SubControlVolume & subControlVolume(unsigned dofIdx) const
Returns the sub-control volume object belonging to a given degree of freedom.
Definition: ecfvstencil.hh:391
const GridView & gridView_
Definition: ecfvstencil.hh:421
Definition: blackoilboundaryratevector.hh:39
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)