Dune Core Modules (2.5.2)
ygrid.hh
Go to the documentation of this file.
106 : _origin(origin), _shift(enclosing.shift()), _coords(enclosing.getCoords()), _size(size), _supersize(enclosing.supersize())
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)
532 inline std::ostream& operator<< (std::ostream& s, typename YGridComponent<Coordinates>::Iterator& e)
534 s << "Printing YGridComponent Iterator:" << std::endl << "Iterator at " << e.coord() << " (index ";
804 inline std::ostream& operator<< (std::ostream& s, const YGrid<Coordinates>& e)
820 class YGridList
826 struct Intersection
829 YGridComponent<Coordinates> grid;
831 int rank;
833 int distance;
835 YGrid<Coordinates> yg;
848 Iterator(const YGridList<Coordinates>& ygl, bool end=false) : _end(ygl.dataEnd()), _which(ygl.dataBegin())
921 Iterator begin() const
927 Iterator end() const
933 void setBegin(typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator begin)
939 DAI dataBegin() const
945 DAI dataEnd() const
951 int size() const
960 void finalize(DAI end, const YGrid<Coordinates>& ygrid)
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
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
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 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
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
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.
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
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 13, 23:29, 2024)