Dune Core Modules (2.5.0)

ygrid.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_GRID_YASPGRID_YGRID_HH
4#define DUNE_GRID_YASPGRID_YGRID_HH
5
6#include <array>
7#include <vector>
8#include <bitset>
9#include <deque>
10
12#include <dune/common/power.hh>
13
18namespace Dune {
19
20 namespace Yasp {
25 template<int d, typename ct>
26 std::array<int,d> sizeArray(const std::array<std::vector<ct>,d>& v)
27 {
28 std::array<int,d> tmp;
29 for (int i=0; i<d; ++i)
30 tmp[i] = v[i].size() - 1;
31 return tmp;
32 }
33 } //namespace Yasp
34
70 template<class Coordinates>
72 {
73 public:
74 //extract coordinate type and dimension from the coordinate container
75 typedef typename Coordinates::ctype ct;
76 static const int d = Coordinates::dimension;
77
78 typedef std::array<int, d> iTupel;
80
82 YGridComponent () : _shift(0ULL)
83 {
84 std::fill(_origin.begin(), _origin.end(), 0);
85 std::fill(_offset.begin(), _offset.end(), 0);
86 std::fill(_size.begin(), _size.end(), 0);
87 }
88
96 YGridComponent(iTupel origin, iTupel size)
97 : _origin(origin), _size(size)
98 {}
99
105 YGridComponent (iTupel origin, iTupel size, const YGridComponent<Coordinates>& enclosing)
106 : _origin(origin), _shift(enclosing.shift()), _coords(enclosing.getCoords()), _size(size), _supersize(enclosing.supersize())
107 {
108 for (int i=0; i<d; i++)
109 _offset[i] = origin[i] - enclosing.origin(i) + enclosing.offset(i);
110
111 // compute superincrements
112 int inc = 1;
113 for (int i=0; i<d; ++i)
114 {
115 _superincrement[i] = inc;
116 inc *= _supersize[i];
117 }
118 }
119
128 YGridComponent (iTupel origin, std::bitset<d> shift, Coordinates* coords, iTupel size, iTupel offset, iTupel supersize)
129 : _origin(origin), _shift(shift), _coords(coords), _size(size), _offset(offset), _supersize(supersize)
130 {
131 // compute superincrements
132 int inc = 1;
133 for (int i=0; i<d; ++i)
134 {
135 _superincrement[i] = inc;
136 inc *= _supersize[i];
137 }
138 }
139
141 int origin (int i) const
142 {
143 return _origin[i];
144 }
145
147 const iTupel& origin () const
148 {
149 return _origin;
150 }
151
153 bool shift (int i) const
154 {
155 return _shift[i];
156 }
157
159 const std::bitset<d>& shift () const
160 {
161 return _shift;
162 }
163
164 Coordinates* getCoords() const
165 {
166 return _coords;
167 }
168
170 int offset (int i) const
171 {
172 return _offset[i];
173 }
174
176 const iTupel & offset () const
177 {
178 return _offset;
179 }
180
182 int supersize (int i) const
183 {
184 return _supersize[i];
185 }
186
188 const iTupel & supersize () const
189 {
190 return _supersize;
191 }
192
194 int size (int i) const
195 {
196 return _size[i];
197 }
198
200 iTupel size () const
201 {
202 return _size;
203 }
204
206 int totalsize () const
207 {
208 int s=1;
209 for (int i=0; i<d; ++i)
210 s *= size(i);
211 return s;
212 }
213
215 int min (int i) const
216 {
217 return _origin[i];
218 }
219
221 int max (int i) const
222 {
223 return _origin[i] + size(i) - 1;
224 }
225
227 bool empty () const
228 {
229 for (int i=0; i<d; ++i)
230 {
231 if (size(i) == 0)
232 return true;
233 }
234 return false;
235 }
236
238 bool inside (const iTupel& coord) const
239 {
240 for (int i=0; i<d; i++)
241 {
242 if ((coord[i]<_origin[i]) || (coord[i]>=_origin[i]+_size[i]))
243 return false;
244 }
245 return true;
246 }
247
249 int index (const iTupel& coord) const
250 {
251 int index = (coord[d-1]-_origin[d-1]);
252
253 for (int i=d-2; i>=0; i--)
254 index = index*_size[i] + (coord[i]-_origin[i]);
255
256 return index;
257 }
258
261 {
262 for (int i=0; i<d; i++)
263 v[i] += _origin[i];
264 return YGridComponent<Coordinates>(v,_size,*this);
265 }
266
269 {
270 for (int i=0; i<d; i++)
271 {
272 //empty coordinate vectors result in empty intersections
273 if (empty() || r.empty())
275 }
276
277 iTupel neworigin;
278 iTupel newsize;
279 for (int i=0; i<d; ++i)
280 {
281 neworigin[i] = std::max(origin(i),r.origin(i));
282 newsize[i] = std::min(max(i),r.max(i)) - neworigin[i] + 1;
283 }
284
285 return YGridComponent<Coordinates>(neworigin,newsize,*this);
286 }
287
288
295 class Iterator {
296 public:
297 // default constructor
298 Iterator () {}
299
302 {
303 iTupel coord(r.origin());
304 reinit(r,coord);
305 }
306
309 {
310 reinit(r,coord);
311 }
312
314 void reinit (const YGridComponent<Coordinates>& r, const iTupel& coord)
315 {
316 // initialize to given position in index set
317 for (int i=0; i<d; ++i)
318 _coord[i] = coord[i];
319
320 // move superindex to first cell in subgrid
321 _superindex = 0;
322 for (int i=0; i<d; ++i)
323 _superindex += (r.offset(i)+coord[i]-r.origin(i))*r.superincrement(i);
324
325 _grid = &r;
326 }
327
329 bool operator== (const Iterator& i) const
330 {
331 return _superindex == i._superindex;
332 }
333
335 bool operator!= (const Iterator& i) const
336 {
337 return _superindex != i._superindex;
338 }
339
341 int superindex () const
342 {
343 return _superindex;
344 }
345
347 int coord (int i) const
348 {
349 return _coord[i];
350 }
351
353 const iTupel& coord () const
354 {
355 return _coord;
356 }
357
359 void move (int i, int dist)
360 {
361 _coord[i] += dist;
362 _superindex += dist*_grid->superincrement(i);
363 }
364
366 void move (const iTupel& dist)
367 {
368 for (int i = 0; i < d; ++i)
369 {
370 _coord[i] += dist[i];
371 _superindex += dist[i]*_grid->superincrement(i);
372 }
373 }
374
377 {
378 for (int i=0; i<d; i++) // check for wrap around
379 {
380 _superindex += _grid->superincrement(i); // move on cell in direction i
381 if (++_coord[i] <= _grid->max(i))
382 return *this;
383 else
384 {
385 _coord[i] = _grid->origin(i); // move back to origin in direction i
386 _superindex -= _grid->size(i) * _grid->superincrement(i);
387 }
388 }
389 // if we wrapped around, back to to begin(), we must put the iterator to end()
390 if (_coord == _grid->origin())
391 {
392 for (int i=0; i<d; i++)
393 _superindex += (_grid->size(i)-1) * _grid->superincrement(i);
394 _superindex += _grid->superincrement(0);
395 }
396 return *this;
397 }
398
400 ct lowerleft(int i) const
401 {
402 return _grid->getCoords()->coordinate(i,_coord[i]);
403 }
404
407 {
408 fTupel ll;
409 for (int i=0; i<d; i++)
410 ll[i] = lowerleft(i);
411 return ll;
412 }
413
415 ct upperright(int i) const
416 {
417 int coord = _coord[i];
418 if (shift(i))
419 coord++;
420 return _grid->getCoords()->coordinate(i,coord);
421 }
422
425 {
426 fTupel ur;
427 for (int i=0; i<d; i++)
428 ur[i] = upperright(i);
429 return ur;
430 }
431
433 ct meshsize (int i) const
434 {
435 return _grid->getCoords()->meshsize(i,_coord[i]);
436 }
437
440 {
441 fTupel h;
442 for (int i=0; i<d; i++)
443 h[i] = meshsize(i);
444 return h;
445 }
446
447 bool shift (int i) const
448 {
449 return _grid->shift(i);
450 }
451
452 std::bitset<d> shift() const
453 {
454 return _grid->shift();
455 }
456
457 Coordinates* coordCont() const
458 {
459 return _grid->getCoords();
460 }
461
462 protected:
463 iTupel _coord;
465 const YGridComponent<Coordinates>* _grid;
466 };
467
468
469 int superindex(iTupel coord) const
470 {
471 // move superindex to first cell in subgrid
472 int si = 0;
473 for (int i=0; i<d; ++i)
474 si += (offset(i)+coord[i]-origin(i))*_superincrement[i];
475 return si;
476 }
477
478 int superincrement(int i) const
479 {
480 return _superincrement[i];
481 }
482
485 {
486 return Iterator(*this);
487 }
488
490 Iterator begin (const iTupel& co) const
491 {
492 return Iterator(*this,co);
493 }
494
496 Iterator end () const
497 {
498 iTupel last;
499 for (int i=0; i<d; i++)
500 last[i] = max(i);
501 last[0] += 1;
502 return Iterator(*this,last);
503 }
504
505 private:
506 iTupel _origin;
507 std::bitset<d> _shift;
508 Coordinates* _coords;
509 iTupel _size;
510 iTupel _offset;
511 iTupel _supersize;
512 iTupel _superincrement;
513
514 };
515
516
518 template <class Coordinates>
519 inline std::ostream& operator<< (std::ostream& s, YGridComponent<Coordinates> e)
520 {
521 s << "Printing YGridComponent structure:" << std::endl;
522 s << "Origin: " << e.origin() << std::endl;
523 s << "Shift: " << e.shift() << std::endl;
524 s << "Size: " << e.size() << std::endl;
525 s << "Offset: " << e.offset() << std::endl;
526 s << "Supersize: " << e.supersize() << std::endl;
527 return s;
528 }
529
531 template <class Coordinates>
532 inline std::ostream& operator<< (std::ostream& s, typename YGridComponent<Coordinates>::Iterator& e)
533 {
534 s << "Printing YGridComponent Iterator:" << std::endl << "Iterator at " << e.coord() << " (index ";
535 s << e.index() << "), position " << e.position();
536 return s;
537 }
538
546 template<class Coordinates>
547 class YGrid
548 {
549 public:
550 static const int dim = Coordinates::dimension;
551
552 // define data array iterator
554
555 typedef typename std::array<int, dim> iTupel;
556
559 {
560 _begin = begin;
561 }
562
564 int shiftmapping(const std::bitset<dim>& shift) const
565 {
566 return _shiftmapping[shift.to_ulong()];
567 }
568
571 {
572 return _begin;
573 }
574
576 DAI dataEnd() const
577 {
578 return _end;
579 }
580
582 bool inside(const iTupel& coord, const std::bitset<dim>& shift = std::bitset<dim>()) const
583 {
584 return (_begin+_shiftmapping[shift.to_ulong()])->inside(coord);
585 }
586
591 {
592 public:
593
596 {}
597
599 Iterator (const YGrid<Coordinates>& yg, const std::array<int,dim>& coords, int which = 0)
600 : _which(which), _yg(&yg)
601 {
602 _it = typename YGridComponent<Coordinates>::Iterator(*(_yg->dataBegin()+which),coords);
603 }
604
606 Iterator (const YGrid<Coordinates>& yg, bool end=false) : _yg(&yg)
607 {
608 if (end)
609 {
610 _it = _yg->_itends.back();
611 _which = _yg->_itends.size() - 1;
612 }
613 else
614 {
615 _it = _yg->_itbegins[0];
616 _which = 0;
617 }
618 }
619
621 void reinit(const YGrid<Coordinates>& yg, const std::array<int,dim>& coords, int which = 0)
622 {
623 _yg = &yg;
624 _which = which;
625 _it = typename YGridComponent<Coordinates>::Iterator(*(_yg->dataBegin()+which),coords);
626 }
627
629 int coord (int i) const
630 {
631 return _it.coord(i);
632 }
633
635 const std::array<int, dim>& coord () const
636 {
637 return _it.coord();
638 }
639
640 typename Coordinates::ctype lowerleft(int i) const
641 {
642 return _it.lowerleft(i);
643 }
644
646 {
647 return _it.lowerleft();
648 }
649
650 typename Coordinates::ctype upperright(int i) const
651 {
652 return _it.upperright(i);
653 }
654
656 {
657 return _it.upperright();
658 }
659
661 typename Coordinates::ctype meshsize (int i) const
662 {
663 return _it.meshsize(i);
664 }
665
668 {
669 return _it.meshsize();
670 }
671
673 bool shift (int i) const
674 {
675 return _it.shift(i);
676 }
677
679 std::bitset<dim> shift () const
680 {
681 return _it.shift();
682 }
683
685 int superindex() const
686 {
687 // the offset of the current component has to be taken into account
688 return _yg->_indexOffset[_which] + _it.superindex();
689 }
690
693 {
694 if ((++_it == _yg->_itends[_which]) && (_which < _yg->_itends.size()-1))
695 _it = _yg->_itbegins[++_which];
696 return *this;
697 }
698
700 bool operator==(const Iterator& i) const
701 {
702 if (_which != i._which)
703 return false;
704 return _it == i._it;
705 }
706
708 bool operator!=(const Iterator& i) const
709 {
710 if (_it != i._it)
711 return true;
712 return _which != i._which;
713 }
714
716 int which() const
717 {
718 return _which;
719 }
720
722 void move(int i, int dist)
723 {
724 _it.move(i,dist);
725 }
726
727 void move(const iTupel& dist)
728 {
729 _it.move(dist);
730 }
731
732 Coordinates* coordCont() const
733 {
734 return _it.coordCont();
735 }
736
737
738 private:
739 unsigned int _which;
740 const YGrid<Coordinates>* _yg;
741 typename YGridComponent<Coordinates>::Iterator _it;
742 };
743
746 {
747 return Iterator(*this);
748 }
749
751 Iterator begin(const std::array<int, dim>& coord, int which = 0) const
752 {
753 return Iterator(*this, coord, which);
754 }
755
757 Iterator end() const
758 {
759 return Iterator(*this,true);
760 }
761
762 int superindex(const iTupel& coord, int which) const
763 {
764 return _indexOffset[which] + (dataBegin()+which)->superindex(coord);
765 }
766
767
768 // finalize the ygrid construction by storing component iterators
769 void finalize(const DAI& end, int artificialOffset = 0)
770 {
771 // set the end iterator in the ygrid component array
772 _end = end;
773
774 _indexOffset.push_back(artificialOffset);
775 int k = 0;
776 for (DAI i=_begin; i != _end; ++i, ++k)
777 {
778 //store begin and end iterators
779 _itbegins.push_back(i->begin());
780 _itends.push_back(i->end());
781
782 // store index offset
783 _indexOffset.push_back(_indexOffset.back() + i->totalsize());
784
785 // store shift to component mapping
786 _shiftmapping[i->shift().to_ulong()] = k;
787 }
788 _indexOffset.resize(_itends.size());
789 }
790
791 private:
792
793 friend class YGrid<Coordinates>::Iterator;
794 DAI _begin;
795 DAI _end;
796 std::array<int,StaticPower<2,dim>::power> _shiftmapping;
797 std::vector<typename YGridComponent<Coordinates>::Iterator> _itbegins;
798 std::vector<typename YGridComponent<Coordinates>::Iterator> _itends;
799 std::vector<int> _indexOffset;
800 };
801
803 template <class Coordinates>
804 inline std::ostream& operator<< (std::ostream& s, const YGrid<Coordinates>& e)
805 {
806 s << "Printing YGrid structure:" << std::endl;
807 for (auto it = e.dataBegin(); it != e.dataEnd(); ++it)
808 s << *it << std::endl;
809 return s;
810 }
811
819 template<class Coordinates>
820 class YGridList
821 {
822 public:
823 static const int dim = Coordinates::dimension;
824
826 struct Intersection
827 {
829 YGridComponent<Coordinates> grid;
831 int rank;
833 int distance;
835 YGrid<Coordinates> yg;
836 };
837
838 // define data array iterator type
839 typedef typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator DAI;
840
841 // iterator that allows to iterate over a concatenation of deques. namely those
842 // that belong to the same codimension.
843 class Iterator
844 {
845 public:
846
848 Iterator(const YGridList<Coordinates>& ygl, bool end=false) : _end(ygl.dataEnd()), _which(ygl.dataBegin())
849 {
850 _it = _which->begin();
851
852 // advance the iterator to the first element that exists.
853 // some deques might be empty and should be skipped
854 while ((_which != _end) && (_it == _which->end()))
855 {
856 ++_which;
857 if (_which != _end)
858 _it = _which->begin();
859 }
860 // the iterator is at the end if and only if _which==_end
861 if (end)
862 {
863 _which = _end;
864 }
865 }
866
868 Iterator& operator++ ()
869 {
870 ++_it;
871 // advance the iterator to the next element that exists.
872 // some deques might be empty and should be skipped
873 while ((_which != _end) && (_it == _which->end()))
874 {
875 ++_which;
876 if (_which != _end)
877 _it = _which->begin();
878 }
879 return *this;
880 }
881
883 typename std::deque<Intersection>::iterator operator->() const
884 {
885 return _it;
886 }
887
889 typename std::deque<Intersection>::iterator operator*() const
890 {
891 return _it;
892 }
893
895 bool operator== (const Iterator& i) const
896 {
897 if (_which != i._which)
898 return false;
899 if (_which == _end)
900 return true;
901 return _it == i._it;
902 }
903
905 bool operator!= (const Iterator& i) const
906 {
907 if (_which != i._which)
908 return true;
909 if (_which == _end)
910 return false;
911 return _it != i._it;
912 }
913
914 private:
915 typename std::deque<Intersection>::iterator _it;
916 DAI _end;
917 DAI _which;
918 };
919
921 Iterator begin() const
922 {
923 return Iterator(*this);
924 }
925
927 Iterator end() const
928 {
929 return Iterator(*this,true);
930 }
931
933 void setBegin(typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator begin)
934 {
935 _begin = begin;
936 }
937
939 DAI dataBegin() const
940 {
941 return _begin;
942 }
943
945 DAI dataEnd() const
946 {
947 return _end;
948 }
949
951 int size() const
952 {
953 int count = 0;
954 for (DAI it = _begin; it != _end; ++it)
955 count += it->size();
956 return count;
957 }
958
960 void finalize(DAI end, const YGrid<Coordinates>& ygrid)
961 {
962 // Instead of directly iterating over the intersection deques, this code
963 // iterates over the components of an associated ygrid and works its way
964 // through the list of intersection deques in parallel.
965 // The reason for this convoluted iteration technique is that there are not
966 // necessarily intersections for all possible shifts, but we have to make
967 // sure that we stop at each shift to update the per-component index shift.
968 // This is ensured by iterating over the ygrid, which is guaranteed to have
969 // a component for each shift vector.
970
971 // set end iterator in the data array
972 _end = end;
973
975 int offset = 0;
976
977 DAI i = _begin;
978
979 // make sure that we have a valid deque (i.e. a non-empty one)
980 while (i != _end && i->begin() == i->end())
981 ++i;
982
983 for (auto yit = ygrid.dataBegin(); yit != ygrid.dataEnd(); ++yit)
984 {
985 if (i == _end)
986 break;
987 auto it = i->begin();
988 if (it->grid.shift() == yit->shift())
989 {
990 // iterate over the intersections in the deque and set the offset
991 for (; it != i->end(); ++it)
992 {
993 it->yg.setBegin(&(it->grid));
994 it->yg.finalize(&(it->grid)+1, offset);
995 }
996
997 // advance to next non-empty deque
998 ++i;
999 while (i != _end && i->begin() == i->end())
1000 ++i;
1001 }
1002
1003 // update the offset from the ygrid component
1004 int add = 1;
1005 for (int j=0; j<dim; j++)
1006 add *= yit->supersize(j);
1007 offset += add;
1008 }
1009 assert (i == end);
1010 }
1011
1012 private:
1013 DAI _begin;
1014 DAI _end;
1015 };
1016
1017} // namespace Dune
1018
1019#endif
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Definition: ygrid.hh:295
Iterator & operator++()
Increment iterator to next cell with position.
Definition: ygrid.hh:376
int superindex() const
Return consecutive index in enclosing grid.
Definition: ygrid.hh:341
bool operator!=(const Iterator &i) const
Return true when two iterators over the same grid are not equal (!).
Definition: ygrid.hh:335
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:400
ct meshsize(int i) const
Return meshsize in direction i.
Definition: ygrid.hh:433
void reinit(const YGridComponent< Coordinates > &r, const iTupel &coord)
reinitialize iterator to given position
Definition: ygrid.hh:314
void move(const iTupel &dist)
move this iterator dist cells in direction i
Definition: ygrid.hh:366
fTupel meshsize() const
Return meshsize of current cell as reference.
Definition: ygrid.hh:439
fTupel lowerleft() const
Return lower left corner of the entity associated with the current coordinates and shift.
Definition: ygrid.hh:406
iTupel _coord
current position in index set
Definition: ygrid.hh:463
const iTupel & coord() const
Return coordinate of the cell as reference (do not modify).
Definition: ygrid.hh:353
void move(int i, int dist)
move this iterator dist cells in direction i
Definition: ygrid.hh:359
Iterator(const YGridComponent< Coordinates > &r)
Make iterator pointing to first cell in a grid.
Definition: ygrid.hh:301
ct upperright(int i) const
Return ith component of upper right corder of the entity associated with the current coordinates and ...
Definition: ygrid.hh:415
int coord(int i) const
Return coordinate of the cell in direction i.
Definition: ygrid.hh:347
fTupel upperright() const
Return upper right corder of the entity associated with the current coordinates and shift.
Definition: ygrid.hh:424
bool operator==(const Iterator &i) const
Return true when two iterators over the same grid are equal (!).
Definition: ygrid.hh:329
Iterator(const YGridComponent< Coordinates > &r, const iTupel &coord)
Make iterator pointing to given cell in a grid.
Definition: ygrid.hh:308
int _superindex
consecutive index in enclosing grid
Definition: ygrid.hh:464
Definition: ygrid.hh:72
int index(const iTupel &coord) const
given a tupel compute its index in the lexicographic numbering
Definition: ygrid.hh:249
int offset(int i) const
Return offset to origin of enclosing grid.
Definition: ygrid.hh:170
YGridComponent< Coordinates > move(iTupel v) const
return grid moved by the vector v
Definition: ygrid.hh:260
iTupel size() const
retrun size
Definition: ygrid.hh:200
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:268
int min(int i) const
Return minimum index in direction i.
Definition: ygrid.hh:215
YGridComponent(iTupel origin, iTupel size)
make ygrid without coordinate information
Definition: ygrid.hh:96
int totalsize() const
Return total size of index set which is the product of all size per direction.
Definition: ygrid.hh:206
YGridComponent(iTupel origin, iTupel size, const YGridComponent< Coordinates > &enclosing)
make a subgrid by taking coordinates from a larger grid
Definition: ygrid.hh:105
const iTupel & supersize() const
return size of enclosing grid
Definition: ygrid.hh:188
YGridComponent(iTupel origin, std::bitset< d > shift, Coordinates *coords, iTupel size, iTupel offset, iTupel supersize)
Make YGridComponent by giving all parameters.
Definition: ygrid.hh:128
Iterator begin() const
return iterator to first element of index set
Definition: ygrid.hh:484
Iterator begin(const iTupel &co) const
return iterator to given element of index set
Definition: ygrid.hh:490
YGridComponent()
make uninitialized ygrid
Definition: ygrid.hh:82
int origin(int i) const
Return origin in direction i.
Definition: ygrid.hh:141
int supersize(int i) const
return size of enclosing grid
Definition: ygrid.hh:182
int size(int i) const
return size in direction i
Definition: ygrid.hh:194
int max(int i) const
Return maximum index in direction i.
Definition: ygrid.hh:221
const std::bitset< d > & shift() const
Return shift tupel.
Definition: ygrid.hh:159
bool shift(int i) const
Return shift in direction i.
Definition: ygrid.hh:153
Iterator end() const
return subiterator to last element of index set
Definition: ygrid.hh:496
bool inside(const iTupel &coord) const
given a coordinate, return true if it is in the grid
Definition: ygrid.hh:238
bool empty() const
Return true if YGrid is empty, i.e. has size 0 in all directions.
Definition: ygrid.hh:227
const iTupel & offset() const
Return offset to origin of enclosing grid.
Definition: ygrid.hh:176
const iTupel & origin() const
return reference to origin
Definition: ygrid.hh:147
Iterator over a collection o YGrids A YGrid::Iterator is the heart of an entity in YaspGrid.
Definition: ygrid.hh:591
int which() const
return the current component number
Definition: ygrid.hh:716
void move(int i, int dist)
move the grid, this is only done and needed for codim 0
Definition: ygrid.hh:722
bool operator==(const Iterator &i) const
compare two iterators: component has to match
Definition: ygrid.hh:700
Dune::FieldVector< typename Coordinates::ctype, dim > meshsize() const
return the current meshsize vector
Definition: ygrid.hh:667
Iterator(const YGrid< Coordinates > &yg, const std::array< int, dim > &coords, int which=0)
construct an iterator from coordinates and component
Definition: ygrid.hh:599
std::bitset< dim > shift() const
return the shift vector
Definition: ygrid.hh:679
bool shift(int i) const
return the shift in direction i
Definition: ygrid.hh:673
bool operator!=(const Iterator &i) const
compare two iterators: component has to match
Definition: ygrid.hh:708
Iterator()
default constructor
Definition: ygrid.hh:595
const std::array< int, dim > & coord() const
return coordinate array at the current position
Definition: ygrid.hh:635
Iterator & operator++()
increment to the next entity jumping to next component if necessary
Definition: ygrid.hh:692
Coordinates::ctype meshsize(int i) const
return the current meshsize in direction i
Definition: ygrid.hh:661
int coord(int i) const
return coordinate at the current position (direction i)
Definition: ygrid.hh:629
Iterator(const YGrid< Coordinates > &yg, bool end=false)
create an iterator to start or end of the codimension
Definition: ygrid.hh:606
int superindex() const
return the superindex
Definition: ygrid.hh:685
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:621
implements a collection of YGridComponents which form a codimension Entities of given codimension c n...
Definition: ygrid.hh:548
int shiftmapping(const std::bitset< dim > &shift) const
get which component belongs to a given shift vector
Definition: ygrid.hh:564
DAI dataBegin() const
get start iterator in the data array
Definition: ygrid.hh:570
Iterator begin() const
return begin iterator for the codimension and partition the ygrid represents
Definition: ygrid.hh:745
DAI dataEnd() const
get end iterator in the data array
Definition: ygrid.hh:576
void setBegin(DAI begin)
set start iterator in the data array
Definition: ygrid.hh:558
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:582
Iterator begin(const std::array< int, dim > &coord, int which=0) const
return iterator pointint to a specified position
Definition: ygrid.hh:751
Iterator end() const
return end iterator for the codimension and partition the ygrid represents
Definition: ygrid.hh:757
Implements a vector constructed from a given type representing a field and a compile-time given size.
Dune namespace.
Definition: alignment.hh:11
Various implementations of the power function for run-time and static arguments.
std::array< int, d > sizeArray(const std::array< std::vector< ct >, d > &v)
Definition: ygrid.hh:26
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)