5#ifndef DUNE_GRID_YASPGRID_HH
6#define DUNE_GRID_YASPGRID_HH
19#include <dune/common/hybridutilities.hh>
20#include <dune/common/bigunsignedint.hh>
21#include <dune/common/math.hh>
22#include <dune/common/typetraits.hh>
23#include <dune/common/reservedvector.hh>
24#include <dune/common/parallel/communication.hh>
25#include <dune/common/parallel/mpihelper.hh>
26#include <dune/geometry/axisalignedcubegeometry.hh>
27#include <dune/geometry/type.hh>
33#include <dune/common/parallel/mpicommunication.hh>
54 template<
int dim,
class Coordinates>
class YaspGrid;
55 template<
int mydim,
int cdim,
class Gr
idImp>
class YaspGeometry;
56 template<
int codim,
int dim,
class Gr
idImp>
class YaspEntity;
58 template<
int codim, PartitionIteratorType pitype,
class Gr
idImp>
class YaspLevelIterator;
62 template<
class Gr
idImp,
bool isLeafIndexSet>
class YaspIndexSet;
90 template<
int dim,
class Coordinates>
109 bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim>,
111 bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim>,
117 std::array<GeometryType,1>>
122 template<
int dim,
int codim>
123 struct YaspCommunicateMeta {
124 template<
class G,
class DataHandle>
127 if (data.contains(dim,codim))
129 g.template communicateCodim<DataHandle,codim>(data,iftype,dir,level);
131 YaspCommunicateMeta<dim,codim-1>::comm(g,data,iftype,dir,level);
136 struct YaspCommunicateMeta<dim,0> {
137 template<
class G,
class DataHandle>
138 static void comm (
const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir,
int level)
140 if (data.contains(dim,0))
141 g.template communicateCodim<DataHandle,0>(data,iftype,dir,level);
163 template<
int dim,
class Coordinates = Equ
idistantCoordinates<
double, dim> >
168 template<
int, PartitionIteratorType,
typename>
179 typedef typename Coordinates::ctype
ctype;
201 std::array<YGrid, dim+1> overlap;
202 std::array<YGridComponent<Coordinates>, Dune::power(2,dim)>
overlap_data;
205 std::array<YGrid, dim+1> interior;
206 std::array<YGridComponent<Coordinates>, Dune::power(2,dim)>
interior_data;
238 typedef std::array<int, dim> iTupel;
239 typedef FieldVector<ctype, dim> fTupel;
266 return _coarseSize[i] * (1 <<
l);
273 for (
int i=0; i<dim; ++i)
295 return _levels.begin();
303 return std::next(_levels.begin(), i);
309 return _levels.end();
330 g.overlapSize = overlap;
334 g.keepOverlap = keep_ovlp;
339 typename std::array<YGridComponent<Coordinates>,
power(2,dim)>::iterator
overlapfront_it = g.overlapfront_data.begin();
340 typename std::array<YGridComponent<Coordinates>,
power(2,dim)>::iterator
overlap_it = g.overlap_data.begin();
341 typename std::array<YGridComponent<Coordinates>,
power(2,dim)>::iterator
interiorborder_it = g.interiorborder_data.begin();
342 typename std::array<YGridComponent<Coordinates>,
power(2,dim)>::iterator
interior_it = g.interior_data.begin();
344 typename std::array<std::deque<Intersection>,
power(2,dim)>::iterator
346 typename std::array<std::deque<Intersection>,
power(2,dim)>::iterator
349 typename std::array<std::deque<Intersection>,
power(2,dim)>::iterator
351 typename std::array<std::deque<Intersection>,
power(2,dim)>::iterator
354 typename std::array<std::deque<Intersection>,
power(2,dim)>::iterator
356 typename std::array<std::deque<Intersection>,
power(2,dim)>::iterator
359 typename std::array<std::deque<Intersection>,
power(2,dim)>::iterator
361 typename std::array<std::deque<Intersection>,
power(2,dim)>::iterator
365 std::array<int,dim> n;
366 std::fill(n.begin(), n.end(), 0);
376 for (
int i=0; i<dim; i++)
406 for (
unsigned int codim = 0; codim < dim + 1; codim++)
423 for (
unsigned int index = 0; index < (1<<dim); index++)
426 std::bitset<dim>
r(index);
427 if (
r.count() != dim-codim)
436 for (
int i=0; i<dim; i++)
442 for (
int i=0; i<dim; i++)
458 for (
int i=0; i<dim; i++)
462 origin[i] += overlap;
480 for (
int i=0; i<dim; i++)
540 struct mpifriendly_ygrid {
543 std::fill(origin.begin(), origin.end(), 0);
544 std::fill(
size.begin(),
size.end(), 0);
546 mpifriendly_ygrid (
const YGridComponent<Coordinates>& grid)
547 : origin(grid.origin()),
size(grid.
size())
585 iTupel coord = _torus.
coord();
586 iTupel delta = i.delta();
588 for (
int k=0;
k<dim;
k++)
nb[
k] += delta[
k];
590 std::fill(v.begin(), v.end(), 0);
592 for (
int k=0;
k<dim;
k++)
687 std::array<int, dim>
sides;
689 for (
int i=0; i<dim; i++)
692 ((
begin()->overlap[0].dataBegin()->origin(i) == 0)+
693 (
begin()->overlap[0].dataBegin()->origin(i) +
begin()->overlap[0].dataBegin()->size(i)
698 for (
int k=0;
k<dim;
k++)
701 for (
int l=0;
l<dim;
l++)
704 offset *=
begin()->overlap[0].dataBegin()->size(
l);
706 nBSegments +=
sides[
k]*offset;
733 std::bitset<dim>
periodic = std::bitset<dim>(0
ULL),
738 , leafIndexSet_(*
this)
748 for (std::size_t i=0; i<dim; i++)
755 std::fill(
o.begin(),
o.end(), 0);
765 for (std::size_t i=0; i<dim; i++)
771 for (
int i=0; i<dim; i++)
777 int mysteryFactor = (std::is_same<Coordinates,EquidistantCoordinates<ctype,dim> >::value) ? 1 : 2;
780 for (
int i=0; i<dim; i++)
789 DUNE_THROW(
Dune::GridError,
"YaspGrid does not support degrees of freedom shared by more than immediately neighboring subdomains."
790 " Note that this also holds for DOFs on subdomain boundaries."
791 " Increase grid elements or decrease overlap accordingly.");
799 for (
int i=0; i<dim; i++)
808 for (
int i=0; i<dim; i++)
811 if constexpr (std::is_same_v<Coordinates,EquidistantCoordinates<ctype,dim>>)
820 if constexpr (std::is_same_v<Coordinates,EquidistantOffsetCoordinates<ctype,dim>>)
822 Dune::FieldVector<ctype,dim> lowerleft;
823 for (
int i=0; i<dim; i++)
834 if constexpr (std::is_same_v<Coordinates,TensorProductCoordinates<ctype,dim>>)
836 std::array<std::vector<ctype>,dim>
newCoords;
840 for (
int i=0; i<dim; ++i)
851 offset[i] -= overlap;
859 for (std::size_t j=
begin; j<
end; j++)
870 for (
int j=0; j<overlap; ++j)
876 offset[i] -= overlap;
903 template<
class C = Coordinates,
904 typename std::enable_if_t< std::is_same_v<C, EquidistantCoordinates<ctype,dim> >,
int> = 0>
906 std::array<
int, std::size_t{dim}> s,
907 std::bitset<std::size_t{dim}> periodic = std::bitset<std::size_t{dim}>{0ULL},
911 : ccobj(
comm), _torus(
comm,tag,s,
overlap,partitioner), leafIndexSet_(*this),
912 _L(L), _periodic(periodic), _coarseSize(s), _overlap(
overlap),
913 keep_ovlp(true), adaptRefCount(0), adaptActive(false)
918 std::fill(o.begin(), o.end(), 0);
919 iTupel o_interior(o);
920 iTupel s_interior(s);
926 for (
int i=0; i<dim; i++)
929 int toosmall = (s_interior[i] < 2*
overlap) &&
930 (periodic[i] || (s_interior[i] != s[i]));
933 MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR,
comm);
935 DUNE_THROW(
Dune::GridError,
"YaspGrid does not support degrees of freedom shared by more than immediately neighboring subdomains."
936 " Note that this also holds for DOFs on subdomain boundaries."
937 " Increase grid elements or decrease overlap accordingly.");
941 iTupel s_overlap(s_interior);
942 for (
int i=0; i<dim; i++)
944 if ((o_interior[i] - overlap > 0) || (periodic[i]))
946 if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
950 FieldVector<ctype,dim> upperRightWithOverlap;
952 for (
int i=0; i<dim; i++)
953 upperRightWithOverlap[i] = (L[i] / s[i]) * s_overlap[i];
956 EquidistantCoordinates<ctype,dim> cc(upperRightWithOverlap,s_overlap);
959 makelevel(cc,periodic,o_interior,overlap);
973 template<
class C = Coordinates,
974 typename std::enable_if_t< std::is_same_v<C, EquidistantOffsetCoordinates<ctype,dim> >,
int> = 0>
976 Dune::FieldVector<ctype, dim> upperright,
977 std::array<
int, std::size_t{dim}> s,
978 std::bitset<std::size_t{dim}> periodic = std::bitset<std::size_t{dim}>(0ULL),
982 : ccobj(
comm), _torus(
comm,tag,s,
overlap,partitioner), leafIndexSet_(*this),
983 _L(upperright - lowerleft),
984 _periodic(periodic), _coarseSize(s), _overlap(
overlap),
985 keep_ovlp(true), adaptRefCount(0), adaptActive(false)
990 std::fill(o.begin(), o.end(), 0);
991 iTupel o_interior(o);
992 iTupel s_interior(s);
998 for (
int i=0; i<dim; i++)
1001 int toosmall = (s_interior[i] < 2*
overlap) &&
1002 (periodic[i] || (s_interior[i] != s[i]));
1005 MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR,
comm);
1007 DUNE_THROW(
Dune::GridError,
"YaspGrid does not support degrees of freedom shared by more than immediately neighboring subdomains."
1008 " Note that this also holds for DOFs on subdomain boundaries."
1009 " Increase grid elements or decrease overlap accordingly.");
1013 iTupel s_overlap(s_interior);
1014 for (
int i=0; i<dim; i++)
1016 if ((o_interior[i] - overlap > 0) || (periodic[i]))
1018 if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
1022 FieldVector<ctype,dim> upperRightWithOverlap;
1023 for (
int i=0; i<dim; i++)
1024 upperRightWithOverlap[i] = lowerleft[i]
1025 + s_overlap[i] * (upperright[i]-lowerleft[i]) / s[i];
1027 EquidistantOffsetCoordinates<ctype,dim> cc(lowerleft,upperRightWithOverlap,s_overlap);
1030 makelevel(cc,periodic,o_interior,overlap);
1042 template<
class C = Coordinates,
1043 typename std::enable_if_t< std::is_same_v<C, TensorProductCoordinates<ctype,dim> >,
int> = 0>
1044 YaspGrid (std::array<std::vector<ctype>, std::size_t{dim}> coords,
1045 std::bitset<std::size_t{dim}> periodic = std::bitset<std::size_t{dim}>(0ULL),
1050 leafIndexSet_(*this), _periodic(periodic), _overlap(
overlap),
1051 keep_ovlp(true), adaptRefCount(0), adaptActive(false)
1053 if (!Dune::Yasp::checkIfMonotonous(coords))
1054 DUNE_THROW(
Dune::GridError,
"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
1059 for (
int i=0; i<dim; i++) {
1060 _coarseSize[i] = coords[i].size() - 1;
1061 _L[i] = coords[i][_coarseSize[i]] - coords[i][0];
1065 std::fill(o.begin(), o.end(), 0);
1066 iTupel o_interior(o);
1067 iTupel s_interior(_coarseSize);
1069 _torus.
partition(_torus.
rank(),o,_coarseSize,o_interior,s_interior);
1073 for (
int i=0; i<dim; i++)
1076 int toosmall = (s_interior[i] < 2*
overlap) &&
1077 (periodic[i] || (s_interior[i] != _coarseSize[i]));
1080 MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR,
comm);
1082 DUNE_THROW(
Dune::GridError,
"YaspGrid does not support degrees of freedom shared by more than immediately neighboring subdomains."
1083 " Note that this also holds for DOFs on subdomain boundaries."
1084 " Increase grid elements or decrease overlap accordingly.");
1089 std::array<std::vector<ctype>,dim> newcoords;
1090 std::array<int, dim> offset(o_interior);
1093 for (
int i=0; i<dim; ++i)
1096 typename std::vector<ctype>::iterator
begin = coords[i].begin() + o_interior[i];
1097 typename std::vector<ctype>::iterator
end =
begin + s_interior[i] + 1;
1101 if (o_interior[i] - overlap > 0)
1106 if (o_interior[i] + s_interior[i] + overlap < _coarseSize[i])
1115 if ((periodic[i]) && (o_interior[i] + s_interior[i] + overlap >= _coarseSize[i]))
1118 typename std::vector<ctype>::iterator it = coords[i].begin();
1120 newcoords[i].push_back(newcoords[i].back() - *it + *(++it));
1123 if ((periodic[i]) && (o_interior[i] - overlap <= 0))
1128 typename std::vector<ctype>::iterator it = coords[i].end() - 1;
1130 newcoords[i].insert(newcoords[i].
begin(), newcoords[i].
front() - *it + *(--it));
1134 TensorProductCoordinates<ctype,dim> cc(newcoords, offset);
1137 makelevel(cc,periodic,o_interior,overlap);
1157 YaspGrid (std::array<std::vector<ctype>, std::size_t{dim}> coords,
1158 std::bitset<std::size_t{dim}> periodic,
1161 std::array<int,dim> coarseSize,
1163 : ccobj(
comm), _torus(
comm,tag,coarseSize,overlap,partitioner), leafIndexSet_(*this),
1164 _periodic(periodic), _coarseSize(coarseSize), _overlap(overlap),
1165 keep_ovlp(true), adaptRefCount(0), adaptActive(false)
1168 static_assert(std::is_same<Coordinates,TensorProductCoordinates<ctype,dim> >::value,
1169 "YaspGrid coordinate container template parameter and given constructor values do not match!");
1171 if (!Dune::Yasp::checkIfMonotonous(coords))
1172 DUNE_THROW(
Dune::GridError,
"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
1174 for (
int i=0; i<dim; i++)
1175 _L[i] = coords[i][coords[i].
size() - 1] - coords[i][0];
1179 std::array<int,dim> o;
1180 std::fill(o.begin(), o.end(), 0);
1181 std::array<int,dim> o_interior(o);
1182 std::array<int,dim> s_interior(coarseSize);
1184 _torus.
partition(_torus.
rank(),o,coarseSize,o_interior,s_interior);
1187 std::array<int,dim> offset(o_interior);
1188 for (
int i=0; i<dim; i++)
1189 if ((periodic[i]) || (o_interior[i] > 0))
1190 offset[i] -= overlap;
1192 TensorProductCoordinates<ctype,dim> cc(coords, offset);
1195 makelevel(cc,periodic,o_interior,overlap);
1201 friend struct BackupRestoreFacility<
YaspGrid<dim,Coordinates> >;
1213 return _levels.size()-1;
1221 "Coarsening " << -refCount <<
" levels requested!");
1224 for (
int k=refCount;
k<0;
k++)
1228 _levels.back() = empty;
1232 indexsets.pop_back();
1236 for (
int k=0;
k<refCount;
k++)
1242 for (
int i=0; i<dim; i++)
1244 if (
cg.overlap[0].dataBegin()->origin(i) > 0 || _periodic[i])
1246 if (
cg.overlap[0].dataBegin()->max(i) + 1 <
globalSize(i) || _periodic[i])
1252 int overlap = (keep_ovlp) ? 2*
cg.overlapSize :
cg.overlapSize;
1256 for (
int i=0; i<dim; i++)
1257 o_interior[i] = 2*
cg.interior[0].dataBegin()->origin(i);
1260 _levels.resize(_levels.size() + 1);
1289 assert(adaptActive ==
false);
1290 if (
e.level() !=
maxLevel())
return false;
1291 adaptRefCount = std::max(adaptRefCount, refCount);
1303 return (
e.level() ==
maxLevel() ) ? adaptRefCount : 0;
1310 return (adaptRefCount > 0);
1317 adaptRefCount =
comm().max(adaptRefCount);
1318 return (adaptRefCount < 0);
1324 adaptActive =
false;
1329 template<
int cd, PartitionIteratorType pitype>
1336 template<
int cd, PartitionIteratorType pitype>
1357 template<
int cd, PartitionIteratorType pitype>
1364 template<
int cd, PartitionIteratorType pitype>
1385 template <
typename Seed>
1386 typename Traits::template Codim<Seed::codimension>::Entity
1389 const int codim = Seed::codimension;
1396 return Entity(EntityImp(g,
YIterator(g->overlapfront[codim],seed.impl().coord(),seed.impl().offset())));
1403 return g->overlapSize;
1410 return g->overlapSize;
1426 int size (
int level,
int codim)
const
1432 typedef typename std::array<YGridComponent<Coordinates>, Dune::power(2,dim)>::iterator DAI;
1433 for (DAI
it = g->overlapfront[codim].dataBegin();
it != g->overlapfront[codim].dataEnd(); ++
it)
1434 count +=
it->totalsize();
1446 int size (
int level, GeometryType type)
const
1448 return (type.isCube()) ?
size(level,dim-type.dim()) : 0;
1472 template<
class DataHandleImp,
class DataType>
1482 template<
class DataHandleImp,
class DataType>
1492 template<
class DataHandle,
int codim>
1496 if (!data.contains(dim,codim))
return;
1499 typedef typename DataHandle::DataType DataType;
1510 sendlist = &g->send_interiorborder_interiorborder[codim];
1511 recvlist = &g->recv_interiorborder_interiorborder[codim];
1515 sendlist = &g->send_interiorborder_overlapfront[codim];
1516 recvlist = &g->recv_overlapfront_interiorborder[codim];
1520 sendlist = &g->send_overlap_overlapfront[codim];
1521 recvlist = &g->recv_overlapfront_overlap[codim];
1525 sendlist = &g->send_overlapfront_overlapfront[codim];
1526 recvlist = &g->recv_overlapfront_overlapfront[codim];
1544 if (data.fixedSize(dim,codim))
1571 size_t *buf =
new size_t[is->grid.totalsize()];
1575 int i=0;
size_t n=0;
1582 buf[i] = data.size(*
it);
1591 torus().send(is->rank,buf,is->grid.totalsize()*
sizeof(
size_t));
1600 size_t *buf =
new size_t[is->grid.totalsize()];
1604 torus().recv(is->rank,buf,is->grid.totalsize()*
sizeof(
size_t));
1629 for (
int i=0; i<is->grid.totalsize(); ++i)
1640 std::vector<DataType*>
sends(
sendlist->size(),
static_cast<DataType*
>(0));
1659 data.gather(
mb,*
it);
1667 std::vector<DataType*>
recvs(
recvlist->size(),
static_cast<DataType*
>(0));
1705 if (data.fixedSize(dim,codim))
1709 size_t n=data.size(*
it);
1713 data.scatter(
mb,*
it,n);
1737 return theglobalidset;
1742 return theglobalidset;
1748 return *(indexsets[level]);
1753 return leafIndexSet_;
1778 template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
1782 class MessageBuffer {
1785 MessageBuffer (DT *p)
1794 void write (
const Y& data)
1796 static_assert(( std::is_same<DT,Y>::value ),
"DataType mismatch");
1802 void read (
Y& data)
const
1804 static_assert(( std::is_same<DT,Y>::value ),
"DataType mismatch");
1815 template<
int cd, PartitionIteratorType pitype>
1819 if (level<0 || level>
maxLevel()) DUNE_THROW(RangeError,
"level out of range");
1830 return levelend <cd, pitype> (level);
1832 DUNE_THROW(
GridError,
"YaspLevelIterator with this codim or partition type not implemented");
1836 template<
int cd, PartitionIteratorType pitype>
1837 YaspLevelIterator<cd,pitype,GridImp> levelend (
int level)
const
1840 if (level<0 || level>
maxLevel()) DUNE_THROW(RangeError,
"level out of range");
1843 return YaspLevelIterator<cd,pitype,GridImp>(g,g->interior[cd].end());
1845 return YaspLevelIterator<cd,pitype,GridImp>(g,g->interiorborder[cd].end());
1847 return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlap[cd].end());
1849 return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlapfront[cd].end());
1851 DUNE_THROW(GridError,
"YaspLevelIterator with this codim or partition type not implemented");
1856 Torus<Communication,dim> _torus;
1858 std::vector< std::shared_ptr< YaspIndexSet<const YaspGrid<dim,Coordinates>,
false > > > indexsets;
1859 YaspIndexSet<const YaspGrid<dim,Coordinates>,
true> leafIndexSet_;
1860 YaspGlobalIdSet<const YaspGrid<dim,Coordinates> > theglobalidset;
1862 Dune::FieldVector<ctype, dim> _L;
1864 std::bitset<dim> _periodic;
1866 ReservedVector<YGridLevel,32> _levels;
1874 template<
typename ctype,
int dim>
1876 std::array<
int, std::size_t{dim}>,
1877 std::bitset<std::size_t{dim}> = std::bitset<std::size_t{dim}>{0ULL},
1880 const Yasp::Partitioning<dim>* = YaspGrid<dim, EquidistantCoordinates<ctype, dim>>::defaultPartitioner())
1881 -> YaspGrid< dim, EquidistantCoordinates<ctype, dim> >;
1883 template<
typename ctype,
int dim>
1885 FieldVector<ctype, dim>,
1886 std::array<
int, std::size_t{dim}>,
1887 std::bitset<std::size_t{dim}> = std::bitset<std::size_t{dim}>{0ULL},
1890 const Yasp::Partitioning<dim>* = YaspGrid<dim, EquidistantOffsetCoordinates<ctype, dim>>::defaultPartitioner())
1891 -> YaspGrid< dim, EquidistantOffsetCoordinates<ctype, dim> >;
1893 template<
typename ctype, std::
size_t dim>
1895 std::bitset<dim> = std::bitset<dim>{0ULL},
1898 const Yasp::Partitioning<dim>* = YaspGrid<
int{dim}, TensorProductCoordinates<ctype,
int{dim}>>::defaultPartitioner())
1899 -> YaspGrid<
int{dim}, TensorProductCoordinates<ctype,
int{dim}> >;
1902 template <
int d,
class CC>
1905 int rank = grid.
torus().rank();
1907 s <<
"[" << rank <<
"]:" <<
" YaspGrid maxlevel=" << grid.
maxLevel() << std::endl;
1909 s <<
"Printing the torus: " <<std::endl;
1910 s << grid.
torus() << std::endl;
1914 s <<
"[" << rank <<
"]: " << std::endl;
1915 s <<
"[" << rank <<
"]: " <<
"==========================================" << std::endl;
1916 s <<
"[" << rank <<
"]: " <<
"level=" << g->level() << std::endl;
1918 for (
int codim = 0; codim < d + 1; ++codim)
1920 s <<
"[" << rank <<
"]: " <<
"overlapfront[" << codim <<
"]: " << g->overlapfront[codim] << std::endl;
1921 s <<
"[" << rank <<
"]: " <<
"overlap[" << codim <<
"]: " << g->overlap[codim] << std::endl;
1922 s <<
"[" << rank <<
"]: " <<
"interiorborder[" << codim <<
"]: " << g->interiorborder[codim] << std::endl;
1923 s <<
"[" << rank <<
"]: " <<
"interior[" << codim <<
"]: " << g->interior[codim] << std::endl;
1926 for (I i=g->send_overlapfront_overlapfront[codim].begin();
1927 i!=g->send_overlapfront_overlapfront[codim].end(); ++i)
1928 s <<
"[" << rank <<
"]: " <<
" s_of_of[" << codim <<
"] to rank "
1929 << i->rank <<
" " << i->grid << std::endl;
1931 for (I i=g->recv_overlapfront_overlapfront[codim].begin();
1932 i!=g->recv_overlapfront_overlapfront[codim].end(); ++i)
1933 s <<
"[" << rank <<
"]: " <<
" r_of_of[" << codim <<
"] to rank "
1934 << i->rank <<
" " << i->grid << std::endl;
1936 for (I i=g->send_overlap_overlapfront[codim].begin();
1937 i!=g->send_overlap_overlapfront[codim].end(); ++i)
1938 s <<
"[" << rank <<
"]: " <<
" s_o_of[" << codim <<
"] to rank "
1939 << i->rank <<
" " << i->grid << std::endl;
1941 for (I i=g->recv_overlapfront_overlap[codim].begin();
1942 i!=g->recv_overlapfront_overlap[codim].end(); ++i)
1943 s <<
"[" << rank <<
"]: " <<
" r_of_o[" << codim <<
"] to rank "
1944 << i->rank <<
" " << i->grid << std::endl;
1946 for (I i=g->send_interiorborder_interiorborder[codim].begin();
1947 i!=g->send_interiorborder_interiorborder[codim].end(); ++i)
1948 s <<
"[" << rank <<
"]: " <<
" s_ib_ib[" << codim <<
"] to rank "
1949 << i->rank <<
" " << i->grid << std::endl;
1951 for (I i=g->recv_interiorborder_interiorborder[codim].begin();
1952 i!=g->recv_interiorborder_interiorborder[codim].end(); ++i)
1953 s <<
"[" << rank <<
"]: " <<
" r_ib_ib[" << codim <<
"] to rank "
1954 << i->rank <<
" " << i->grid << std::endl;
1956 for (I i=g->send_interiorborder_overlapfront[codim].begin();
1957 i!=g->send_interiorborder_overlapfront[codim].end(); ++i)
1958 s <<
"[" << rank <<
"]: " <<
" s_ib_of[" << codim <<
"] to rank "
1959 << i->rank <<
" " << i->grid << std::endl;
1961 for (I i=g->recv_overlapfront_interiorborder[codim].begin();
1962 i!=g->recv_overlapfront_interiorborder[codim].end(); ++i)
1963 s <<
"[" << rank <<
"]: " <<
" r_of_ib[" << codim <<
"] to rank "
1964 << i->rank <<
" " << i->grid << std::endl;
1973 namespace Capabilities
1983 template<
int dim,
class Coordinates>
1984 struct hasBackupRestoreFacilities< YaspGrid<dim, Coordinates> >
1986 static const bool v =
true;
1992 template<
int dim,
class Coordinates>
1993 struct hasSingleGeometryType<
YaspGrid<dim, Coordinates> >
1995 static const bool v =
true;
1996 static const unsigned int topologyId = GeometryTypes::cube(dim).id();
2002 template<
int dim,
class Coordinates>
2003 struct isCartesian<
YaspGrid<dim, Coordinates> >
2005 static const bool v =
true;
2011 template<
int dim,
class Coordinates,
int codim>
2012 struct hasEntity<
YaspGrid<dim, Coordinates>, codim>
2014 static const bool v =
true;
2021 template<
int dim,
class Coordinates,
int codim>
2022 struct hasEntityIterator<
YaspGrid<dim, Coordinates>, codim>
2024 static const bool v =
true;
2030 template<
int dim,
int codim,
class Coordinates>
2031 struct canCommunicate<
YaspGrid< dim, Coordinates>, codim >
2033 static const bool v =
true;
2039 template<
int dim,
class Coordinates>
2040 struct isLevelwiseConforming<
YaspGrid<dim, Coordinates> >
2042 static const bool v =
true;
2048 template<
int dim,
class Coordinates>
2049 struct isLeafwiseConforming<
YaspGrid<dim, Coordinates> >
2051 static const bool v =
true;
2057 template<
int dim,
class Coordinates>
2058 struct viewThreadSafe<
YaspGrid<dim, Coordinates> >
2060 static const bool v =
true;
This provides a YGrid, the elemental component of the yaspgrid implementation.
Specialization of the PersistentContainer for YaspGrid.
The YaspLevelIterator class.
The YaspIntersectionIterator class.
The YaspIntersection class.
level-wise, non-persistent, consecutive indices for YaspGrid
The YaspGeometry class and its specializations.
The YaspEntitySeed class.
the YaspEntity class and its specializations
This file provides the infrastructure for toroidal communication in YaspGrid.
Specialization of the StructuredGridFactory class for YaspGrid.
This provides container classes for the coordinates to be used in YaspGrid Upon implementation of the...
Describes the parallel communication interface class for MessageBuffers and DataHandles.
std::ostream & operator<<(std::ostream &out, const PartitionType &type)
write a PartitionType to a stream
Definition gridenums.hh:72
CommunicationDirection
Define a type for communication direction parameter.
Definition gridenums.hh:170
InterfaceType
Parameter to be used for the communication functions.
Definition gridenums.hh:86
@ All_Partition
all entities
Definition gridenums.hh:141
@ Interior_Partition
only interior entities
Definition gridenums.hh:137
@ InteriorBorder_Partition
interior and border entities
Definition gridenums.hh:138
@ Overlap_Partition
interior, border, and overlap entities
Definition gridenums.hh:139
@ Ghost_Partition
only ghost entities
Definition gridenums.hh:142
@ BackwardCommunication
reverse communication direction
Definition gridenums.hh:172
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition gridenums.hh:88
@ All_All_Interface
send all and receive all entities
Definition gridenums.hh:91
@ Overlap_All_Interface
send overlap, receive all entities
Definition gridenums.hh:90
@ Overlap_OverlapFront_Interface
send overlap, receive overlap and front entities
Definition gridenums.hh:89
@ InteriorBorder_InteriorBorder_Interface
send/receive interior and border entities
Definition gridenums.hh:87
void swap(Dune::PersistentContainer< G, T > &a, Dune::PersistentContainer< G, T > &b)
Definition utility/persistentcontainer.hh:83
Include standard header files.
Definition agrid.hh:60
const int yaspgrid_level_bits
Definition yaspgrid.hh:48
Communication< MPI_Comm > YaspCommunication
Definition yaspgrid.hh:85
const int yaspgrid_dim_bits
Definition yaspgrid.hh:47
YaspGrid(FieldVector< ctype, dim >, std::array< int, std::size_t{dim}>, std::bitset< std::size_t{dim}>=std::bitset< std::size_t{dim}>{0ULL}, int=1, YaspCommunication=YaspCommunication(), const Yasp::Partitioning< dim > *=YaspGrid< dim, EquidistantCoordinates< ctype, dim > >::defaultPartitioner()) -> YaspGrid< dim, EquidistantCoordinates< ctype, dim > >
constexpr Overlap overlap
PartitionSet for the overlap partition.
Definition partitionset.hh:277
constexpr Front front
PartitionSet for the front partition.
Definition partitionset.hh:280
std::array< int, d > sizeArray(const std::array< std::vector< ct >, d > &v)
Definition ygrid.hh:29
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition common/intersection.hh:164
static const bool v
Definition common/capabilities.hh:28
static const unsigned int topologyId
Definition common/capabilities.hh:31
static const bool v
Definition common/capabilities.hh:50
static const bool v
Definition common/capabilities.hh:59
static const bool v
Definition common/capabilities.hh:75
static const bool v
Definition common/capabilities.hh:98
static const bool v
Definition common/capabilities.hh:107
static const bool v
Definition common/capabilities.hh:116
static const bool v
Definition common/capabilities.hh:125
static const bool v
Definition common/capabilities.hh:170
Definition defaultgridview.hh:26
Definition defaultgridview.hh:211
Wrapper class for entities.
Definition common/entity.hh:66
Base class for exceptions in Dune grid modules.
Definition exceptions.hh:20
Definition common/grid.hh:848
int size(int level, int codim) const
Return number of grid entities of a given codim on a given level in this process.
Definition common/grid.hh:538
A Traits struct that collects all associated types of one implementation.
Definition common/grid.hh:411
A traits struct that collects all associated types of one grid model.
Definition common/grid.hh:1013
Traits associated with a specific codim.
Definition common/grid.hh:1035
[ provides Dune::Grid ]
Definition yaspgrid.hh:166
YaspGridFamily< dim, Coordinates > GridFamily
the GridFamily of this grid
Definition yaspgrid.hh:716
YaspIndexSet< YaspGrid< dim, Coordinates >, false > LevelIndexSetType
Definition yaspgrid.hh:721
int ghostSize(int level, int codim) const
return size (= distance in graph) of ghost region
Definition yaspgrid.hh:1414
const Torus< Communication, dim > & torus() const
return reference to torus
Definition yaspgrid.hh:246
typename Base::Communication Communication
Definition yaspgrid.hh:181
const Traits::LeafIndexSet & leafIndexSet() const
Definition yaspgrid.hh:1751
void init()
Definition yaspgrid.hh:678
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lbegin(int level) const
version without second template parameter for convenience
Definition yaspgrid.hh:1344
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition yaspgrid.hh:1446
Traits::template Codim< Seed::codimension >::Entity entity(const Seed &seed) const
Definition yaspgrid.hh:1387
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity in maxLevel
Definition yaspgrid.hh:1372
const Traits::GlobalIdSet & globalIdSet() const
Definition yaspgrid.hh:1735
const Traits::LocalIdSet & localIdSet() const
Definition yaspgrid.hh:1740
YaspGridFamily< dim, Coordinates >::Traits Traits
Definition yaspgrid.hh:718
bool getRefineOption() const
Definition yaspgrid.hh:284
void globalRefine(int refCount)
refine the grid refCount times.
Definition yaspgrid.hh:1217
int getMark(const typename Traits::template Codim< 0 >::Entity &e) const
returns adaptation mark for given entity
Definition yaspgrid.hh:1301
void boundarysegmentssize()
Definition yaspgrid.hh:684
YaspGrid(std::array< std::vector< ctype >, std::size_t{dim}> coords, std::bitset< std::size_t{dim}> periodic=std::bitset< std::size_t{dim}>(0ULL), int overlap=1, Communication comm=Communication(), const Yasp::Partitioning< dim > *partitioner=defaultPartitioner())
Standard constructor for a tensorproduct YaspGrid.
Definition yaspgrid.hh:1044
YaspIndexSet< YaspGrid< dim, Coordinates >, true > LeafIndexSetType
Definition yaspgrid.hh:722
int ghostSize(int codim) const
return size (= distance in graph) of ghost region
Definition yaspgrid.hh:1420
int globalSize(int i) const
return number of cells on finest level in given direction on all processors
Definition yaspgrid.hh:252
YaspGrid(Dune::FieldVector< ctype, dim > lowerleft, Dune::FieldVector< ctype, dim > upperright, std::array< int, std::size_t{dim}> s, std::bitset< std::size_t{dim}> periodic=std::bitset< std::size_t{dim}>(0ULL), int overlap=1, Communication comm=Communication(), const Yasp::Partitioning< dim > *partitioner=defaultPartitioner())
Definition yaspgrid.hh:975
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafend() const
return LeafIterator which points behind the last entity in maxLevel
Definition yaspgrid.hh:1365
void postAdapt()
clean up some markers
Definition yaspgrid.hh:1322
const Dune::FieldVector< ctype, dim > & domainSize() const
returns the size of the physical domain
Definition yaspgrid.hh:1464
YGridLevelIterator end() const
return iterator pointing to one past the finest level
Definition yaspgrid.hh:307
friend class Entity
Definition yaspgrid.hh:1779
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition yaspgrid.hh:1452
int maxLevel() const
Definition yaspgrid.hh:1211
YaspGrid(Dune::FieldVector< ctype, dim > L, std::array< int, std::size_t{dim}> s, std::bitset< std::size_t{dim}> periodic=std::bitset< std::size_t{dim}>{0ULL}, int overlap=1, Communication comm=Communication(), const Yasp::Partitioning< dim > *partitioner=defaultPartitioner())
Definition yaspgrid.hh:905
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafend() const
return LeafIterator which points behind the last entity in maxLevel
Definition yaspgrid.hh:1379
void communicate(CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir, int level) const
Definition yaspgrid.hh:1473
void intersections(const YGridComponent< Coordinates > &sendgrid, const YGridComponent< Coordinates > &recvgrid, std::deque< Intersection > &sendlist, std::deque< Intersection > &recvlist)
Construct list of intersections with neighboring processors.
Definition yaspgrid.hh:562
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity in maxLevel
Definition yaspgrid.hh:1358
bool preAdapt()
returns true, if the grid will be coarsened
Definition yaspgrid.hh:1314
iTupel levelSize(int l) const
return size vector of the grid (in cells) on level l
Definition yaspgrid.hh:270
bool mark(int refCount, const typename Traits::template Codim< 0 >::Entity &e)
Marks an entity to be refined/coarsened in a subsequent adapt.
Definition yaspgrid.hh:1287
int overlapSize(int level, int codim) const
return size (= distance in graph) of overlap region
Definition yaspgrid.hh:1400
void communicate(CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir) const
Definition yaspgrid.hh:1483
int size(int codim) const
number of leaf entities per codim in this process
Definition yaspgrid.hh:1440
bool isPeriodic(int i) const
return whether the grid is periodic in direction i
Definition yaspgrid.hh:279
const Traits::LevelIndexSet & levelIndexSet(int level) const
Definition yaspgrid.hh:1745
void refineOptions(bool keepPhysicalOverlap)
set options for refinement
Definition yaspgrid.hh:1271
int size(int level, int codim) const
number of entities per level and codim in this process
Definition yaspgrid.hh:1426
YaspGlobalIdSet< YaspGrid< dim, Coordinates > > GlobalIdSetType
Definition yaspgrid.hh:723
YaspGrid(const Coordinates &coordinates, std::bitset< dim > periodic=std::bitset< dim >(0ULL), int overlap=1, Communication comm=Communication(), const Yasp::Partitioning< dim > *partitioner=defaultPartitioner())
Definition yaspgrid.hh:732
bool adapt()
map adapt to global refine
Definition yaspgrid.hh:1307
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lend(int level) const
version without second template parameter for convenience
Definition yaspgrid.hh:1351
ReservedVector< YGridLevel, 32 >::const_iterator YGridLevelIterator
Iterator over the grid levels.
Definition yaspgrid.hh:290
size_t numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition yaspgrid.hh:1458
void communicateCodim(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int level) const
Definition yaspgrid.hh:1493
void makelevel(const Coordinates &coords, std::bitset< dim > periodic, iTupel o_interior, int overlap)
Make a new YGridLevel structure.
Definition yaspgrid.hh:327
bigunsignedint< dim *yaspgrid_dim_bits+yaspgrid_level_bits+dim > PersistentIndexType
Definition yaspgrid.hh:713
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator lend(int level) const
Iterator to one past the last entity of given codim on level for partition type.
Definition yaspgrid.hh:1337
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator lbegin(int level) const
one past the end on this level
Definition yaspgrid.hh:1330
static const Yasp::Partitioning< dim > * defaultPartitioner()
Definition yaspgrid.hh:313
iTupel globalSize() const
return number of cells on finest level on all processors
Definition yaspgrid.hh:258
const Communication & comm() const
return a communication object
Definition yaspgrid.hh:1758
int overlapSize(int odim) const
return size (= distance in graph) of overlap region
Definition yaspgrid.hh:1407
const YaspGrid< dim, Coordinates > GridImp
Definition yaspgrid.hh:676
YGridLevelIterator begin(int i) const
return iterator pointing to given level
Definition yaspgrid.hh:299
Coordinates::ctype ctype
Type used for coordinates.
Definition yaspgrid.hh:179
int levelSize(int l, int i) const
return size of the grid (in cells) on level l in direction i
Definition yaspgrid.hh:264
YGridLevelIterator begin() const
return iterator pointing to coarsest level
Definition yaspgrid.hh:293
The general version that handles all codimensions but 0 and dim.
Definition yaspgridgeometry.hh:31
Definition yaspgridentity.hh:242
Describes the minimal information necessary to create a fully functional YaspEntity.
Definition yaspgridentityseed.hh:18
Iterates over entities of one grid level.
Definition yaspgridleveliterator.hh:19
YaspIntersectionIterator enables iteration over intersections with neighboring codim 0 entities.
Definition yaspgridintersectioniterator.hh:22
YaspIntersection provides data about intersection with neighboring codim 0 entities.
Definition yaspgridintersection.hh:22
YaspHierarchicIterator enables iteration over son entities of codim 0.
Definition yaspgridhierarchiciterator.hh:20
Implementation of Level- and LeafIndexSets for YaspGrid.
Definition yaspgridindexsets.hh:25
persistent, globally unique Ids
Definition yaspgrididset.hh:25
Definition yaspgridpersistentcontainer.hh:35
Definition yaspgrid.hh:92
YaspCommunication CCType
Definition yaspgrid.hh:93
GridTraits< dim, dim, Dune::YaspGrid< dim, Coordinates >, YaspGeometry, YaspEntity, YaspLevelIterator, YaspIntersection, YaspIntersection, YaspIntersectionIterator, YaspIntersectionIterator, YaspHierarchicIterator, YaspLevelIterator, YaspIndexSet< const YaspGrid< dim, Coordinates >, false >, YaspIndexSet< const YaspGrid< dim, Coordinates >, true >, YaspGlobalIdSet< const YaspGrid< dim, Coordinates > >, bigunsignedint< dim *yaspgrid_dim_bits+yaspgrid_level_bits+dim >, YaspGlobalIdSet< const YaspGrid< dim, Coordinates > >, bigunsignedint< dim *yaspgrid_dim_bits+yaspgrid_level_bits+dim >, CCType, DefaultLevelGridViewTraits, DefaultLeafGridViewTraits, YaspEntitySeed, YaspGeometry, unsigned int, std::array< GeometryType, 1 > > Traits
Definition yaspgrid.hh:118
a base class for the yaspgrid partitioning strategy
Definition partitioning.hh:38
Definition partitioning.hh:47
double partition(int rank, iTupel origin_in, iTupel size_in, iTupel &origin_out, iTupel &size_out) const
partition the given grid onto the torus and return the piece of the process with given rank; returns ...
Definition torus.hh:239
iTupel coord() const
return own coordinates
Definition torus.hh:100
int rank() const
return own rank
Definition torus.hh:94
const iTupel & dims() const
return dimensions of torus
Definition torus.hh:112
void recv(int rank, void *buffer, int size) const
store a receive request; buffers are received in order; handles also local requests with memcpy
Definition torus.hh:374
void send(int rank, void *buffer, int size) const
store a send request; buffers are sent in order; handles also local requests with memcpy
Definition torus.hh:361
int neighbors() const
return the number of neighbors, which is
Definition torus.hh:203
ProcListIterator recvend() const
last process in receive list
Definition torus.hh:355
ProcListIterator sendend() const
end of send list
Definition torus.hh:343
ProcListIterator sendbegin() const
first process in send list
Definition torus.hh:337
void exchange() const
exchange messages stored in request buffers; clear request buffers afterwards
Definition torus.hh:387
ProcListIterator recvbegin() const
first process in receive list
Definition torus.hh:349
implements a collection of YGridComponents which form a codimension Entities of given codimension c n...
Definition ygrid.hh:551
Iterator over a collection o YGrids A YGrid::Iterator is the heart of an entity in YaspGrid.
Definition ygrid.hh:594
type describing an intersection with a neighboring processor
Definition ygrid.hh:829
A set of traits classes to store static information about grid implementation.
Different resources needed by all grid implementations.
Provides base classes for index and id sets.