Dune Core Modules (2.9.0)

ygrid.hh
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4 // vi: set et ts=4 sw=2 sts=2:
5 #ifndef DUNE_GRID_YASPGRID_YGRID_HH
6 #define DUNE_GRID_YASPGRID_YGRID_HH
7 
8 #include <array>
9 #include <vector>
10 #include <bitset>
11 #include <deque>
12 
13 #include <dune/common/fvector.hh>
14 #include <dune/common/math.hh>
16 
21 namespace Dune {
22 
23  namespace Yasp {
28  template<int d, typename ct>
29  std::array<int,d> sizeArray(const std::array<std::vector<ct>,d>& v)
30  {
31  std::array<int,d> tmp;
32  for (int i=0; i<d; ++i)
33  tmp[i] = v[i].size() - 1;
34  return tmp;
35  }
36  } //namespace Yasp
37 
73  template<class Coordinates>
75  {
76  public:
77  //extract coordinate type and dimension from the coordinate container
78  typedef typename Coordinates::ctype ct;
79  static const int d = Coordinates::dimension;
80 
81  typedef std::array<int, d> iTupel;
82  typedef FieldVector<ct,d> fTupel;
83 
85  YGridComponent () : _shift(0ULL)
86  {
87  std::fill(_origin.begin(), _origin.end(), 0);
88  std::fill(_offset.begin(), _offset.end(), 0);
89  std::fill(_size.begin(), _size.end(), 0);
90  }
91 
99  YGridComponent(iTupel origin, iTupel size)
100  : _origin(origin), _size(size)
101  {}
102 
108  YGridComponent (iTupel origin, iTupel size, const YGridComponent<Coordinates>& enclosing)
109  : _origin(origin), _shift(enclosing.shift()), _coords(enclosing.getCoords()), _size(size), _supersize(enclosing.supersize())
110  {
111  for (int i=0; i<d; i++)
112  _offset[i] = origin[i] - enclosing.origin(i) + enclosing.offset(i);
113 
114  // compute superincrements
115  int inc = 1;
116  for (int i=0; i<d; ++i)
117  {
118  _superincrement[i] = inc;
119  inc *= _supersize[i];
120  }
121  }
122 
131  YGridComponent (iTupel origin, std::bitset<d> shift, Coordinates* coords, iTupel size, iTupel offset, iTupel supersize)
132  : _origin(origin), _shift(shift), _coords(coords), _size(size), _offset(offset), _supersize(supersize)
133  {
134  // compute superincrements
135  int inc = 1;
136  for (int i=0; i<d; ++i)
137  {
138  _superincrement[i] = inc;
139  inc *= _supersize[i];
140  }
141  }
142 
144  int origin (int i) const
145  {
146  return _origin[i];
147  }
148 
150  const iTupel& origin () const
151  {
152  return _origin;
153  }
154 
156  bool shift (int i) const
157  {
158  return _shift[i];
159  }
160 
162  const std::bitset<d>& shift () const
163  {
164  return _shift;
165  }
166 
167  Coordinates* getCoords() const
168  {
169  return _coords;
170  }
171 
173  int offset (int i) const
174  {
175  return _offset[i];
176  }
177 
179  const iTupel & offset () const
180  {
181  return _offset;
182  }
183 
185  int supersize (int i) const
186  {
187  return _supersize[i];
188  }
189 
191  const iTupel & supersize () const
192  {
193  return _supersize;
194  }
195 
197  int size (int i) const
198  {
199  return _size[i];
200  }
201 
203  iTupel size () const
204  {
205  return _size;
206  }
207 
209  int totalsize () const
210  {
211  int s=1;
212  for (int i=0; i<d; ++i)
213  s *= size(i);
214  return s;
215  }
216 
218  int min (int i) const
219  {
220  return _origin[i];
221  }
222 
224  int max (int i) const
225  {
226  return _origin[i] + size(i) - 1;
227  }
228 
230  bool empty () const
231  {
232  for (int i=0; i<d; ++i)
233  {
234  if (size(i) == 0)
235  return true;
236  }
237  return false;
238  }
239 
241  bool inside (const iTupel& coord) const
242  {
243  for (int i=0; i<d; i++)
244  {
245  if ((coord[i]<_origin[i]) || (coord[i]>=_origin[i]+_size[i]))
246  return false;
247  }
248  return true;
249  }
250 
252  int index (const iTupel& coord) const
253  {
254  int index = (coord[d-1]-_origin[d-1]);
255 
256  for (int i=d-2; i>=0; i--)
257  index = index*_size[i] + (coord[i]-_origin[i]);
258 
259  return index;
260  }
261 
264  {
265  for (int i=0; i<d; i++)
266  v[i] += _origin[i];
267  return YGridComponent<Coordinates>(v,_size,*this);
268  }
269 
272  {
273  for (int i=0; i<d; i++)
274  {
275  //empty coordinate vectors result in empty intersections
276  if (empty() || r.empty())
278  }
279 
280  iTupel neworigin;
281  iTupel newsize;
282  for (int i=0; i<d; ++i)
283  {
284  neworigin[i] = std::max(origin(i),r.origin(i));
285  newsize[i] = std::min(max(i),r.max(i)) - neworigin[i] + 1;
286  }
287 
288  return YGridComponent<Coordinates>(neworigin,newsize,*this);
289  }
290 
291 
298  class Iterator {
299  public:
300  // default constructor
301  Iterator () = default;
302 
304  Iterator (const YGridComponent<Coordinates>& r) : _grid(&r)
305  {
306  iTupel coord(r.origin());
307  reinit(r,coord);
308  }
309 
311  Iterator (const YGridComponent<Coordinates>& r, const iTupel& coord)
312  {
313  reinit(r,coord);
314  }
315 
317  void reinit (const YGridComponent<Coordinates>& r, const iTupel& coord)
318  {
319  // initialize to given position in index set
320  for (int i=0; i<d; ++i)
321  _coord[i] = coord[i];
322 
323  // move superindex to first cell in subgrid
324  _superindex = 0;
325  for (int i=0; i<d; ++i)
326  _superindex += (r.offset(i)+coord[i]-r.origin(i))*r.superincrement(i);
327 
328  _grid = &r;
329  }
330 
332  bool operator== (const Iterator& i) const
333  {
334  return _superindex == i._superindex;
335  }
336 
338  bool operator!= (const Iterator& i) const
339  {
340  return _superindex != i._superindex;
341  }
342 
344  int superindex () const
345  {
346  return _superindex;
347  }
348 
350  int coord (int i) const
351  {
352  return _coord[i];
353  }
354 
356  const iTupel& coord () const
357  {
358  return _coord;
359  }
360 
362  void move (int i, int dist)
363  {
364  _coord[i] += dist;
365  _superindex += dist*_grid->superincrement(i);
366  }
367 
369  void move (const iTupel& dist)
370  {
371  for (int i = 0; i < d; ++i)
372  {
373  _coord[i] += dist[i];
374  _superindex += dist[i]*_grid->superincrement(i);
375  }
376  }
377 
380  {
381  for (int i=0; i<d; i++) // check for wrap around
382  {
383  _superindex += _grid->superincrement(i); // move on cell in direction i
384  if (++_coord[i] <= _grid->max(i))
385  return *this;
386  else
387  {
388  _coord[i] = _grid->origin(i); // move back to origin in direction i
389  _superindex -= _grid->size(i) * _grid->superincrement(i);
390  }
391  }
392  // if we wrapped around, back to to begin(), we must put the iterator to end()
393  if (_coord == _grid->origin())
394  {
395  for (int i=0; i<d; i++)
396  _superindex += (_grid->size(i)-1) * _grid->superincrement(i);
397  _superindex += _grid->superincrement(0);
398  }
399  return *this;
400  }
401 
403  ct lowerleft(int i) const
404  {
405  return _grid->getCoords()->coordinate(i,_coord[i]);
406  }
407 
410  {
411  fTupel ll;
412  for (int i=0; i<d; i++)
413  ll[i] = lowerleft(i);
414  return ll;
415  }
416 
418  ct upperright(int i) const
419  {
420  int coord = _coord[i];
421  if (shift(i))
422  coord++;
423  return _grid->getCoords()->coordinate(i,coord);
424  }
425 
428  {
429  fTupel ur;
430  for (int i=0; i<d; i++)
431  ur[i] = upperright(i);
432  return ur;
433  }
434 
436  ct meshsize (int i) const
437  {
438  return _grid->getCoords()->meshsize(i,_coord[i]);
439  }
440 
442  fTupel meshsize () const
443  {
444  fTupel h;
445  for (int i=0; i<d; i++)
446  h[i] = meshsize(i);
447  return h;
448  }
449 
450  bool shift (int i) const
451  {
452  return _grid->shift(i);
453  }
454 
455  std::bitset<d> shift() const
456  {
457  return _grid->shift();
458  }
459 
460  Coordinates* coordCont() const
461  {
462  return _grid->getCoords();
463  }
464 
465  protected:
466  iTupel _coord;
467  int _superindex = 0;
468  const YGridComponent<Coordinates>* _grid = nullptr;
469  };
470 
471 
472  int superindex(iTupel coord) const
473  {
474  // move superindex to first cell in subgrid
475  int si = 0;
476  for (int i=0; i<d; ++i)
477  si += (offset(i)+coord[i]-origin(i))*_superincrement[i];
478  return si;
479  }
480 
481  int superincrement(int i) const
482  {
483  return _superincrement[i];
484  }
485 
487  Iterator begin () const
488  {
489  return Iterator(*this);
490  }
491 
493  Iterator begin (const iTupel& co) const
494  {
495  return Iterator(*this,co);
496  }
497 
499  Iterator end () const
500  {
501  iTupel last;
502  for (int i=0; i<d; i++)
503  last[i] = max(i);
504  last[0] += 1;
505  return Iterator(*this,last);
506  }
507 
508  private:
509  iTupel _origin;
510  std::bitset<d> _shift;
511  Coordinates* _coords;
512  iTupel _size;
513  iTupel _offset;
514  iTupel _supersize;
515  iTupel _superincrement;
516 
517  };
518 
519 
521  template <class Coordinates>
522  inline std::ostream& operator<< (std::ostream& s, YGridComponent<Coordinates> e)
523  {
524  s << "Printing YGridComponent structure:" << std::endl;
525  s << "Origin: " << e.origin() << std::endl;
526  s << "Shift: " << e.shift() << std::endl;
527  s << "Size: " << e.size() << std::endl;
528  s << "Offset: " << e.offset() << std::endl;
529  s << "Supersize: " << e.supersize() << std::endl;
530  return s;
531  }
532 
534  template <class Coordinates>
535  inline std::ostream& operator<< (std::ostream& s, typename YGridComponent<Coordinates>::Iterator& e)
536  {
537  s << "Printing YGridComponent Iterator:" << std::endl << "Iterator at " << e.coord() << " (index ";
538  s << e.index() << "), position " << e.position();
539  return s;
540  }
541 
549  template<class Coordinates>
550  class YGrid
551  {
552  public:
553  static const int dim = Coordinates::dimension;
554 
555  // define data array iterator
557 
558  typedef typename std::array<int, dim> iTupel;
559 
562  {
563  _begin = begin;
564  }
565 
567  int shiftmapping(const std::bitset<dim>& shift) const
568  {
569  return _shiftmapping[shift.to_ulong()];
570  }
571 
573  DAI dataBegin() const
574  {
575  return _begin;
576  }
577 
579  DAI dataEnd() const
580  {
581  return _end;
582  }
583 
585  bool inside(const iTupel& coord, const std::bitset<dim>& shift = std::bitset<dim>()) const
586  {
587  return (_begin+_shiftmapping[shift.to_ulong()])->inside(coord);
588  }
589 
593  class Iterator
594  {
595  public:
596 
598  Iterator () = default;
599 
601  Iterator (const YGrid<Coordinates>& yg, const std::array<int,dim>& coords, int which = 0)
602  : _which(which), _yg(&yg)
603  {
604  _it = typename YGridComponent<Coordinates>::Iterator(*(_yg->dataBegin()+which),coords);
605  }
606 
608  Iterator (const YGrid<Coordinates>& yg, bool end=false) : _yg(&yg)
609  {
610  if (end)
611  {
612  _it = _yg->_itends.back();
613  _which = _yg->_itends.size() - 1;
614  }
615  else
616  {
617  _it = _yg->_itbegins[0];
618  _which = 0;
619  }
620  }
621 
623  void reinit(const YGrid<Coordinates>& yg, const std::array<int,dim>& coords, int which = 0)
624  {
625  _yg = &yg;
626  _which = which;
627  _it = typename YGridComponent<Coordinates>::Iterator(*(_yg->dataBegin()+which),coords);
628  }
629 
631  int coord (int i) const
632  {
633  return _it.coord(i);
634  }
635 
637  const std::array<int, dim>& coord () const
638  {
639  return _it.coord();
640  }
641 
642  typename Coordinates::ctype lowerleft(int i) const
643  {
644  return _it.lowerleft(i);
645  }
646 
648  {
649  return _it.lowerleft();
650  }
651 
652  typename Coordinates::ctype upperright(int i) const
653  {
654  return _it.upperright(i);
655  }
656 
658  {
659  return _it.upperright();
660  }
661 
663  typename Coordinates::ctype meshsize (int i) const
664  {
665  return _it.meshsize(i);
666  }
667 
670  {
671  return _it.meshsize();
672  }
673 
675  bool shift (int i) const
676  {
677  return _it.shift(i);
678  }
679 
681  std::bitset<dim> shift () const
682  {
683  return _it.shift();
684  }
685 
687  int superindex() const
688  {
689  // the offset of the current component has to be taken into account
690  return _yg->_indexOffset[_which] + _it.superindex();
691  }
692 
695  {
696  if ((++_it == _yg->_itends[_which]) && (_which < _yg->_itends.size()-1))
697  _it = _yg->_itbegins[++_which];
698  return *this;
699  }
700 
702  bool operator==(const Iterator& i) const
703  {
704  if (_which != i._which)
705  return false;
706  return _it == i._it;
707  }
708 
710  bool operator!=(const Iterator& i) const
711  {
712  if (_it != i._it)
713  return true;
714  return _which != i._which;
715  }
716 
718  int which() const
719  {
720  return _which;
721  }
722 
724  void move(int i, int dist)
725  {
726  _it.move(i,dist);
727  }
728 
729  void move(const iTupel& dist)
730  {
731  _it.move(dist);
732  }
733 
734  Coordinates* coordCont() const
735  {
736  return _it.coordCont();
737  }
738 
739 
740  private:
741  unsigned int _which = 0;
742  const YGrid<Coordinates>* _yg = nullptr;
743  typename YGridComponent<Coordinates>::Iterator _it;
744  };
745 
747  Iterator begin() const
748  {
749  return Iterator(*this);
750  }
751 
753  Iterator begin(const std::array<int, dim>& coord, int which = 0) const
754  {
755  return Iterator(*this, coord, which);
756  }
757 
759  Iterator end() const
760  {
761  return Iterator(*this,true);
762  }
763 
764  int superindex(const iTupel& coord, int which) const
765  {
766  return _indexOffset[which] + (dataBegin()+which)->superindex(coord);
767  }
768 
769 
770  // finalize the ygrid construction by storing component iterators
771  void finalize(const DAI& end, int artificialOffset = 0)
772  {
773  // set the end iterator in the ygrid component array
774  _end = end;
775 
776  _indexOffset.push_back(artificialOffset);
777  int k = 0;
778  for (DAI i=_begin; i != _end; ++i, ++k)
779  {
780  //store begin and end iterators
781  _itbegins.push_back(i->begin());
782  _itends.push_back(i->end());
783 
784  // store index offset
785  _indexOffset.push_back(_indexOffset.back() + i->totalsize());
786 
787  // store shift to component mapping
788  _shiftmapping[i->shift().to_ulong()] = k;
789  }
790  _indexOffset.resize(_itends.size());
791  }
792 
793  private:
794 
795  friend class YGrid<Coordinates>::Iterator;
796  DAI _begin;
797  DAI _end;
798  std::array<int,Dune::power(2,dim)> _shiftmapping;
799  std::vector<typename YGridComponent<Coordinates>::Iterator> _itbegins;
800  std::vector<typename YGridComponent<Coordinates>::Iterator> _itends;
801  std::vector<int> _indexOffset;
802  };
803 
805  template <class Coordinates>
806  inline std::ostream& operator<< (std::ostream& s, const YGrid<Coordinates>& e)
807  {
808  s << "Printing YGrid structure:" << std::endl;
809  for (auto it = e.dataBegin(); it != e.dataEnd(); ++it)
810  s << *it << std::endl;
811  return s;
812  }
813 
821  template<class Coordinates>
822  class YGridList
823  {
824  public:
825  static const int dim = Coordinates::dimension;
826 
828  struct Intersection
829  {
831  YGridComponent<Coordinates> grid;
833  int rank;
835  int distance;
837  YGrid<Coordinates> yg;
838  };
839 
840  // define data array iterator type
841  typedef typename std::array<std::deque<Intersection>, Dune::power(2,dim)>::iterator DAI;
842 
843  // iterator that allows to iterate over a concatenation of deques. namely those
844  // that belong to the same codimension.
845  class Iterator
846  {
847  public:
848 
850  Iterator(const YGridList<Coordinates>& ygl, bool end=false) : _end(ygl.dataEnd()), _which(ygl.dataBegin())
851  {
852  _it = _which->begin();
853 
854  // advance the iterator to the first element that exists.
855  // some deques might be empty and should be skipped
856  while ((_which != _end) && (_it == _which->end()))
857  {
858  ++_which;
859  if (_which != _end)
860  _it = _which->begin();
861  }
862  // the iterator is at the end if and only if _which==_end
863  if (end)
864  {
865  _which = _end;
866  }
867  }
868 
870  Iterator& operator++ ()
871  {
872  ++_it;
873  // advance the iterator to the next element that exists.
874  // some deques might be empty and should be skipped
875  while ((_which != _end) && (_it == _which->end()))
876  {
877  ++_which;
878  if (_which != _end)
879  _it = _which->begin();
880  }
881  return *this;
882  }
883 
885  typename std::deque<Intersection>::iterator operator->() const
886  {
887  return _it;
888  }
889 
891  typename std::deque<Intersection>::iterator operator*() const
892  {
893  return _it;
894  }
895 
897  bool operator== (const Iterator& i) const
898  {
899  if (_which != i._which)
900  return false;
901  if (_which == _end)
902  return true;
903  return _it == i._it;
904  }
905 
907  bool operator!= (const Iterator& i) const
908  {
909  if (_which != i._which)
910  return true;
911  if (_which == _end)
912  return false;
913  return _it != i._it;
914  }
915 
916  private:
917  typename std::deque<Intersection>::iterator _it;
918  DAI _end;
919  DAI _which;
920  };
921 
923  Iterator begin() const
924  {
925  return Iterator(*this);
926  }
927 
929  Iterator end() const
930  {
931  return Iterator(*this,true);
932  }
933 
935  void setBegin(typename std::array<std::deque<Intersection>, Dune::power(2,dim)>::iterator begin)
936  {
937  _begin = begin;
938  }
939 
941  DAI dataBegin() const
942  {
943  return _begin;
944  }
945 
947  DAI dataEnd() const
948  {
949  return _end;
950  }
951 
953  int size() const
954  {
955  int count = 0;
956  for (DAI it = _begin; it != _end; ++it)
957  count += it->size();
958  return count;
959  }
960 
962  void finalize(DAI end, const YGrid<Coordinates>& ygrid)
963  {
964  // Instead of directly iterating over the intersection deques, this code
965  // iterates over the components of an associated ygrid and works its way
966  // through the list of intersection deques in parallel.
967  // The reason for this convoluted iteration technique is that there are not
968  // necessarily intersections for all possible shifts, but we have to make
969  // sure that we stop at each shift to update the per-component index shift.
970  // This is ensured by iterating over the ygrid, which is guaranteed to have
971  // a component for each shift vector.
972 
973  // set end iterator in the data array
974  _end = end;
975 
977  int offset = 0;
978 
979  DAI i = _begin;
980 
981  // make sure that we have a valid deque (i.e. a non-empty one)
982  while (i != _end && i->begin() == i->end())
983  ++i;
984 
985  for (auto yit = ygrid.dataBegin(); yit != ygrid.dataEnd(); ++yit)
986  {
987  if (i == _end)
988  break;
989  auto it = i->begin();
990  if (it->grid.shift() == yit->shift())
991  {
992  // iterate over the intersections in the deque and set the offset
993  for (; it != i->end(); ++it)
994  {
995  it->yg.setBegin(&(it->grid));
996  it->yg.finalize(&(it->grid)+1, offset);
997  }
998 
999  // advance to next non-empty deque
1000  ++i;
1001  while (i != _end && i->begin() == i->end())
1002  ++i;
1003  }
1004 
1005  // update the offset from the ygrid component
1006  int add = 1;
1007  for (int j=0; j<dim; j++)
1008  add *= yit->supersize(j);
1009  offset += add;
1010  }
1011  assert (i == end);
1012  }
1013 
1014  private:
1015  DAI _begin;
1016  DAI _end;
1017  };
1018 
1019 } // namespace Dune
1020 
1021 #endif
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Definition: ygrid.hh:298
const iTupel & coord() const
Return coordinate of the cell as reference (do not modify).
Definition: ygrid.hh:356
int superindex() const
Return consecutive index in enclosing grid.
Definition: ygrid.hh:344
bool operator!=(const Iterator &i) const
Return true when two iterators over the same grid are not equal (!).
Definition: ygrid.hh:338
ct lowerleft(int i) const
Return ith component of lower left corner of the entity associated with the current coordinates and s...
Definition: ygrid.hh:403
ct meshsize(int i) const
Return meshsize in direction i.
Definition: ygrid.hh:436
void reinit(const YGridComponent< Coordinates > &r, const iTupel &coord)
reinitialize iterator to given position
Definition: ygrid.hh:317
void move(const iTupel &dist)
move this iterator dist cells in direction i
Definition: ygrid.hh:369
fTupel meshsize() const
Return meshsize of current cell as reference.
Definition: ygrid.hh:442
fTupel lowerleft() const
Return lower left corner of the entity associated with the current coordinates and shift.
Definition: ygrid.hh:409
iTupel _coord
current position in index set
Definition: ygrid.hh:466
void move(int i, int dist)
move this iterator dist cells in direction i
Definition: ygrid.hh:362
Iterator(const YGridComponent< Coordinates > &r)
Make iterator pointing to first cell in a grid.
Definition: ygrid.hh:304
ct upperright(int i) const
Return ith component of upper right corder of the entity associated with the current coordinates and ...
Definition: ygrid.hh:418
int coord(int i) const
Return coordinate of the cell in direction i.
Definition: ygrid.hh:350
fTupel upperright() const
Return upper right corder of the entity associated with the current coordinates and shift.
Definition: ygrid.hh:427
bool operator==(const Iterator &i) const
Return true when two iterators over the same grid are equal (!).
Definition: ygrid.hh:332
Iterator & operator++()
Increment iterator to next cell with position.
Definition: ygrid.hh:379
Iterator(const YGridComponent< Coordinates > &r, const iTupel &coord)
Make iterator pointing to given cell in a grid.
Definition: ygrid.hh:311
int _superindex
consecutive index in enclosing grid
Definition: ygrid.hh:467
Definition: ygrid.hh:75
int index(const iTupel &coord) const
given a tupel compute its index in the lexicographic numbering
Definition: ygrid.hh:252
int offset(int i) const
Return offset to origin of enclosing grid.
Definition: ygrid.hh:173
iTupel size() const
retrun size
Definition: ygrid.hh:203
int min(int i) const
Return minimum index in direction i.
Definition: ygrid.hh:218
const iTupel & supersize() const
return size of enclosing grid
Definition: ygrid.hh:191
YGridComponent(iTupel origin, iTupel size)
make ygrid without coordinate information
Definition: ygrid.hh:99
int totalsize() const
Return total size of index set which is the product of all size per direction.
Definition: ygrid.hh:209
YGridComponent(iTupel origin, iTupel size, const YGridComponent< Coordinates > &enclosing)
make a subgrid by taking coordinates from a larger grid
Definition: ygrid.hh:108
YGridComponent(iTupel origin, std::bitset< d > shift, Coordinates *coords, iTupel size, iTupel offset, iTupel supersize)
Make YGridComponent by giving all parameters.
Definition: ygrid.hh:131
Iterator begin() const
return iterator to first element of index set
Definition: ygrid.hh:487
YGridComponent< Coordinates > move(iTupel v) const
return grid moved by the vector v
Definition: ygrid.hh:263
Iterator begin(const iTupel &co) const
return iterator to given element of index set
Definition: ygrid.hh:493
const iTupel & origin() const
return reference to origin
Definition: ygrid.hh:150
YGridComponent()
make uninitialized ygrid
Definition: ygrid.hh:85
int origin(int i) const
Return origin in direction i.
Definition: ygrid.hh:144
int supersize(int i) const
return size of enclosing grid
Definition: ygrid.hh:185
const std::bitset< d > & shift() const
Return shift tupel.
Definition: ygrid.hh:162
int size(int i) const
return size in direction i
Definition: ygrid.hh:197
int max(int i) const
Return maximum index in direction i.
Definition: ygrid.hh:224
bool shift(int i) const
Return shift in direction i.
Definition: ygrid.hh:156
YGridComponent< Coordinates > intersection(const YGridComponent< Coordinates > &r) const
Return YGridComponent of supergrid of self which is the intersection of self and another YGridCompone...
Definition: ygrid.hh:271
Iterator end() const
return subiterator to last element of index set
Definition: ygrid.hh:499
bool inside(const iTupel &coord) const
given a coordinate, return true if it is in the grid
Definition: ygrid.hh:241
bool empty() const
Return true if YGrid is empty, i.e. has size 0 in all directions.
Definition: ygrid.hh:230
const iTupel & offset() const
Return offset to origin of enclosing grid.
Definition: ygrid.hh:179
Iterator over a collection o YGrids A YGrid::Iterator is the heart of an entity in YaspGrid.
Definition: ygrid.hh:594
Iterator()=default
default constructor
int which() const
return the current component number
Definition: ygrid.hh:718
void move(int i, int dist)
move the grid, this is only done and needed for codim 0
Definition: ygrid.hh:724
bool operator==(const Iterator &i) const
compare two iterators: component has to match
Definition: ygrid.hh:702
Iterator(const YGrid< Coordinates > &yg, const std::array< int, dim > &coords, int which=0)
construct an iterator from coordinates and component
Definition: ygrid.hh:601
std::bitset< dim > shift() const
return the shift vector
Definition: ygrid.hh:681
Dune::FieldVector< typename Coordinates::ctype, dim > meshsize() const
return the current meshsize vector
Definition: ygrid.hh:669
bool shift(int i) const
return the shift in direction i
Definition: ygrid.hh:675
bool operator!=(const Iterator &i) const
compare two iterators: component has to match
Definition: ygrid.hh:710
Iterator & operator++()
increment to the next entity jumping to next component if necessary
Definition: ygrid.hh:694
Coordinates::ctype meshsize(int i) const
return the current meshsize in direction i
Definition: ygrid.hh:663
int coord(int i) const
return coordinate at the current position (direction i)
Definition: ygrid.hh:631
Iterator(const YGrid< Coordinates > &yg, bool end=false)
create an iterator to start or end of the codimension
Definition: ygrid.hh:608
const std::array< int, dim > & coord() const
return coordinate array at the current position
Definition: ygrid.hh:637
int superindex() const
return the superindex
Definition: ygrid.hh:687
void reinit(const YGrid< Coordinates > &yg, const std::array< int, dim > &coords, int which=0)
reinitializes an iterator, as if it was just constructed.
Definition: ygrid.hh:623
implements a collection of YGridComponents which form a codimension Entities of given codimension c n...
Definition: ygrid.hh:551
int shiftmapping(const std::bitset< dim > &shift) const
get which component belongs to a given shift vector
Definition: ygrid.hh:567
DAI dataBegin() const
get start iterator in the data array
Definition: ygrid.hh:573
Iterator begin() const
return begin iterator for the codimension and partition the ygrid represents
Definition: ygrid.hh:747
DAI dataEnd() const
get end iterator in the data array
Definition: ygrid.hh:579
void setBegin(DAI begin)
set start iterator in the data array
Definition: ygrid.hh:561
bool inside(const iTupel &coord, const std::bitset< dim > &shift=std::bitset< dim >()) const
decide whether a coordinate is in the grid (depending on the component)
Definition: ygrid.hh:585
Iterator begin(const std::array< int, dim > &coord, int which=0) const
return iterator pointint to a specified position
Definition: ygrid.hh:753
Iterator end() const
return end iterator for the codimension and partition the ygrid represents
Definition: ygrid.hh:759
Implements a vector constructed from a given type representing a field and a compile-time given size.
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:89
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
Implementation of stream operators for std::array and std::tuple.
std::array< int, d > sizeArray(const std::array< std::vector< ct >, d > &v)
Definition: ygrid.hh:29
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 27, 22:29, 2024)