My Project
GridInterfaceEuler.hpp
1 //===========================================================================
2 //
3 // File: GridInterfaceEuler.hpp
4 //
5 // Created: Mon Jun 15 12:53:38 2009
6 //
7 // Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8 // Bård Skaflestad <bard.skaflestad@sintef.no>
9 //
10 // $Date$
11 //
12 // $Revision$
13 //
14 //===========================================================================
15 
16 /*
17  Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
18  Copyright 2009, 2010 Statoil ASA.
19 
20  This file is part of The Open Reservoir Simulator Project (OpenRS).
21 
22  OpenRS is free software: you can redistribute it and/or modify
23  it under the terms of the GNU General Public License as published by
24  the Free Software Foundation, either version 3 of the License, or
25  (at your option) any later version.
26 
27  OpenRS is distributed in the hope that it will be useful,
28  but WITHOUT ANY WARRANTY; without even the implied warranty of
29  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  GNU General Public License for more details.
31 
32  You should have received a copy of the GNU General Public License
33  along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
34 */
35 
36 #ifndef OPENRS_GRIDINTERFACEEULER_HEADER
37 #define OPENRS_GRIDINTERFACEEULER_HEADER
38 
39 #include <opm/grid/utility/SparseTable.hpp>
40 #include <opm/grid/utility/StopWatch.hpp>
41 #include <opm/grid/CpGrid.hpp> // How to avoid this? Needed for the explicit mapper specialization below.
42 
43 #include <opm/common/utility/platform_dependent/disable_warnings.h>
44 
45 #include <dune/common/fvector.hh>
46 #include <dune/grid/common/mcmgmapper.hh>
47 #include <dune/grid/common/defaultgridview.hh>
48 
49 #include <boost/iterator/iterator_facade.hpp>
50 
51 #include <opm/common/utility/platform_dependent/reenable_warnings.h>
52 
53 #include <climits>
54 #include <iostream>
55 #include <memory>
56 
57 
58 
59 namespace Opm
60 {
61 
63  template <class DuneGrid>
64  struct GICellMapper
65  {
66  using Type = Dune::MultipleCodimMultipleGeomTypeMapper<Dune::GridView<Dune::DefaultLeafGridViewTraits<DuneGrid>>>;
67  };
68 
71  {
72  explicit CpGridCellMapper(const Dune::CpGrid&)
73  {
74  }
75  template<class EntityType>
76  int map (const EntityType& e) const
77  {
78  return e.index();
79  }
80  };
81 
82  namespace GIE
83  {
84  template <class GridInterface, class EntityPointerType>
85  class Cell;
86 
87  template <class GI>
88  class Face {
89  public:
90  typedef typename GI::DuneIntersectionIterator DuneIntersectionIter;
91  typedef typename GI::GridType::Traits::template Codim<0>::EntityPointer CellPtr;
92  typedef typename GI::GridType::ctype Scalar;
93  typedef Dune::FieldVector<Scalar, GI::GridType::dimension> Vector;
94  typedef Dune::FieldVector<Scalar, GI::GridType::dimension-1> LocalVector;
95  typedef int Index;
97 
98  enum { BoundaryMarkerIndex = -999,
99  LocalEndIndex = INT_MAX };
100 
101  Face()
102  : pgrid_(0), iter_(), local_index_(-1)
103  {
104  }
105 
106  Face(const GI& grid,
107  const DuneIntersectionIter& it,
108  const Index loc_ind)
109  : pgrid_(&grid), iter_(it), local_index_(loc_ind)
110  {
111  }
112 
116  Scalar area() const
117  {
118  return iter_->geometry().volume();
119  }
120 
124  Vector centroid() const
125  {
126  //return iter_->geometry().global(localCentroid());
127  return iter_->geometry().center();
128  }
129 
133  Vector normal() const
134  {
135  //return iter_->unitOuterNormal(localCentroid());
136  return iter_->centerUnitOuterNormal();
137  }
138 
142  bool boundary() const
143  {
144  return iter_->boundary();
145  }
146 
150  int boundaryId() const
151  {
152  return iter_->boundaryId();
153  }
154 
158  Cell cell() const
159  {
160  return Cell(*pgrid_, iter_->inside());
161  }
162 
166  Index cellIndex() const
167  {
168  return pgrid_->mapper().index(*iter_->inside());
169  }
170 
175  {
176  return Cell(*pgrid_, iter_->outside());
177  }
178 
182  Index neighbourCellIndex() const
183  {
184  if (iter_->boundary()) {
185  return BoundaryMarkerIndex;
186  } else {
187  return pgrid_->mapper().index(*iter_->outside());
188  }
189  }
190 
194  Index index() const
195  {
196  return pgrid_->faceIndex(cellIndex(), localIndex());
197  }
198 
202  Index localIndex() const
203  {
204  return local_index_;
205  }
206 
210  Scalar neighbourCellVolume() const
211  {
212  return iter_->outside()->geometry().volume();
213  }
214 
215  protected:
216  const GI* pgrid_;
217  DuneIntersectionIter iter_;
218  Index local_index_;
219 
220  private:
221 // LocalVector localCentroid() const
222 // {
223 // typedef Dune::ReferenceElements<Scalar, GI::GridType::dimension-1> RefElems;
224 // return RefElems::general(iter_->type()).position(0,0);
225 // }
226  };
227 
228 
234  template <class GridInterface>
235  class FaceIterator :
236  public boost::iterator_facade<FaceIterator<GridInterface>,
237  const Face<GridInterface>,
238  boost::forward_traversal_tag>,
239  public Face<GridInterface>
240  {
241  private:
243 
244  public:
249 
252  : FaceType()
253  {
254  }
255 
267  FaceIterator(const GridInterface& grid,
268  const DuneIntersectionIter& it,
269  const int local_index)
270  : FaceType(grid, it, local_index)
271  {
272  }
273 
275  const FaceIterator& dereference() const
276  {
277  return *this;
278  }
279 
281  bool equal(const FaceIterator& other) const
282  {
283  // Note that we do not compare the local_index_ members,
284  // since they may or may not be equal for end iterators.
285  return FaceType::iter_ == other.FaceType::iter_;
286  }
287 
289  void increment()
290  {
291  ++FaceType::iter_;
292  ++FaceType::local_index_;
293  }
294 
296  bool operator<(const FaceIterator& other) const
297  {
298  if (FaceType::cellIndex() == other.FaceType::cellIndex()) {
299  return FaceType::localIndex() < other.FaceType::localIndex();
300  } else {
301  return FaceType::cellIndex() < other.FaceType::cellIndex();
302  }
303  }
304  };
305 
306 
307  template <class GridInterface, class EntityPointerType>
308  class Cell
309  {
310  public:
311  Cell()
312  : pgrid_(0), iter_()
313  {
314  }
315  Cell(const GridInterface& grid,
316  const EntityPointerType& it)
317  : pgrid_(&grid), iter_(it)
318  {
319  }
321  typedef typename FaceIterator::Vector Vector;
322  typedef typename FaceIterator::Scalar Scalar;
323  typedef typename FaceIterator::Index Index;
324 
325  FaceIterator facebegin() const
326  {
327  return FaceIterator(*pgrid_, iter_->ileafbegin(), 0);
328  }
329 
330  FaceIterator faceend() const
331  {
332  return FaceIterator(*pgrid_, iter_->ileafend(),
333  FaceIterator::LocalEndIndex);
334  }
335 
336  Scalar volume() const
337  {
338  return iter_->geometry().volume();
339  }
340 
341  Vector centroid() const
342  {
343 // typedef Dune::ReferenceElements<Scalar, GridInterface::GridType::dimension> RefElems;
344 // Vector localpt
345 // = RefElems::general(iter_->type()).position(0,0);
346 // return iter_->geometry().global(localpt);
347  return iter_->geometry().center();
348  }
349 
350  Index index() const
351  {
352  return pgrid_->mapper().index(*iter_);
353  }
354  protected:
355  const GridInterface* pgrid_;
356  EntityPointerType iter_;
357  };
358 
359 
360  template <class GridInterface>
362  : public boost::iterator_facade<CellIterator<GridInterface>,
363  const Cell<GridInterface,
364  typename GridInterface::GridType::template Codim<0>::LeafIterator>,
365  boost::forward_traversal_tag>,
366  public Cell<GridInterface, typename GridInterface::GridType::template Codim<0>::LeafIterator>
367  {
368  private:
369  typedef typename GridInterface::GridType::template Codim<0>::LeafIterator DuneCellIter;
371  public:
372  typedef typename CellType::Vector Vector;
373  typedef typename CellType::Scalar Scalar;
374  typedef typename CellType::Index Index;
375 
376  CellIterator(const GridInterface& grid, DuneCellIter it)
377  : CellType(grid, it)
378  {
379  }
381  const CellIterator& dereference() const
382  {
383  return *this;
384  }
386  bool equal(const CellIterator& other) const
387  {
388  return CellType::iter_ == other.CellType::iter_;
389  }
391  void increment()
392  {
393  ++CellType::iter_;
394  }
395  };
396 
397  } // namespace GIE
398 
399 
400  template <class DuneGrid>
402  {
403  public:
404  typedef typename DuneGrid::LeafIntersectionIterator DuneIntersectionIterator;
405  typedef DuneGrid GridType;
408  typedef typename CellIterator::Vector Vector;
409  typedef typename CellIterator::Scalar Scalar;
410  typedef typename CellIterator::Index Index;
411 
412  typedef typename GICellMapper<DuneGrid>::Type Mapper;
413 
414  enum { Dimension = DuneGrid::dimension };
415 
417  : pgrid_(0), num_faces_(0), max_faces_per_cell_(0)
418  {
419  }
420  explicit GridInterfaceEuler(const DuneGrid& grid, bool build_facemap = true)
421  : pgrid_(&grid), pmapper_(new Mapper(grid.leafGridView(), Dune::mcmgElementLayout())), num_faces_(0), max_faces_per_cell_(0)
422  {
423  if (build_facemap) {
424  buildFaceIndices();
425  }
426  }
427  void init(const DuneGrid& grid, bool build_facemap = true)
428  {
429  pgrid_ = &grid;
430  pmapper_.reset(new Mapper(grid.leafGridView(), Dune::mcmgElementLayout()));
431  if (build_facemap) {
432  buildFaceIndices();
433  }
434  }
435  CellIterator cellbegin() const
436  {
437  return CellIterator(*this, grid().template leafbegin<0>());
438  }
439  CellIterator cellend() const
440  {
441  return CellIterator(*this, grid().template leafend<0>());
442  }
443  int numberOfCells() const
444  {
445  return grid().size(0);
446  }
447  int numberOfFaces() const
448  {
449  assert(num_faces_ != 0);
450  return num_faces_;
451  }
452  int maxFacesPerCell() const
453  {
454  assert(max_faces_per_cell_ != 0);
455  return max_faces_per_cell_;
456  }
457  const DuneGrid& grid() const
458  {
459  assert(pgrid_);
460  return *pgrid_;
461  }
462 
463  // The following are primarily helpers for the implementation,
464  // perhaps they should be private?
465  const Mapper& mapper() const
466  {
467  assert (pmapper_);
468  return *pmapper_;
469  }
470  Index faceIndex(int cell_index, int local_face_index) const
471  {
472  assert(num_faces_ != 0);
473  return face_indices_[cell_index][local_face_index];
474  }
475  typedef Opm::SparseTable<int>::row_type Indices;
476  Indices faceIndices(int cell_index) const
477  {
478  assert(num_faces_ != 0);
479  return face_indices_[cell_index];
480  }
481  private:
482  const DuneGrid* pgrid_;
483  std::unique_ptr<Mapper> pmapper_;
484  int num_faces_;
485  int max_faces_per_cell_;
486  Opm::SparseTable<int> face_indices_;
487 
488  void buildFaceIndices()
489  {
490 #ifdef VERBOSE
491  std::cout << "Building unique face indices... " << std::flush;
492  Opm::time::StopWatch clock;
493  clock.start();
494 #endif
495  typedef CellIterator CI;
496  typedef typename CI::FaceIterator FI;
497 
498  // We build the actual cell to face mapping in two passes.
499  // [code mostly lifted from IncompFlowSolverHybrid::enumerateGridDof(),
500  // but with a twist: This code builds a mapping from cells in index
501  // order to unique face numbers, while the mapping built in the
502  // enumerateGridDof() method was ordered by cell iterator order]
503 
504  // Allocate and reserve structures.
505  const int nc = numberOfCells();
506  std::vector<int> cell(nc, -1);
507  std::vector<int> num_faces(nc); // In index order.
508  std::vector<int> fpos; fpos .reserve(nc + 1);
509  std::vector<int> num_cf; num_cf.reserve(nc); // In iterator order.
510  std::vector<int> faces ;
511 
512  // First pass: enumerate internal faces.
513  int cellno = 0; fpos.push_back(0);
514  int tot_ncf = 0, tot_ncf2 = 0, max_ncf = 0;
515  for (CI c = cellbegin(); c != cellend(); ++c, ++cellno) {
516  const int c0 = c->index();
517  assert((0 <= c0) && (c0 < nc) && (cell[c0] == -1));
518  cell[c0] = cellno;
519  num_cf.push_back(0);
520  int& ncf = num_cf.back();
521  for (FI f = c->facebegin(); f != c-> faceend(); ++f) {
522  if (!f->boundary()) {
523  const int c1 = f->neighbourCellIndex();
524  assert((0 <= c1) && (c1 < nc) && (c1 != c0));
525 
526  if (cell[c1] == -1) {
527  // Previously undiscovered internal face.
528  faces.push_back(c1);
529  }
530  }
531  ++ncf;
532  }
533  num_faces[c0] = ncf;
534  fpos.push_back(int(faces.size()));
535  max_ncf = std::max(max_ncf, ncf);
536  tot_ncf += ncf;
537  tot_ncf2 += ncf * ncf;
538  }
539  assert(cellno == nc);
540 
541  // Build cumulative face sizes enabling direct insertion of
542  // face indices into cfdata later.
543  std::vector<int> cumul_num_faces(numberOfCells() + 1);
544  cumul_num_faces[0] = 0;
545  std::partial_sum(num_faces.begin(), num_faces.end(), cumul_num_faces.begin() + 1);
546 
547  // Avoid (most) allocation(s) inside 'c' loop.
548  std::vector<int> l2g;
549  l2g.reserve(max_ncf);
550  std::vector<double> cfdata(tot_ncf);
551  int total_num_faces = int(faces.size());
552 
553  // Second pass: build cell-to-face mapping, including boundary.
554  typedef std::vector<int>::iterator VII;
555  for (CI c = cellbegin(); c != cellend(); ++c) {
556  const int c0 = c->index();
557  assert ((0 <= c0 ) && ( c0 < nc) &&
558  (0 <= cell[c0]) && (cell[c0] < nc));
559  const int ncf = num_cf[cell[c0]];
560  l2g.resize(ncf, 0);
561  for (FI f = c->facebegin(); f != c->faceend(); ++f) {
562  if (f->boundary()) {
563  // External, not counted before. Add new face...
564  l2g[f->localIndex()] = total_num_faces++;
565  } else {
566  // Internal face. Need to determine during
567  // traversal of which cell we discovered this
568  // face first, and extract the face number
569  // from the 'faces' table range of that cell.
570 
571  // Note: std::find() below is potentially
572  // *VERY* expensive (e.g., large number of
573  // seeks in moderately sized data in case of
574  // faulted cells).
575  const int c1 = f->neighbourCellIndex();
576  assert ((0 <= c1 ) && ( c1 < nc) &&
577  (0 <= cell[c1]) && (cell[c1] < nc));
578 
579  int t = c0, seek = c1;
580  if (cell[seek] < cell[t])
581  std::swap(t, seek);
582  int s = fpos[cell[t]], e = fpos[cell[t] + 1];
583  VII p = std::find(faces.begin() + s, faces.begin() + e, seek);
584  assert(p != faces.begin() + e);
585  l2g[f->localIndex()] = p - faces.begin();
586  }
587  }
588  assert(int(l2g.size()) == num_faces[c0]);
589  std::copy(l2g.begin(), l2g.end(), cfdata.begin() + cumul_num_faces[c0]);
590  }
591  num_faces_ = total_num_faces;
592  max_faces_per_cell_ = max_ncf;
593  face_indices_.assign(cfdata.begin(), cfdata.end(),
594  num_faces.begin(), num_faces.end());
595 
596 #ifdef VERBOSE
597  clock.stop();
598  double elapsed = clock.secsSinceStart();
599  std::cout << "done. Time elapsed: " << elapsed << std::endl;
600 #endif
601  }
602 
603  };
604 
605 
606 
607 } // namespace Opm
608 
609 
610 #endif // OPENRS_GRIDINTERFACEEULER_HEADER
Definition: GridInterfaceEuler.hpp:367
const CellIterator & dereference() const
Used by iterator facade.
Definition: GridInterfaceEuler.hpp:381
bool equal(const CellIterator &other) const
Used by iterator facade.
Definition: GridInterfaceEuler.hpp:386
void increment()
Used by iterator facade.
Definition: GridInterfaceEuler.hpp:391
Definition: GridInterfaceEuler.hpp:309
Intersection (face) iterator for solver-near grid interface.
Definition: GridInterfaceEuler.hpp:240
const FaceIterator & dereference() const
Used by iterator facade.
Definition: GridInterfaceEuler.hpp:275
bool equal(const FaceIterator &other) const
Used by iterator facade.
Definition: GridInterfaceEuler.hpp:281
void increment()
Used by iterator facade.
Definition: GridInterfaceEuler.hpp:289
Face< GridInterface >::DuneIntersectionIter DuneIntersectionIter
Type of low-level intersection iterator.
Definition: GridInterfaceEuler.hpp:248
bool operator<(const FaceIterator &other) const
Gives an ordering of intersectionIterators.
Definition: GridInterfaceEuler.hpp:296
FaceIterator(const GridInterface &grid, const DuneIntersectionIter &it, const int local_index)
Constructor.
Definition: GridInterfaceEuler.hpp:267
FaceIterator()
Default constructor.
Definition: GridInterfaceEuler.hpp:251
Definition: GridInterfaceEuler.hpp:88
Vector normal() const
Definition: GridInterfaceEuler.hpp:133
Vector centroid() const
Definition: GridInterfaceEuler.hpp:124
Index neighbourCellIndex() const
Definition: GridInterfaceEuler.hpp:182
bool boundary() const
Definition: GridInterfaceEuler.hpp:142
Cell cell() const
Definition: GridInterfaceEuler.hpp:158
Index cellIndex() const
Definition: GridInterfaceEuler.hpp:166
Cell neighbourCell() const
Definition: GridInterfaceEuler.hpp:174
int boundaryId() const
Definition: GridInterfaceEuler.hpp:150
Index localIndex() const
Definition: GridInterfaceEuler.hpp:202
Scalar neighbourCellVolume() const
Definition: GridInterfaceEuler.hpp:210
Scalar area() const
Definition: GridInterfaceEuler.hpp:116
Index index() const
Definition: GridInterfaceEuler.hpp:194
Definition: GridInterfaceEuler.hpp:402
Inverting small matrices.
Definition: ImplicitAssembly.hpp:43
A mapper for Dune::CpGrid cells only.
Definition: GridInterfaceEuler.hpp:71
Mapper for general grids.
Definition: GridInterfaceEuler.hpp:65