dune-grid  2.4
yaspgrid.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_HH
4 #define DUNE_GRID_YASPGRID_HH
5 
6 #include <iostream>
7 #include <vector>
8 #include <algorithm>
9 #include <stack>
10 
11 // either include stdint.h or provide fallback for uint8_t
12 #if HAVE_STDINT_H
13 #include <stdint.h>
14 #else
15 typedef unsigned char uint8_t;
16 #endif
17 
19 #include <dune/grid/common/grid.hh> // the grid base classes
20 #include <dune/grid/common/capabilities.hh> // the capabilities
21 #include <dune/common/power.hh>
22 #include <dune/common/bigunsignedint.hh>
23 #include <dune/common/typetraits.hh>
24 #include <dune/common/reservedvector.hh>
25 #include <dune/common/parallel/collectivecommunication.hh>
26 #include <dune/common/parallel/mpihelper.hh>
27 #include <dune/common/deprecated.hh>
28 #include <dune/geometry/genericgeometry/topologytypes.hh>
29 #include <dune/geometry/axisalignedcubegeometry.hh>
32 
33 
34 #if HAVE_MPI
35 #include <dune/common/parallel/mpicollectivecommunication.hh>
36 #endif
37 
45 namespace Dune {
46 
47  /* some sizes for building global ids
48  */
49  const int yaspgrid_dim_bits = 24; // bits for encoding each dimension
50  const int yaspgrid_level_bits = 5; // bits for encoding level number
51 
52 
53  //************************************************************************
54  // forward declaration of templates
55 
56  template<int dim, class Coordinates> class YaspGrid;
57  template<int mydim, int cdim, class GridImp> class YaspGeometry;
58  template<int codim, int dim, class GridImp> class YaspEntity;
59  template<int codim, class GridImp> class YaspEntityPointer;
60  template<int codim, class GridImp> class YaspEntitySeed;
61  template<int codim, PartitionIteratorType pitype, class GridImp> class YaspLevelIterator;
62  template<class GridImp> class YaspIntersectionIterator;
63  template<class GridImp> class YaspIntersection;
64  template<class GridImp> class YaspHierarchicIterator;
65  template<class GridImp, bool isLeafIndexSet> class YaspIndexSet;
66  template<class GridImp> class YaspGlobalIdSet;
67  template<class GridImp> class YaspPersistentContainerIndex;
68 
69 } // namespace Dune
70 
85 
86 namespace Dune {
87 
88  template<int dim, class Coordinates>
90  {
91 #if HAVE_MPI
92  typedef CollectiveCommunication<MPI_Comm> CCType;
93 #else
94  typedef CollectiveCommunication<No_Comm> CCType;
95 #endif
96 
97  typedef GridTraits<dim, // dimension of the grid
98  dim, // dimension of the world space
102  YaspLevelIterator, // type used for the level iterator
103  YaspIntersection, // leaf intersection
104  YaspIntersection, // level intersection
105  YaspIntersectionIterator, // leaf intersection iter
106  YaspIntersectionIterator, // level intersection iter
108  YaspLevelIterator, // type used for the leaf(!) iterator
109  YaspIndexSet< const YaspGrid< dim, Coordinates >, false >, // level index set
110  YaspIndexSet< const YaspGrid< dim, Coordinates >, true >, // leaf index set
112  bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim>,
113  YaspGlobalIdSet<const YaspGrid<dim, Coordinates> >,
114  bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim>,
115  CCType,
119  };
120 
121 #ifndef DOXYGEN
122  template<int dim, int codim>
123  struct YaspCommunicateMeta {
124  template<class G, class DataHandle>
125  static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
126  {
127  if (data.contains(dim,codim))
128  {
129  g.template communicateCodim<DataHandle,codim>(data,iftype,dir,level);
130  }
131  YaspCommunicateMeta<dim,codim-1>::comm(g,data,iftype,dir,level);
132  }
133  };
134 
135  template<int dim>
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)
139  {
140  if (data.contains(dim,0))
141  g.template communicateCodim<DataHandle,0>(data,iftype,dir,level);
142  }
143  };
144 #endif
145 
146  //************************************************************************
163  template<int dim, class Coordinates = EquidistantCoordinates<double, dim> >
164  class YaspGrid
165  : public GridDefaultImplementation<dim,dim,typename Coordinates::ctype,YaspGridFamily<dim, Coordinates> >
166  {
167 
168  template<int, PartitionIteratorType, typename>
169  friend class YaspLevelIterator;
170 
171  template<typename>
173 
174  protected:
175 
177 
178  public:
180  typedef typename Coordinates::ctype ctype;
181 #if HAVE_MPI
182  typedef CollectiveCommunication<MPI_Comm> CollectiveCommunicationType;
183 #else
184  typedef CollectiveCommunication<No_Comm> CollectiveCommunicationType;
185 #endif
186 
187 #ifndef DOXYGEN
188  typedef typename Dune::YGrid<Coordinates> YGrid;
190 
193  struct YGridLevel {
194 
196  int level() const
197  {
198  return level_;
199  }
200 
201  Coordinates coords;
202 
203  std::array<YGrid, dim+1> overlapfront;
204  std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power> overlapfront_data;
205  std::array<YGrid, dim+1> overlap;
206  std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power> overlap_data;
207  std::array<YGrid, dim+1> interiorborder;
208  std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power> interiorborder_data;
209  std::array<YGrid, dim+1> interior;
210  std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power> interior_data;
211 
212  std::array<YGridList<Coordinates>,dim+1> send_overlapfront_overlapfront;
213  std::array<std::deque<Intersection>, StaticPower<2,dim>::power> send_overlapfront_overlapfront_data;
214  std::array<YGridList<Coordinates>,dim+1> recv_overlapfront_overlapfront;
215  std::array<std::deque<Intersection>, StaticPower<2,dim>::power> recv_overlapfront_overlapfront_data;
216 
217  std::array<YGridList<Coordinates>,dim+1> send_overlap_overlapfront;
218  std::array<std::deque<Intersection>, StaticPower<2,dim>::power> send_overlap_overlapfront_data;
219  std::array<YGridList<Coordinates>,dim+1> recv_overlapfront_overlap;
220  std::array<std::deque<Intersection>, StaticPower<2,dim>::power> recv_overlapfront_overlap_data;
221 
222  std::array<YGridList<Coordinates>,dim+1> send_interiorborder_interiorborder;
223  std::array<std::deque<Intersection>, StaticPower<2,dim>::power> send_interiorborder_interiorborder_data;
224  std::array<YGridList<Coordinates>,dim+1> recv_interiorborder_interiorborder;
225  std::array<std::deque<Intersection>, StaticPower<2,dim>::power> recv_interiorborder_interiorborder_data;
226 
227  std::array<YGridList<Coordinates>,dim+1> send_interiorborder_overlapfront;
228  std::array<std::deque<Intersection>, StaticPower<2,dim>::power> send_interiorborder_overlapfront_data;
229  std::array<YGridList<Coordinates>,dim+1> recv_overlapfront_interiorborder;
230  std::array<std::deque<Intersection>, StaticPower<2,dim>::power> recv_overlapfront_interiorborder_data;
231 
232  // general
233  YaspGrid<dim,Coordinates>* mg; // each grid level knows its multigrid
234  int overlapSize; // in mesh cells on this level
235  bool keepOverlap;
236 
238  int level_;
239  };
240 
242  typedef std::array<int, dim> iTupel;
243  typedef FieldVector<ctype, dim> fTupel;
244 
245  // communication tag used by multigrid
246  enum { tag = 17 };
247 #endif
248 
251  {
252  return _torus;
253  }
254 
256  int globalSize(int i) const
257  {
258  return levelSize(maxLevel(),i);
259  }
260 
262  iTupel globalSize() const
263  {
264  return levelSize(maxLevel());
265  }
266 
268  int levelSize(int l, int i) const
269  {
270  return _coarseSize[i] * (1 << l);
271  }
272 
274  iTupel levelSize(int l) const
275  {
276  iTupel s;
277  for (int i=0; i<dim; ++i)
278  s[i] = levelSize(l,i);
279  return s;
280  }
281 
283  bool isPeriodic(int i) const
284  {
285  return _periodic[i];
286  }
287 
288  bool getRefineOption() const
289  {
290  return keep_ovlp;
291  }
292 
294  typedef typename ReservedVector<YGridLevel,32>::const_iterator YGridLevelIterator;
295 
297  YGridLevelIterator begin () const
298  {
299  return YGridLevelIterator(_levels,0);
300  }
301 
303  YGridLevelIterator begin (int i) const
304  {
305  if (i<0 || i>maxLevel())
306  DUNE_THROW(GridError, "level not existing");
307  return YGridLevelIterator(_levels,i);
308  }
309 
311  YGridLevelIterator end () const
312  {
313  return YGridLevelIterator(_levels,maxLevel()+1);
314  }
315 
316  // static method to create the default load balance strategy
318  {
319  static YLoadBalanceDefault<dim> lb;
320  return & lb;
321  }
322 
323  protected:
331  void makelevel (const Coordinates& coords, std::bitset<dim> periodic, iTupel o_interior, int overlap)
332  {
333  YGridLevel& g = _levels.back();
334  g.overlapSize = overlap;
335  g.mg = this;
336  g.level_ = maxLevel();
337  g.coords = coords;
338  g.keepOverlap = keep_ovlp;
339 
340  // set the inserting positions in the corresponding arrays of YGridLevelStructure
341  typename std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power>::iterator overlapfront_it = g.overlapfront_data.begin();
342  typename std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power>::iterator overlap_it = g.overlap_data.begin();
343  typename std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power>::iterator interiorborder_it = g.interiorborder_data.begin();
344  typename std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power>::iterator interior_it = g.interior_data.begin();
345 
346  typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator
347  send_overlapfront_overlapfront_it = g.send_overlapfront_overlapfront_data.begin();
348  typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator
349  recv_overlapfront_overlapfront_it = g.recv_overlapfront_overlapfront_data.begin();
350 
351  typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator
352  send_overlap_overlapfront_it = g.send_overlap_overlapfront_data.begin();
353  typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator
354  recv_overlapfront_overlap_it = g.recv_overlapfront_overlap_data.begin();
355 
356  typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator
357  send_interiorborder_interiorborder_it = g.send_interiorborder_interiorborder_data.begin();
358  typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator
359  recv_interiorborder_interiorborder_it = g.recv_interiorborder_interiorborder_data.begin();
360 
361  typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator
362  send_interiorborder_overlapfront_it = g.send_interiorborder_overlapfront_data.begin();
363  typename std::array<std::deque<Intersection>, StaticPower<2,dim>::power>::iterator
364  recv_overlapfront_interiorborder_it = g.recv_overlapfront_interiorborder_data.begin();
365 
366  // have a null array for constructor calls around
367  std::array<int,dim> n;
368  std::fill(n.begin(), n.end(), 0);
369 
370  // determine origin of the grid with overlap and store whether an overlap area exists in direction i.
371  std::bitset<dim> ovlp_low(0ULL);
372  std::bitset<dim> ovlp_up(0ULL);
373 
374  iTupel o_overlap;
375  iTupel s_overlap;
376 
377  // determine at where we have overlap and how big the size of the overlap partition is
378  for (int i=0; i<dim; i++)
379  {
380  // the coordinate container has been contructed to hold the entire grid on
381  // this processor, including overlap. this is the element size.
382  s_overlap[i] = g.coords.size(i);
383 
384  //in the periodic case there is always overlap
385  if (periodic[i])
386  {
387  o_overlap[i] = o_interior[i]-overlap;
388  ovlp_low[i] = true;
389  ovlp_up[i] = true;
390  }
391  else
392  {
393  //check lower boundary
394  if (o_interior[i] - overlap < 0)
395  o_overlap[i] = 0;
396  else
397  {
398  o_overlap[i] = o_interior[i] - overlap;
399  ovlp_low[i] = true;
400  }
401 
402  //check upper boundary
403  if (o_overlap[i] + g.coords.size(i) < globalSize(i))
404  ovlp_up[i] = true;
405  }
406  }
407 
408  for (unsigned int codim = 0; codim < dim + 1; codim++)
409  {
410  // set the begin iterator for the corresponding ygrids
411  g.overlapfront[codim].setBegin(overlapfront_it);
412  g.overlap[codim].setBegin(overlap_it);
413  g.interiorborder[codim].setBegin(interiorborder_it);
414  g.interior[codim].setBegin(interior_it);
415  g.send_overlapfront_overlapfront[codim].setBegin(send_overlapfront_overlapfront_it);
416  g.recv_overlapfront_overlapfront[codim].setBegin(recv_overlapfront_overlapfront_it);
417  g.send_overlap_overlapfront[codim].setBegin(send_overlap_overlapfront_it);
418  g.recv_overlapfront_overlap[codim].setBegin(recv_overlapfront_overlap_it);
419  g.send_interiorborder_interiorborder[codim].setBegin(send_interiorborder_interiorborder_it);
420  g.recv_interiorborder_interiorborder[codim].setBegin(recv_interiorborder_interiorborder_it);
421  g.send_interiorborder_overlapfront[codim].setBegin(send_interiorborder_overlapfront_it);
422  g.recv_overlapfront_interiorborder[codim].setBegin(recv_overlapfront_interiorborder_it);
423 
424  // find all combinations of unit vectors that span entities of the given codimension
425  for (unsigned int index = 0; index < (1<<dim); index++)
426  {
427  // check whether the given shift is of our codimension
428  std::bitset<dim> r(index);
429  if (r.count() != dim-codim)
430  continue;
431 
432  // get an origin and a size array for subsequent modification
433  std::array<int,dim> origin(o_overlap);
434  std::array<int,dim> size(s_overlap);
435 
436  // build overlapfront
437  // we have to extend the element size by one in all directions without shift.
438  for (int i=0; i<dim; i++)
439  if (!r[i])
440  size[i]++;
441  *overlapfront_it = YGridComponent<Coordinates>(origin, r, &g.coords, size, n, size);
442 
443  // build overlap
444  for (int i=0; i<dim; i++)
445  {
446  if (!r[i])
447  {
448  if (ovlp_low[i])
449  {
450  origin[i]++;
451  size[i]--;
452  }
453  if (ovlp_up[i])
454  size[i]--;
455  }
456  }
457  *overlap_it = YGridComponent<Coordinates>(origin,size,*overlapfront_it);
458 
459  // build interiorborder
460  for (int i=0; i<dim; i++)
461  {
462  if (ovlp_low[i])
463  {
464  origin[i] += overlap;
465  size[i] -= overlap;
466  if (!r[i])
467  {
468  origin[i]--;
469  size[i]++;
470  }
471  }
472  if (ovlp_up[i])
473  {
474  size[i] -= overlap;
475  if (!r[i])
476  size[i]++;
477  }
478  }
479  *interiorborder_it = YGridComponent<Coordinates>(origin,size,*overlapfront_it);
480 
481  // build interior
482  for (int i=0; i<dim; i++)
483  {
484  if (!r[i])
485  {
486  if (ovlp_low[i])
487  {
488  origin[i]++;
489  size[i]--;
490  }
491  if (ovlp_up[i])
492  size[i]--;
493  }
494  }
495  *interior_it = YGridComponent<Coordinates>(origin, size, *overlapfront_it);
496 
497  intersections(*overlapfront_it,*overlapfront_it,*send_overlapfront_overlapfront_it, *recv_overlapfront_overlapfront_it);
498  intersections(*overlap_it,*overlapfront_it,*send_overlap_overlapfront_it, *recv_overlapfront_overlap_it);
499  intersections(*interiorborder_it,*interiorborder_it,*send_interiorborder_interiorborder_it,*recv_interiorborder_interiorborder_it);
500  intersections(*interiorborder_it,*overlapfront_it,*send_interiorborder_overlapfront_it,*recv_overlapfront_interiorborder_it);
501 
502  // advance all iterators pointing to the next insertion point
503  ++overlapfront_it;
504  ++overlap_it;
505  ++interiorborder_it;
506  ++interior_it;
507  ++send_overlapfront_overlapfront_it;
508  ++recv_overlapfront_overlapfront_it;
509  ++send_overlap_overlapfront_it;
510  ++recv_overlapfront_overlap_it;
511  ++send_interiorborder_interiorborder_it;
512  ++recv_interiorborder_interiorborder_it;
513  ++send_interiorborder_overlapfront_it;
514  ++recv_overlapfront_interiorborder_it;
515  }
516 
517  // set end iterators in the corresonding ygrids
518  g.overlapfront[codim].finalize(overlapfront_it);
519  g.overlap[codim].finalize(overlap_it);
520  g.interiorborder[codim].finalize(interiorborder_it);
521  g.interior[codim].finalize(interior_it);
522  g.send_overlapfront_overlapfront[codim].finalize(send_overlapfront_overlapfront_it,g.overlapfront[codim]);
523  g.recv_overlapfront_overlapfront[codim].finalize(recv_overlapfront_overlapfront_it,g.overlapfront[codim]);
524  g.send_overlap_overlapfront[codim].finalize(send_overlap_overlapfront_it,g.overlapfront[codim]);
525  g.recv_overlapfront_overlap[codim].finalize(recv_overlapfront_overlap_it,g.overlapfront[codim]);
526  g.send_interiorborder_interiorborder[codim].finalize(send_interiorborder_interiorborder_it,g.overlapfront[codim]);
527  g.recv_interiorborder_interiorborder[codim].finalize(recv_interiorborder_interiorborder_it,g.overlapfront[codim]);
528  g.send_interiorborder_overlapfront[codim].finalize(send_interiorborder_overlapfront_it,g.overlapfront[codim]);
529  g.recv_overlapfront_interiorborder[codim].finalize(recv_overlapfront_interiorborder_it,g.overlapfront[codim]);
530  }
531  }
532 
533 #ifndef DOXYGEN
534 
542  struct mpifriendly_ygrid {
543  mpifriendly_ygrid ()
544  {
545  std::fill(origin.begin(), origin.end(), 0);
546  std::fill(size.begin(), size.end(), 0);
547  }
548  mpifriendly_ygrid (const YGridComponent<Coordinates>& grid)
549  : origin(grid.origin()), size(grid.size())
550  {}
551  iTupel origin;
552  iTupel size;
553  };
554 #endif
555 
565  std::deque<Intersection>& sendlist, std::deque<Intersection>& recvlist)
566  {
567  iTupel size = globalSize();
568 
569  // the exchange buffers
570  std::vector<YGridComponent<Coordinates> > send_recvgrid(_torus.neighbors());
571  std::vector<YGridComponent<Coordinates> > recv_recvgrid(_torus.neighbors());
572  std::vector<YGridComponent<Coordinates> > send_sendgrid(_torus.neighbors());
573  std::vector<YGridComponent<Coordinates> > recv_sendgrid(_torus.neighbors());
574 
575  // new exchange buffers to send simple struct without virtual functions
576  std::vector<mpifriendly_ygrid> mpifriendly_send_recvgrid(_torus.neighbors());
577  std::vector<mpifriendly_ygrid> mpifriendly_recv_recvgrid(_torus.neighbors());
578  std::vector<mpifriendly_ygrid> mpifriendly_send_sendgrid(_torus.neighbors());
579  std::vector<mpifriendly_ygrid> mpifriendly_recv_sendgrid(_torus.neighbors());
580 
581  // fill send buffers; iterate over all neighboring processes
582  // non-periodic case is handled automatically because intersection will be zero
583  for (typename Torus<CollectiveCommunicationType,dim>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
584  {
585  // determine if we communicate with this neighbor (and what)
586  bool skip = false;
587  iTupel coord = _torus.coord(); // my coordinates
588  iTupel delta = i.delta(); // delta to neighbor
589  iTupel nb = coord; // the neighbor
590  for (int k=0; k<dim; k++) nb[k] += delta[k];
591  iTupel v; // grid movement
592  std::fill(v.begin(), v.end(), 0);
593 
594  for (int k=0; k<dim; k++)
595  {
596  if (nb[k]<0)
597  {
598  if (_periodic[k])
599  v[k] += size[k];
600  else
601  skip = true;
602  }
603  if (nb[k]>=_torus.dims(k))
604  {
605  if (_periodic[k])
606  v[k] -= size[k];
607  else
608  skip = true;
609  }
610  // neither might be true, then v=0
611  }
612 
613  // store moved grids in send buffers
614  if (!skip)
615  {
616  send_sendgrid[i.index()] = sendgrid.move(v);
617  send_recvgrid[i.index()] = recvgrid.move(v);
618  }
619  else
620  {
621  send_sendgrid[i.index()] = YGridComponent<Coordinates>();
622  send_recvgrid[i.index()] = YGridComponent<Coordinates>();
623  }
624  }
625 
626  // issue send requests for sendgrid being sent to all neighbors
627  for (typename Torus<CollectiveCommunicationType,dim>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
628  {
629  mpifriendly_send_sendgrid[i.index()] = mpifriendly_ygrid(send_sendgrid[i.index()]);
630  _torus.send(i.rank(), &mpifriendly_send_sendgrid[i.index()], sizeof(mpifriendly_ygrid));
631  }
632 
633  // issue recv requests for sendgrids of neighbors
634  for (typename Torus<CollectiveCommunicationType,dim>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
635  _torus.recv(i.rank(), &mpifriendly_recv_sendgrid[i.index()], sizeof(mpifriendly_ygrid));
636 
637  // exchange the sendgrids
638  _torus.exchange();
639 
640  // issue send requests for recvgrid being sent to all neighbors
641  for (typename Torus<CollectiveCommunicationType,dim>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
642  {
643  mpifriendly_send_recvgrid[i.index()] = mpifriendly_ygrid(send_recvgrid[i.index()]);
644  _torus.send(i.rank(), &mpifriendly_send_recvgrid[i.index()], sizeof(mpifriendly_ygrid));
645  }
646 
647  // issue recv requests for recvgrid of neighbors
648  for (typename Torus<CollectiveCommunicationType,dim>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
649  _torus.recv(i.rank(), &mpifriendly_recv_recvgrid[i.index()], sizeof(mpifriendly_ygrid));
650 
651  // exchange the recvgrid
652  _torus.exchange();
653 
654  // process receive buffers and compute intersections
655  for (typename Torus<CollectiveCommunicationType,dim>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
656  {
657  // what must be sent to this neighbor
658  Intersection send_intersection;
659  mpifriendly_ygrid yg = mpifriendly_recv_recvgrid[i.index()];
660  recv_recvgrid[i.index()] = YGridComponent<Coordinates>(yg.origin,yg.size);
661  send_intersection.grid = sendgrid.intersection(recv_recvgrid[i.index()]);
662  send_intersection.rank = i.rank();
663  send_intersection.distance = i.distance();
664  if (!send_intersection.grid.empty()) sendlist.push_front(send_intersection);
665 
666  Intersection recv_intersection;
667  yg = mpifriendly_recv_sendgrid[i.index()];
668  recv_sendgrid[i.index()] = YGridComponent<Coordinates>(yg.origin,yg.size);
669  recv_intersection.grid = recvgrid.intersection(recv_sendgrid[i.index()]);
670  recv_intersection.rank = i.rank();
671  recv_intersection.distance = i.distance();
672  if(!recv_intersection.grid.empty()) recvlist.push_back(recv_intersection);
673  }
674  }
675 
676  protected:
677 
679 
680  void init()
681  {
682  Yasp::BinomialTable<dim>::init();
683  Yasp::EntityShiftTable<Yasp::calculate_entity_shift<dim>,dim>::init();
684  Yasp::EntityShiftTable<Yasp::calculate_entity_move<dim>,dim>::init();
685  indexsets.push_back( std::make_shared< YaspIndexSet<const YaspGrid<dim, Coordinates>, false > >(*this,0) );
687  }
688 
690  {
691  // sizes of local macro grid
692  std::array<int, dim> sides;
693  {
694  for (int i=0; i<dim; i++)
695  {
696  sides[i] =
697  ((begin()->overlap[0].dataBegin()->origin(i) == 0)+
698  (begin()->overlap[0].dataBegin()->origin(i) + begin()->overlap[0].dataBegin()->size(i)
699  == levelSize(0,i)));
700  }
701  }
702  nBSegments = 0;
703  for (int k=0; k<dim; k++)
704  {
705  int offset = 1;
706  for (int l=0; l<dim; l++)
707  {
708  if (l==k) continue;
709  offset *= begin()->overlap[0].dataBegin()->size(l);
710  }
711  nBSegments += sides[k]*offset;
712  }
713  }
714 
715  public:
716 
717  // define the persistent index type
718  typedef bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim> PersistentIndexType;
719 
722  // the Traits
724 
725  // need for friend declarations in entity
727  typedef YaspIndexSet<YaspGrid<dim, Coordinates>, true > LeafIndexSetType;
729 
738  YaspGrid (Dune::FieldVector<ctype, dim> L,
739  std::array<int, dim> s,
740  std::bitset<dim> periodic = std::bitset<dim>(0ULL),
741  int overlap = 1,
742  CollectiveCommunicationType comm = CollectiveCommunicationType(),
744  : ccobj(comm), _torus(comm,tag,s,lb), leafIndexSet_(*this),
745  _L(L), _periodic(periodic), _coarseSize(s), _overlap(overlap),
746  keep_ovlp(true), adaptRefCount(0), adaptActive(false)
747  {
748  // check whether YaspGrid has been given the correct template parameter
749  static_assert(is_same<Coordinates,EquidistantCoordinates<ctype,dim> >::value,
750  "YaspGrid coordinate container template parameter and given constructor values do not match!");
751 
752  _levels.resize(1);
753 
754  iTupel o;
755  std::fill(o.begin(), o.end(), 0);
756  iTupel o_interior(o);
757  iTupel s_interior(s);
758 
759  _torus.partition(_torus.rank(),o,s,o_interior,s_interior);
760 
761  // check whether the grid is large enough to be overlapping
762  for (int i=0; i<dim; i++)
763  if ((s_interior[i] <= overlap) && // interior is very small
764  (periodic[i] || (s_interior[i] != s[i]))) // there is an overlap in that direction
765  DUNE_THROW(Dune::GridError,"YaspGrid is too small to be overlapping");
766 
767  fTupel h(L);
768  for (int i=0; i<dim; i++)
769  h[i] /= s[i];
770 
771  iTupel s_overlap(s_interior);
772  for (int i=0; i<dim; i++)
773  {
774  if ((o_interior[i] - overlap > 0) || (periodic[i]))
775  s_overlap[i] += overlap;
776  if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
777  s_overlap[i] += overlap;
778  }
779 
780  EquidistantCoordinates<ctype,dim> cc(h,s_overlap);
781 
782  // add level
783  makelevel(cc,periodic,o_interior,overlap);
784 
785  init();
786  }
787 
797  YaspGrid (Dune::FieldVector<ctype, dim> lowerleft,
798  Dune::FieldVector<ctype, dim> upperright,
799  std::array<int, dim> s,
800  std::bitset<dim> periodic = std::bitset<dim>(0ULL),
801  int overlap = 1,
802  CollectiveCommunicationType comm = CollectiveCommunicationType(),
804  : ccobj(comm), _torus(comm,tag,s,lb), leafIndexSet_(*this),
805  _L(upperright - lowerleft),
806  _periodic(periodic), _coarseSize(s), _overlap(overlap),
807  keep_ovlp(true), adaptRefCount(0), adaptActive(false)
808  {
809  // check whether YaspGrid has been given the correct template parameter
810  static_assert(is_same<Coordinates,EquidistantOffsetCoordinates<ctype,dim> >::value,
811  "YaspGrid coordinate container template parameter and given constructor values do not match!");
812 
813  _levels.resize(1);
814 
815  iTupel o;
816  std::fill(o.begin(), o.end(), 0);
817  iTupel o_interior(o);
818  iTupel s_interior(s);
819 
820  _torus.partition(_torus.rank(),o,s,o_interior,s_interior);
821 
822  // check whether the grid is large enough to be overlapping
823  for (int i=0; i<dim; i++)
824  if ((s_interior[i] <= overlap) && // interior is very small
825  (periodic[i] || (s_interior[i] != s[i]))) // there is an overlap in that direction
826  DUNE_THROW(Dune::GridError,"YaspGrid is too small to be overlapping");
827 
828  Dune::FieldVector<ctype,dim> extension(upperright);
829  Dune::FieldVector<ctype,dim> h;
830  for (int i=0; i<dim; i++)
831  {
832  extension[i] -= lowerleft[i];
833  h[i] = extension[i] / s[i];
834  }
835 
836  iTupel s_overlap(s_interior);
837  for (int i=0; i<dim; i++)
838  {
839  if ((o_interior[i] - overlap > 0) || (periodic[i]))
840  s_overlap[i] += overlap;
841  if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
842  s_overlap[i] += overlap;
843  }
844 
845  EquidistantOffsetCoordinates<ctype,dim> cc(lowerleft,h,s_overlap);
846 
847  // add level
848  makelevel(cc,periodic,o_interior,overlap);
849 
850  init();
851  }
852 
860  YaspGrid (std::array<std::vector<ctype>, dim> coords,
861  std::bitset<dim> periodic = std::bitset<dim>(0ULL),
862  int overlap = 1,
863  CollectiveCommunicationType comm = CollectiveCommunicationType(),
865  : ccobj(comm), _torus(comm,tag,Dune::Yasp::sizeArray<dim>(coords),lb),
866  leafIndexSet_(*this), _periodic(periodic), _overlap(overlap),
867  keep_ovlp(true), adaptRefCount(0), adaptActive(false)
868  {
869  if (!Dune::Yasp::checkIfMonotonous(coords))
870  DUNE_THROW(Dune::GridError,"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
871 
872  // check whether YaspGrid has been given the correct template parameter
873  static_assert(is_same<Coordinates,TensorProductCoordinates<ctype,dim> >::value,
874  "YaspGrid coordinate container template parameter and given constructor values do not match!");
875 
876  _levels.resize(1);
877 
878  //determine sizes of vector to correctly construct torus structure and store for later size requests
879  for (int i=0; i<dim; i++) {
880  _coarseSize[i] = coords[i].size() - 1;
881  _L[i] = coords[i][_coarseSize[i]] - coords[i][0];
882  }
883 
884  iTupel o;
885  std::fill(o.begin(), o.end(), 0);
886  iTupel o_interior(o);
887  iTupel s_interior(_coarseSize);
888 
889  _torus.partition(_torus.rank(),o,_coarseSize,o_interior,s_interior);
890 
891  // check whether the grid is large enough to be overlapping
892  for (int i=0; i<dim; i++)
893  if ((s_interior[i] <= overlap) && // interior is very small
894  (periodic[i] || (s_interior[i] != _coarseSize[i]))) // there is an overlap in that direction
895  DUNE_THROW(Dune::GridError,"YaspGrid is too small to be overlapping");
896 
897  std::array<std::vector<ctype>,dim> newcoords;
898  std::array<int, dim> offset(o_interior);
899 
900  // find the relevant part of the coords vector for this processor and copy it to newcoords
901  for (int i=0; i<dim; ++i)
902  {
903  //define iterators on coords that specify the coordinate range to be used
904  typename std::vector<ctype>::iterator begin = coords[i].begin() + o_interior[i];
905  typename std::vector<ctype>::iterator end = begin + s_interior[i] + 1;
906 
907  // check whether we are not at the physical boundary. In that case overlap is a simple
908  // extension of the coordinate range to be used
909  if (o_interior[i] - overlap > 0)
910  {
911  begin = begin - overlap;
912  offset[i] -= overlap;
913  }
914  if (o_interior[i] + s_interior[i] + overlap < _coarseSize[i])
915  end = end + overlap;
916 
917  //copy the selected part in the new coord vector
918  newcoords[i].resize(end-begin);
919  std::copy(begin, end, newcoords[i].begin());
920 
921  // check whether we are at the physical boundary and a have a periodic grid.
922  // In this case the coordinate vector has to be tweaked manually.
923  if ((periodic[i]) && (o_interior[i] + s_interior[i] + overlap >= _coarseSize[i]))
924  {
925  // we need to add the first <overlap> cells to the end of newcoords
926  typename std::vector<ctype>::iterator it = coords[i].begin();
927  for (int j=0; j<overlap; ++j)
928  newcoords[i].push_back(newcoords[i].back() - *it + *(++it));
929  }
930 
931  if ((periodic[i]) && (o_interior[i] - overlap <= 0))
932  {
933  offset[i] -= overlap;
934 
935  // we need to add the last <overlap> cells to the begin of newcoords
936  typename std::vector<ctype>::iterator it = coords[i].end() - 1;
937  for (int j=0; j<overlap; ++j)
938  newcoords[i].insert(newcoords[i].begin(), newcoords[i].front() - *it + *(--it));
939  }
940  }
941 
942  TensorProductCoordinates<ctype,dim> cc(newcoords, offset);
943 
944  // add level
945  makelevel(cc,periodic,o_interior,overlap);
946  init();
947  }
948 
960  DUNE_DEPRECATED_MSG("This Yaspgrid constructor is deprecated.")
961  YaspGrid (Dune::MPIHelper::MPICommunicator comm,
962  Dune::FieldVector<ctype, dim> L,
963  std::array<int, dim> s,
964  std::bitset<dim> periodic,
965  int overlap,
966  const YLoadBalance<dim>* lb = defaultLoadbalancer())
967  : ccobj(comm), _torus(comm,tag,s,lb), leafIndexSet_(*this),
968  _L(L), keep_ovlp(true), adaptRefCount(0), adaptActive(false)
969  {
970  _periodic = periodic;
971  _levels.resize(1);
972  _overlap = overlap;
973  _coarseSize = s;
974 
975  iTupel o;
976  std::fill(o.begin(), o.end(), 0);
977  iTupel o_interior(o);
978  iTupel s_interior(s);
979 
980  _torus.partition(_torus.rank(),o,s,o_interior,s_interior);
981 
982  fTupel h(L);
983  for (int i=0; i<dim; i++)
984  h[i] /= s[i];
985 
986  iTupel s_overlap(s_interior);
987  for (int i=0; i<dim; i++)
988  {
989  if ((o_interior[i] - overlap > 0) || (periodic[i]))
990  s_overlap[i] += overlap;
991  if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
992  s_overlap[i] += overlap;
993  }
994 
995  // check whether YaspGrid has been given the correct template parameter
996  static_assert(is_same<Coordinates,EquidistantCoordinates<ctype,dim> >::value,
997  "YaspGrid coordinate container template parameter and given constructor values do not match!");
998 
999  EquidistantCoordinates<ctype,dim> cc(h,s_overlap);
1000 
1001  // add level
1002  makelevel(cc,periodic,o_interior,overlap);
1003  init();
1004  }
1005 
1006 
1017  DUNE_DEPRECATED_MSG("This Yaspgrid constructor is deprecated.")
1018  YaspGrid (Dune::MPIHelper::MPICommunicator comm,
1019  std::array<std::vector<ctype>, dim> coords,
1020  std::bitset<dim> periodic, int overlap,
1021  const YLoadBalance<dim>* lb = defaultLoadbalancer())
1022  : ccobj(comm), _torus(comm,tag,Dune::Yasp::sizeArray<dim>(coords),lb),
1023  leafIndexSet_(*this),
1024  _periodic(std::bitset<dim>(0)),
1025  _overlap(overlap),
1026  keep_ovlp(true),
1027  adaptRefCount(0), adaptActive(false)
1028  {
1029  if (!Dune::Yasp::checkIfMonotonous(coords))
1030  DUNE_THROW(Dune::GridError,"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
1031  _periodic = periodic;
1032  _levels.resize(1);
1033  _overlap = overlap;
1034 
1035  //determine sizes of vector to correctly construct torus structure and store for later size requests
1036  for (int i=0; i<dim; i++) {
1037  _coarseSize[i] = coords[i].size() - 1;
1038  _L[i] = coords[i][_coarseSize[i]] - coords[i][0];
1039  }
1040 
1041  iTupel o;
1042  std::fill(o.begin(), o.end(), 0);
1043  iTupel o_interior(o);
1044  iTupel s_interior(_coarseSize);
1045 
1046  _torus.partition(_torus.rank(),o,_coarseSize,o_interior,s_interior);
1047 
1048  std::array<std::vector<ctype>,dim> newcoords;
1049  std::array<int, dim> offset(o_interior);
1050 
1051  // find the relevant part of the coords vector for this processor and copy it to newcoords
1052  for (int i=0; i<dim; ++i)
1053  {
1054  //define iterators on coords that specify the coordinate range to be used
1055  typename std::vector<ctype>::iterator begin = coords[i].begin() + o_interior[i];
1056  typename std::vector<ctype>::iterator end = begin + s_interior[i] + 1;
1057 
1058  // check whether we are not at the physical boundary. In that case overlap is a simple
1059  // extension of the coordinate range to be used
1060  if (o_interior[i] - overlap > 0)
1061  {
1062  begin = begin - overlap;
1063  offset[i] -= overlap;
1064  }
1065  if (o_interior[i] + s_interior[i] + overlap < _coarseSize[i])
1066  end = end + overlap;
1067 
1068  //copy the selected part in the new coord vector
1069  newcoords[i].resize(end-begin);
1070  std::copy(begin, end, newcoords[i].begin());
1071 
1072  // check whether we are at the physical boundary and a have a periodic grid.
1073  // In this case the coordinate vector has to be tweaked manually.
1074  if ((periodic[i]) && (o_interior[i] + s_interior[i] + overlap >= _coarseSize[i]))
1075  {
1076  // we need to add the first <overlap> cells to the end of newcoords
1077  typename std::vector<ctype>::iterator it = coords[i].begin();
1078  for (int j=0; j<overlap; ++j)
1079  newcoords[i].push_back(newcoords[i].back() - *it + *(++it));
1080  }
1081 
1082  if ((periodic[i]) && (o_interior[i] - overlap <= 0))
1083  {
1084  offset[i] -= overlap;
1085 
1086  // we need to add the last <overlap> cells to the begin of newcoords
1087  typename std::vector<ctype>::iterator it = coords[i].end() - 1;
1088  for (int j=0; j<overlap; ++j)
1089  newcoords[i].insert(newcoords[i].begin(), newcoords[i].front() - *it + *(--it));
1090  }
1091  }
1092 
1093  // check whether YaspGrid has been given the correct template parameter
1094  static_assert(is_same<Coordinates,TensorProductCoordinates<ctype,dim> >::value,
1095  "YaspGrid coordinate container template parameter and given constructor values do not match!");
1096 
1097  TensorProductCoordinates<ctype,dim> cc(newcoords, offset);
1098 
1099  // add level
1100  makelevel(cc,periodic,o_interior,overlap);
1101  init();
1102  }
1103 
1104  private:
1105 
1120  YaspGrid (std::array<std::vector<ctype>, dim> coords,
1121  std::bitset<dim> periodic,
1122  int overlap,
1123  CollectiveCommunicationType comm,
1124  std::array<int,dim> coarseSize,
1125  const YLoadBalance<dim>* lb = defaultLoadbalancer())
1126  : ccobj(comm), _torus(comm,tag,coarseSize,lb), leafIndexSet_(*this),
1127  _periodic(periodic), _coarseSize(coarseSize), _overlap(overlap),
1128  keep_ovlp(true), adaptRefCount(0), adaptActive(false)
1129  {
1130  // check whether YaspGrid has been given the correct template parameter
1131  static_assert(is_same<Coordinates,TensorProductCoordinates<ctype,dim> >::value,
1132  "YaspGrid coordinate container template parameter and given constructor values do not match!");
1133 
1134  if (!Dune::Yasp::checkIfMonotonous(coords))
1135  DUNE_THROW(Dune::GridError,"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
1136 
1137  for (int i=0; i<dim; i++)
1138  _L[i] = coords[i][coords[i].size() - 1] - coords[i][0];
1139 
1140  _levels.resize(1);
1141 
1142  std::array<int,dim> o;
1143  std::fill(o.begin(), o.end(), 0);
1144  std::array<int,dim> o_interior(o);
1145  std::array<int,dim> s_interior(coarseSize);
1146 
1147  _torus.partition(_torus.rank(),o,coarseSize,o_interior,s_interior);
1148 
1149  // get offset by modifying o_interior accoring to overlap
1150  std::array<int,dim> offset(o_interior);
1151  for (int i=0; i<dim; i++)
1152  if ((periodic[i]) || (o_interior[i] > 0))
1153  offset[i] -= overlap;
1154 
1155  TensorProductCoordinates<ctype,dim> cc(coords, offset);
1156 
1157  // add level
1158  makelevel(cc,periodic,o_interior,overlap);
1159 
1160  init();
1161  }
1162 
1163  // the backup restore facility needs to be able to use above constructor
1164  friend struct BackupRestoreFacility<YaspGrid<dim,Coordinates> >;
1165 
1166  // do not copy this class
1167  YaspGrid(const YaspGrid&);
1168 
1169  public:
1170 
1174  int maxLevel() const
1175  {
1176  return _levels.size()-1;
1177  }
1178 
1180  void globalRefine (int refCount)
1181  {
1182  if (refCount < -maxLevel())
1183  DUNE_THROW(GridError, "Only " << maxLevel() << " levels left. " <<
1184  "Coarsening " << -refCount << " levels requested!");
1185 
1186  // If refCount is negative then coarsen the grid
1187  for (int k=refCount; k<0; k++)
1188  {
1189  // create an empty grid level
1190  YGridLevel empty;
1191  _levels.back() = empty;
1192  // reduce maxlevel
1193  _levels.pop_back();
1194 
1195  indexsets.pop_back();
1196  }
1197 
1198  // If refCount is positive refine the grid
1199  for (int k=0; k<refCount; k++)
1200  {
1201  // access to coarser grid level
1202  YGridLevel& cg = _levels[maxLevel()];
1203 
1204  std::bitset<dim> ovlp_low(0ULL), ovlp_up(0ULL);
1205  for (int i=0; i<dim; i++)
1206  {
1207  if (cg.overlap[0].dataBegin()->origin(i) > 0 || _periodic[i])
1208  ovlp_low[i] = true;
1209  if (cg.overlap[0].dataBegin()->max(i) + 1 < globalSize(i) || _periodic[i])
1210  ovlp_up[i] = true;
1211  }
1212 
1213  Coordinates newcont(cg.coords.refine(ovlp_low, ovlp_up, cg.overlapSize, keep_ovlp));
1214 
1215  int overlap = (keep_ovlp) ? 2*cg.overlapSize : cg.overlapSize;
1216 
1217  //determine new origin
1218  iTupel o_interior;
1219  for (int i=0; i<dim; i++)
1220  o_interior[i] = 2*cg.interior[0].dataBegin()->origin(i);
1221 
1222  // add level
1223  _levels.resize(_levels.size() + 1);
1224  makelevel(newcont,_periodic,o_interior,overlap);
1225 
1226  indexsets.push_back( std::make_shared<YaspIndexSet<const YaspGrid<dim,Coordinates>, false > >(*this,maxLevel()) );
1227  }
1228  }
1229 
1234  void refineOptions (bool keepPhysicalOverlap)
1235  {
1236  keep_ovlp = keepPhysicalOverlap;
1237  }
1238 
1250  bool mark( int refCount, const typename Traits::template Codim<0>::Entity & e )
1251  {
1252  assert(adaptActive == false);
1253  if (e.level() != maxLevel()) return false;
1254  adaptRefCount = std::max(adaptRefCount, refCount);
1255  return true;
1256  }
1257 
1264  int getMark ( const typename Traits::template Codim<0>::Entity &e ) const
1265  {
1266  return ( e.level() == maxLevel() ) ? adaptRefCount : 0;
1267  }
1268 
1270  bool adapt ()
1271  {
1272  globalRefine(adaptRefCount);
1273  return (adaptRefCount > 0);
1274  }
1275 
1277  bool preAdapt ()
1278  {
1279  adaptActive = true;
1280  adaptRefCount = comm().max(adaptRefCount);
1281  return (adaptRefCount < 0);
1282  }
1283 
1285  void postAdapt()
1286  {
1287  adaptActive = false;
1288  adaptRefCount = 0;
1289  }
1290 
1292  template<int cd, PartitionIteratorType pitype>
1293  typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lbegin (int level) const
1294  {
1295  return levelbegin<cd,pitype>(level);
1296  }
1297 
1299  template<int cd, PartitionIteratorType pitype>
1300  typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lend (int level) const
1301  {
1302  return levelend<cd,pitype>(level);
1303  }
1304 
1306  template<int cd>
1307  typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lbegin (int level) const
1308  {
1309  return levelbegin<cd,All_Partition>(level);
1310  }
1311 
1313  template<int cd>
1314  typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lend (int level) const
1315  {
1316  return levelend<cd,All_Partition>(level);
1317  }
1318 
1320  template<int cd, PartitionIteratorType pitype>
1322  {
1323  return levelbegin<cd,pitype>(maxLevel());
1324  }
1325 
1327  template<int cd, PartitionIteratorType pitype>
1329  {
1330  return levelend<cd,pitype>(maxLevel());
1331  }
1332 
1334  template<int cd>
1336  {
1337  return levelbegin<cd,All_Partition>(maxLevel());
1338  }
1339 
1341  template<int cd>
1343  {
1344  return levelend<cd,All_Partition>(maxLevel());
1345  }
1346 
1355  template <typename Seed>
1356  DUNE_DEPRECATED_MSG("entityPointer() is deprecated and will be removed after the release of dune-grid 2.4. Use entity() instead to directly obtain an Entity object.")
1357  typename Traits::template Codim<Seed::codimension>::EntityPointer
1358  entityPointer(const Seed& seed) const
1359  {
1360  const int codim = Seed::codimension;
1361  YGridLevelIterator g = begin(this->getRealImplementation(seed).level());
1362 
1364  typename YGrid::Iterator(g->overlapfront[codim], this->getRealImplementation(seed).coord(),this->getRealImplementation(seed).offset()));
1365  }
1366 
1367  // \brief obtain Entity from EntitySeed. */
1368  template <typename Seed>
1369  typename Traits::template Codim<Seed::codimension>::Entity
1370  entity(const Seed& seed) const
1371  {
1372  const int codim = Seed::codimension;
1373  YGridLevelIterator g = begin(this->getRealImplementation(seed).level());
1374 
1375  typedef typename Traits::template Codim<Seed::codimension>::Entity Entity;
1376  typedef YaspEntity<codim,dim,const YaspGrid> EntityImp;
1377  typedef typename YGrid::Iterator YIterator;
1378 
1379  return Entity(EntityImp(g,YIterator(g->overlapfront[codim],this->getRealImplementation(seed).coord(),this->getRealImplementation(seed).offset())));
1380  }
1381 
1383  int overlapSize (int level, int codim) const
1384  {
1385  YGridLevelIterator g = begin(level);
1386  return g->overlapSize;
1387  }
1388 
1390  int overlapSize (int codim) const
1391  {
1392  YGridLevelIterator g = begin(maxLevel());
1393  return g->overlapSize;
1394  }
1395 
1397  int ghostSize (int level, int codim) const
1398  {
1399  return 0;
1400  }
1401 
1403  int ghostSize (int codim) const
1404  {
1405  return 0;
1406  }
1407 
1409  int size (int level, int codim) const
1410  {
1411  YGridLevelIterator g = begin(level);
1412 
1413  // sum over all components of the codimension
1414  int count = 0;
1415  typedef typename std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power>::iterator DAI;
1416  for (DAI it = g->overlapfront[codim].dataBegin(); it != g->overlapfront[codim].dataEnd(); ++it)
1417  count += it->totalsize();
1418 
1419  return count;
1420  }
1421 
1423  int size (int codim) const
1424  {
1425  return size(maxLevel(),codim);
1426  }
1427 
1429  int size (int level, GeometryType type) const
1430  {
1431  return (type.isCube()) ? size(level,dim-type.dim()) : 0;
1432  }
1433 
1435  int size (GeometryType type) const
1436  {
1437  return size(maxLevel(),type);
1438  }
1439 
1441  size_t numBoundarySegments () const
1442  {
1443  return nBSegments;
1444  }
1445 
1447  const Dune::FieldVector<ctype, dim>& domainSize () const {
1448  return _L;
1449  }
1450 
1455  template<class DataHandleImp, class DataType>
1457  {
1458  YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,level);
1459  }
1460 
1465  template<class DataHandleImp, class DataType>
1467  {
1468  YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,this->maxLevel());
1469  }
1470 
1475  template<class DataHandle, int codim>
1476  void communicateCodim (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level) const
1477  {
1478  // check input
1479  if (!data.contains(dim,codim)) return; // should have been checked outside
1480 
1481  // data types
1482  typedef typename DataHandle::DataType DataType;
1483 
1484  // access to grid level
1485  YGridLevelIterator g = begin(level);
1486 
1487  // find send/recv lists or throw error
1488  const YGridList<Coordinates>* sendlist = 0;
1489  const YGridList<Coordinates>* recvlist = 0;
1490 
1492  {
1493  sendlist = &g->send_interiorborder_interiorborder[codim];
1494  recvlist = &g->recv_interiorborder_interiorborder[codim];
1495  }
1496  if (iftype==InteriorBorder_All_Interface)
1497  {
1498  sendlist = &g->send_interiorborder_overlapfront[codim];
1499  recvlist = &g->recv_overlapfront_interiorborder[codim];
1500  }
1502  {
1503  sendlist = &g->send_overlap_overlapfront[codim];
1504  recvlist = &g->recv_overlapfront_overlap[codim];
1505  }
1506  if (iftype==All_All_Interface)
1507  {
1508  sendlist = &g->send_overlapfront_overlapfront[codim];
1509  recvlist = &g->recv_overlapfront_overlapfront[codim];
1510  }
1511 
1512  // change communication direction?
1513  if (dir==BackwardCommunication)
1514  std::swap(sendlist,recvlist);
1515 
1516  int cnt;
1517 
1518  // Size computation (requires communication if variable size)
1519  std::vector<int> send_size(sendlist->size(),-1); // map rank to total number of objects (of type DataType) to be sent
1520  std::vector<int> recv_size(recvlist->size(),-1); // map rank to total number of objects (of type DataType) to be recvd
1521  std::vector<size_t*> send_sizes(sendlist->size(),static_cast<size_t*>(0)); // map rank to array giving number of objects per entity to be sent
1522  std::vector<size_t*> recv_sizes(recvlist->size(),static_cast<size_t*>(0)); // map rank to array giving number of objects per entity to be recvd
1523 
1524  // define type to iterate over send and recv lists
1525  typedef typename YGridList<Coordinates>::Iterator ListIt;
1526 
1527  if (data.fixedsize(dim,codim))
1528  {
1529  // fixed size: just take a dummy entity, size can be computed without communication
1530  cnt=0;
1531  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1532  {
1535  send_size[cnt] = is->grid.totalsize() * data.size(*it);
1536  cnt++;
1537  }
1538  cnt=0;
1539  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1540  {
1543  recv_size[cnt] = is->grid.totalsize() * data.size(*it);
1544  cnt++;
1545  }
1546  }
1547  else
1548  {
1549  // variable size case: sender side determines the size
1550  cnt=0;
1551  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1552  {
1553  // allocate send buffer for sizes per entitiy
1554  size_t *buf = new size_t[is->grid.totalsize()];
1555  send_sizes[cnt] = buf;
1556 
1557  // loop over entities and ask for size
1558  int i=0; size_t n=0;
1562  itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1563  for ( ; it!=itend; ++it)
1564  {
1565  buf[i] = data.size(*it);
1566  n += buf[i];
1567  i++;
1568  }
1569 
1570  // now we know the size for this rank
1571  send_size[cnt] = n;
1572 
1573  // hand over send request to torus class
1574  torus().send(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
1575  cnt++;
1576  }
1577 
1578  // allocate recv buffers for sizes and store receive request
1579  cnt=0;
1580  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1581  {
1582  // allocate recv buffer
1583  size_t *buf = new size_t[is->grid.totalsize()];
1584  recv_sizes[cnt] = buf;
1585 
1586  // hand over recv request to torus class
1587  torus().recv(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
1588  cnt++;
1589  }
1590 
1591  // exchange all size buffers now
1592  torus().exchange();
1593 
1594  // release send size buffers
1595  cnt=0;
1596  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1597  {
1598  delete[] send_sizes[cnt];
1599  send_sizes[cnt] = 0;
1600  cnt++;
1601  }
1602 
1603  // process receive size buffers
1604  cnt=0;
1605  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1606  {
1607  // get recv buffer
1608  size_t *buf = recv_sizes[cnt];
1609 
1610  // compute total size
1611  size_t n=0;
1612  for (int i=0; i<is->grid.totalsize(); ++i)
1613  n += buf[i];
1614 
1615  // ... and store it
1616  recv_size[cnt] = n;
1617  ++cnt;
1618  }
1619  }
1620 
1621 
1622  // allocate & fill the send buffers & store send request
1623  std::vector<DataType*> sends(sendlist->size(), static_cast<DataType*>(0)); // store pointers to send buffers
1624  cnt=0;
1625  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1626  {
1627  // allocate send buffer
1628  DataType *buf = new DataType[send_size[cnt]];
1629 
1630  // remember send buffer
1631  sends[cnt] = buf;
1632 
1633  // make a message buffer
1634  MessageBuffer<DataType> mb(buf);
1635 
1636  // fill send buffer; iterate over cells in intersection
1640  itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1641  for ( ; it!=itend; ++it)
1642  data.gather(mb,*it);
1643 
1644  // hand over send request to torus class
1645  torus().send(is->rank,buf,send_size[cnt]*sizeof(DataType));
1646  cnt++;
1647  }
1648 
1649  // allocate recv buffers and store receive request
1650  std::vector<DataType*> recvs(recvlist->size(),static_cast<DataType*>(0)); // store pointers to send buffers
1651  cnt=0;
1652  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1653  {
1654  // allocate recv buffer
1655  DataType *buf = new DataType[recv_size[cnt]];
1656 
1657  // remember recv buffer
1658  recvs[cnt] = buf;
1659 
1660  // hand over recv request to torus class
1661  torus().recv(is->rank,buf,recv_size[cnt]*sizeof(DataType));
1662  cnt++;
1663  }
1664 
1665  // exchange all buffers now
1666  torus().exchange();
1667 
1668  // release send buffers
1669  cnt=0;
1670  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1671  {
1672  delete[] sends[cnt];
1673  sends[cnt] = 0;
1674  cnt++;
1675  }
1676 
1677  // process receive buffers and delete them
1678  cnt=0;
1679  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1680  {
1681  // get recv buffer
1682  DataType *buf = recvs[cnt];
1683 
1684  // make a message buffer
1685  MessageBuffer<DataType> mb(buf);
1686 
1687  // copy data from receive buffer; iterate over cells in intersection
1688  if (data.fixedsize(dim,codim))
1689  {
1692  size_t n=data.size(*it);
1694  itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1695  for ( ; it!=itend; ++it)
1696  data.scatter(mb,*it,n);
1697  }
1698  else
1699  {
1700  int i=0;
1701  size_t *sbuf = recv_sizes[cnt];
1705  itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1706  for ( ; it!=itend; ++it)
1707  data.scatter(mb,*it,sbuf[i++]);
1708  delete[] sbuf;
1709  }
1710 
1711  // delete buffer
1712  delete[] buf; // hier krachts !
1713  cnt++;
1714  }
1715  }
1716 
1717  // The new index sets from DDM 11.07.2005
1718  const typename Traits::GlobalIdSet& globalIdSet() const
1719  {
1720  return theglobalidset;
1721  }
1722 
1723  const typename Traits::LocalIdSet& localIdSet() const
1724  {
1725  return theglobalidset;
1726  }
1727 
1728  const typename Traits::LevelIndexSet& levelIndexSet(int level) const
1729  {
1730  if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1731  return *(indexsets[level]);
1732  }
1733 
1734  const typename Traits::LeafIndexSet& leafIndexSet() const
1735  {
1736  return leafIndexSet_;
1737  }
1738 
1741  const CollectiveCommunicationType& comm () const
1742  {
1743  return ccobj;
1744  }
1745 
1746  private:
1747 
1748  // number of boundary segments of the level 0 grid
1749  int nBSegments;
1750 
1751  // Index classes need access to the real entity
1752  friend class Dune::YaspIndexSet<const Dune::YaspGrid<dim, Coordinates>, true >;
1753  friend class Dune::YaspIndexSet<const Dune::YaspGrid<dim, Coordinates>, false >;
1754  friend class Dune::YaspGlobalIdSet<const Dune::YaspGrid<dim, Coordinates> >;
1755  friend class Dune::YaspPersistentContainerIndex<const Dune::YaspGrid<dim, Coordinates> >;
1756 
1757  friend class Dune::YaspIntersectionIterator<const Dune::YaspGrid<dim, Coordinates> >;
1758  friend class Dune::YaspIntersection<const Dune::YaspGrid<dim, Coordinates> >;
1759  friend class Dune::YaspEntity<0, dim, const Dune::YaspGrid<dim, Coordinates> >;
1760 
1761  template <int codim_, class GridImp_>
1763 
1764  template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
1765  friend class Entity;
1766 
1767  template<class DT>
1768  class MessageBuffer {
1769  public:
1770  // Constructor
1771  MessageBuffer (DT *p)
1772  {
1773  a=p;
1774  i=0;
1775  j=0;
1776  }
1777 
1778  // write data to message buffer, acts like a stream !
1779  template<class Y>
1780  void write (const Y& data)
1781  {
1782  static_assert(( is_same<DT,Y>::value ), "DataType mismatch");
1783  a[i++] = data;
1784  }
1785 
1786  // read data from message buffer, acts like a stream !
1787  template<class Y>
1788  void read (Y& data) const
1789  {
1790  static_assert(( is_same<DT,Y>::value ), "DataType mismatch");
1791  data = a[j++];
1792  }
1793 
1794  private:
1795  DT *a;
1796  int i;
1797  mutable int j;
1798  };
1799 
1801  template<int cd, PartitionIteratorType pitype>
1802  YaspLevelIterator<cd,pitype,GridImp> levelbegin (int level) const
1803  {
1804  YGridLevelIterator g = begin(level);
1805  if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1806 
1807  if (pitype==Interior_Partition)
1808  return YaspLevelIterator<cd,pitype,GridImp>(g,g->interior[cd].begin());
1809  if (pitype==InteriorBorder_Partition)
1810  return YaspLevelIterator<cd,pitype,GridImp>(g,g->interiorborder[cd].begin());
1811  if (pitype==Overlap_Partition)
1812  return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlap[cd].begin());
1813  if (pitype<=All_Partition)
1814  return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlapfront[cd].begin());
1815  if (pitype==Ghost_Partition)
1816  return levelend <cd, pitype> (level);
1817 
1818  DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
1819  }
1820 
1822  template<int cd, PartitionIteratorType pitype>
1823  YaspLevelIterator<cd,pitype,GridImp> levelend (int level) const
1824  {
1825  YGridLevelIterator g = begin(level);
1826  if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1827 
1828  if (pitype==Interior_Partition)
1829  return YaspLevelIterator<cd,pitype,GridImp>(g,g->interior[cd].end());
1830  if (pitype==InteriorBorder_Partition)
1831  return YaspLevelIterator<cd,pitype,GridImp>(g,g->interiorborder[cd].end());
1832  if (pitype==Overlap_Partition)
1833  return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlap[cd].end());
1834  if (pitype<=All_Partition || pitype == Ghost_Partition)
1835  return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlapfront[cd].end());
1836 
1837  DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
1838  }
1839 
1840  CollectiveCommunicationType ccobj;
1841 
1842  Torus<CollectiveCommunicationType,dim> _torus;
1843 
1844  std::vector< std::shared_ptr< YaspIndexSet<const YaspGrid<dim,Coordinates>, false > > > indexsets;
1845  YaspIndexSet<const YaspGrid<dim,Coordinates>, true> leafIndexSet_;
1846  YaspGlobalIdSet<const YaspGrid<dim,Coordinates> > theglobalidset;
1847 
1848  Dune::FieldVector<ctype, dim> _L;
1849  iTupel _s;
1850  std::bitset<dim> _periodic;
1851  iTupel _coarseSize;
1852  ReservedVector<YGridLevel,32> _levels;
1853  int _overlap;
1854  bool keep_ovlp;
1855  int adaptRefCount;
1856  bool adaptActive;
1857  };
1858 
1860 
1861  template <int d, class CC>
1862  std::ostream& operator<< (std::ostream& s, const YaspGrid<d,CC>& grid)
1863  {
1864  int rank = grid.torus().rank();
1865 
1866  s << "[" << rank << "]:" << " YaspGrid maxlevel=" << grid.maxLevel() << std::endl;
1867 
1868  s << "Printing the torus: " <<std::endl;
1869  s << grid.torus() << std::endl;
1870 
1871  for (typename YaspGrid<d,CC>::YGridLevelIterator g=grid.begin(); g!=grid.end(); ++g)
1872  {
1873  s << "[" << rank << "]: " << std::endl;
1874  s << "[" << rank << "]: " << "==========================================" << std::endl;
1875  s << "[" << rank << "]: " << "level=" << g->level() << std::endl;
1876 
1877  for (int codim = 0; codim < d + 1; ++codim)
1878  {
1879  s << "[" << rank << "]: " << "overlapfront[" << codim << "]: " << g->overlapfront[codim] << std::endl;
1880  s << "[" << rank << "]: " << "overlap[" << codim << "]: " << g->overlap[codim] << std::endl;
1881  s << "[" << rank << "]: " << "interiorborder[" << codim << "]: " << g->interiorborder[codim] << std::endl;
1882  s << "[" << rank << "]: " << "interior[" << codim << "]: " << g->interior[codim] << std::endl;
1883 
1884  typedef typename YGridList<CC>::Iterator I;
1885  for (I i=g->send_overlapfront_overlapfront[codim].begin();
1886  i!=g->send_overlapfront_overlapfront[codim].end(); ++i)
1887  s << "[" << rank << "]: " << " s_of_of[" << codim << "] to rank "
1888  << i->rank << " " << i->grid << std::endl;
1889 
1890  for (I i=g->recv_overlapfront_overlapfront[codim].begin();
1891  i!=g->recv_overlapfront_overlapfront[codim].end(); ++i)
1892  s << "[" << rank << "]: " << " r_of_of[" << codim << "] to rank "
1893  << i->rank << " " << i->grid << std::endl;
1894 
1895  for (I i=g->send_overlap_overlapfront[codim].begin();
1896  i!=g->send_overlap_overlapfront[codim].end(); ++i)
1897  s << "[" << rank << "]: " << " s_o_of[" << codim << "] to rank "
1898  << i->rank << " " << i->grid << std::endl;
1899 
1900  for (I i=g->recv_overlapfront_overlap[codim].begin();
1901  i!=g->recv_overlapfront_overlap[codim].end(); ++i)
1902  s << "[" << rank << "]: " << " r_of_o[" << codim << "] to rank "
1903  << i->rank << " " << i->grid << std::endl;
1904 
1905  for (I i=g->send_interiorborder_interiorborder[codim].begin();
1906  i!=g->send_interiorborder_interiorborder[codim].end(); ++i)
1907  s << "[" << rank << "]: " << " s_ib_ib[" << codim << "] to rank "
1908  << i->rank << " " << i->grid << std::endl;
1909 
1910  for (I i=g->recv_interiorborder_interiorborder[codim].begin();
1911  i!=g->recv_interiorborder_interiorborder[codim].end(); ++i)
1912  s << "[" << rank << "]: " << " r_ib_ib[" << codim << "] to rank "
1913  << i->rank << " " << i->grid << std::endl;
1914 
1915  for (I i=g->send_interiorborder_overlapfront[codim].begin();
1916  i!=g->send_interiorborder_overlapfront[codim].end(); ++i)
1917  s << "[" << rank << "]: " << " s_ib_of[" << codim << "] to rank "
1918  << i->rank << " " << i->grid << std::endl;
1919 
1920  for (I i=g->recv_overlapfront_interiorborder[codim].begin();
1921  i!=g->recv_overlapfront_interiorborder[codim].end(); ++i)
1922  s << "[" << rank << "]: " << " r_of_ib[" << codim << "] to rank "
1923  << i->rank << " " << i->grid << std::endl;
1924  }
1925  }
1926 
1927  s << std::endl;
1928 
1929  return s;
1930  }
1931 
1932  namespace Capabilities
1933  {
1934 
1942  template<int dim, class Coordinates>
1943  struct hasBackupRestoreFacilities< YaspGrid<dim, Coordinates> >
1944  {
1945  static const bool v = true;
1946  };
1947 
1951  template<int dim, class Coordinates>
1952  struct hasSingleGeometryType< YaspGrid<dim, Coordinates> >
1953  {
1954  static const bool v = true;
1955  static const unsigned int topologyId = GenericGeometry :: CubeTopology< dim > :: type :: id ;
1956  };
1957 
1961  template<int dim, class Coordinates>
1962  struct isCartesian< YaspGrid<dim, Coordinates> >
1963  {
1964  static const bool v = true;
1965  };
1966 
1970  template<int dim, class Coordinates, int codim>
1971  struct hasEntity< YaspGrid<dim, Coordinates>, codim>
1972  {
1973  static const bool v = true;
1974  };
1975 
1979  template<int dim, int codim, class Coordinates>
1980  struct canCommunicate< YaspGrid< dim, Coordinates>, codim >
1981  {
1982  static const bool v = true;
1983  };
1984 
1989  template<int dim, class Coordinates>
1990  struct DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") isParallel< YaspGrid<dim, Coordinates> >
1991  {
1992  static const bool DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") v = true;
1993  };
1994 
1998  template<int dim, class Coordinates>
1999  struct isLevelwiseConforming< YaspGrid<dim, Coordinates> >
2000  {
2001  static const bool v = true;
2002  };
2003 
2007  template<int dim, class Coordinates>
2008  struct isLeafwiseConforming< YaspGrid<dim, Coordinates> >
2009  {
2010  static const bool v = true;
2011  };
2012 
2013  }
2014 
2015 } // end namespace
2016 
2017 
2018 #endif
type describing an intersection with a neighboring processor
Definition: ygrid.hh:826
void communicateCodim(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int level) const
Definition: yaspgrid.hh:1476
CollectiveCommunication< MPI_Comm > CCType
Definition: yaspgrid.hh:92
int size() const
return the size of the container, this is the sum of the sizes of all deques
Definition: ygrid.hh:951
bool adapt()
map adapt to global refine
Definition: yaspgrid.hh:1270
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lbegin(int level) const
version without second template parameter for convenience
Definition: yaspgrid.hh:1307
YGridLevelIterator end() const
return iterator pointing to one past the finest level
Definition: yaspgrid.hh:311
YaspGridFamily< dim, Coordinates > GridFamily
the GridFamily of this grid
Definition: yaspgrid.hh:721
Definition: defaultgridview.hh:223
Wrapper class for pointers to entities.
Definition: common/entitypointer.hh:112
Definition: common/geometry.hh:24
int overlapSize(int level, int codim) const
return size (= distance in graph) of overlap region
Definition: yaspgrid.hh:1383
a base class for the yaspgrid partitioning strategy The name might be irritating. It will probably ch...
Definition: partitioning.hh:23
YGridComponent< Coordinates > move(iTupel v) const
return grid moved by the vector v
Definition: ygrid.hh:260
interior, border, and overlap entities
Definition: gridenums.hh:137
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity in maxLevel
Definition: yaspgrid.hh:1321
const YaspGrid< dim, Coordinates > GridImp
Definition: yaspgrid.hh:678
send interior and border, receive all entities
Definition: gridenums.hh:86
bool getRefineOption() const
Definition: yaspgrid.hh:288
[ provides Dune::Grid ]
Definition: yaspgrid.hh:56
int rank() const
return own rank
Definition: torus.hh:96
A Traits struct that collects all associated types of one implementation.
Definition: common/grid.hh:437
int neighbors() const
return the number of neighbors, which is
Definition: torus.hh:205
InterfaceType
Parameter to be used for the communication functions.
Definition: gridenums.hh:84
The YaspIntersection class.
unsigned char uint8_t
Definition: yaspgrid.hh:15
size_t numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition: yaspgrid.hh:1441
bool preAdapt()
returns true, if the grid will be coarsened
Definition: yaspgrid.hh:1277
int levelSize(int l, int i) const
return size of the grid (in cells) on level l in direction i
Definition: yaspgrid.hh:268
bool checkIfMonotonous(const Dune::array< std::vector< ctype >, dim > &coords)
Definition: coordinates.hh:361
YaspHierarchicIterator enables iteration over son entities of codim 0.
Definition: yaspgrid.hh:64
Interior interior
PartitionSet for the interior partition.
Definition: partitionset.hh:226
Dune::array< int, d > sizeArray(const Dune::array< std::vector< ct >, d > &v)
Definition: ygrid.hh:26
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:178
Specialize with 'true' if the grid is a Cartesian grid. Cartesian grids satisfy the following propert...
Definition: common/capabilities.hh:47
void exchange() const
exchange messages stored in request buffers; clear request buffers afterwards
Definition: torus.hh:389
Specialization of the PersistentContainer for YaspGrid.
Definition: yaspgrid.hh:67
The general version that handles all codimensions but 0 and dim.
Definition: yaspgrid.hh:57
int max(const DofVectorPointer< int > &dofVector)
Definition: dofvector.hh:335
Include standard header files.
Definition: agrid.hh:59
const Torus< CollectiveCommunicationType, dim > & torus() const
return reference to torus
Definition: yaspgrid.hh:250
Describes the minimal information necessary to create a fully functional YaspEntity.
Definition: yaspgrid.hh:60
void init()
Definition: yaspgrid.hh:680
Coordinates::ctype ctype
Type used for coordinates.
Definition: yaspgrid.hh:180
Coordinate container for a tensor product YaspGrid.
Definition: coordinates.hh:233
YaspGrid(std::array< std::vector< ctype >, dim > coords, std::bitset< dim > periodic=std::bitset< dim >(0ULL), int overlap=1, CollectiveCommunicationType comm=CollectiveCommunicationType(), const YLoadBalance< dim > *lb=defaultLoadbalancer())
Standard constructor for a tensorproduct YaspGrid.
Definition: yaspgrid.hh:860
The YaspEntityPointer class.
send/receive interior and border entities
Definition: gridenums.hh:85
Definition: defaultgridview.hh:23
static const YLoadBalanceDefault< dim > * defaultLoadbalancer()
Definition: yaspgrid.hh:317
Front front
PartitionSet for the front partition.
Definition: partitionset.hh:235
ReservedVector< YGridLevel, 32 >::const_iterator YGridLevelIterator
Iterator over the grid levels.
Definition: yaspgrid.hh:294
send overlap, receive all entities
Definition: gridenums.hh:88
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: yaspgrid.hh:1429
send all and receive all entities
Definition: gridenums.hh:89
int maxLevel() const
Definition: yaspgrid.hh:1174
level-wise, non-persistent, consecutive indices for YaspGrid
int overlapSize(int codim) const
return size (= distance in graph) of overlap region
Definition: yaspgrid.hh:1390
friend class Entity
Definition: yaspgrid.hh:1765
Traits::template Codim< Seed::codimension >::EntityPointer entityPointer(const Seed &seed) const
obtain EntityPointer from EntitySeed.
Definition: yaspgrid.hh:1358
send overlap, receive overlap and front entities
Definition: gridenums.hh:87
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:1250
Iterator over a collection o YGrids A YGrid::Iterator is the heart of an entity in YaspGrid...
Definition: ygrid.hh:590
bool isPeriodic(int i) const
return whether the grid is periodic in direction i
Definition: yaspgrid.hh:283
Container for equidistant coordinates in a YaspGrid with non-trivial origin.
Definition: coordinates.hh:124
ProcListIterator sendend() const
end of send list
Definition: torus.hh:345
facility for writing and reading grids
Definition: common/backuprestore.hh:40
Specialize with 'true' for all codims that a grid implements entities for. (default=false) ...
Definition: common/capabilities.hh:57
only ghost entities
Definition: gridenums.hh:140
ProcListIterator sendbegin() const
first process in send list
Definition: torus.hh:339
void communicate(CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir) const
Definition: yaspgrid.hh:1466
int ghostSize(int codim) const
return size (= distance in graph) of ghost region
Definition: yaspgrid.hh:1403
implements a collection of multiple std::deque Intersections with neighboring processor...
Definition: ygrid.hh:820
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: yaspgrid.hh:1435
CommDataHandleIF describes the features of a data handle for communication in parallel runs using the...
Definition: datahandleif.hh:72
CollectiveCommunication< MPI_Comm > CollectiveCommunicationType
Definition: yaspgrid.hh:182
const int yaspgrid_level_bits
Definition: yaspgrid.hh:50
persistent, globally unique Ids
Definition: yaspgrid.hh:66
Different resources needed by all grid implementations.
static std::conditional< std::is_reference< InterfaceType >::value, typename std::add_lvalue_reference< typename ReturnImplementationType< typename std::remove_reference< InterfaceType >::type >::ImplementationType >::type, typename std::remove_const< typename ReturnImplementationType< typename std::remove_reference< InterfaceType >::type >::ImplementationType >::type >::type getRealImplementation(InterfaceType &&i)
return real implementation of interface class
Definition: common/grid.hh:1305
ProcListIterator recvbegin() const
first process in receive list
Definition: torus.hh:351
const Dune::FieldVector< ctype, dim > & domainSize() const
returns the size of the physical domain
Definition: yaspgrid.hh:1447
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity in maxLevel
Definition: yaspgrid.hh:1335
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:564
const Traits::LocalIdSet & localIdSet() const
Definition: yaspgrid.hh:1723
const Traits::LeafIndexSet & leafIndexSet() const
Definition: yaspgrid.hh:1734
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafend() const
return LeafIterator which points behind the last entity in maxLevel
Definition: yaspgrid.hh:1342
all entities
Definition: gridenums.hh:139
Definition: yaspgrid.hh:58
int getMark(const typename Traits::template Codim< 0 >::Entity &e) const
returns adaptation mark for given entity
Definition: yaspgrid.hh:1264
void globalRefine(int refCount)
refine the grid refCount times.
Definition: yaspgrid.hh:1180
Iterator begin() const
return iterator pointing to the begin of the container
Definition: ygrid.hh:921
void postAdapt()
clean up some markers
Definition: yaspgrid.hh:1285
The YaspLevelIterator class.
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:376
This file provides the infrastructure for toroidal communication in YaspGrid.
Iterates over entities of one grid level.
Definition: yaspgrid.hh:61
int ghostSize(int level, int codim) const
return size (= distance in graph) of ghost region
Definition: yaspgrid.hh:1397
YaspGridFamily< dim, Coordinates >::Traits Traits
Definition: yaspgrid.hh:723
int size(int codim) const
number of leaf entities per codim in this process
Definition: yaspgrid.hh:1423
Describes the parallel communication interface class for MessageBuffers and DataHandles.
YaspIntersection provides data about intersection with neighboring codim 0 entities.
Definition: yaspgrid.hh:63
int globalSize(int i) const
return number of cells on finest level in given direction on all processors
Definition: yaspgrid.hh:256
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:241
YaspIndexSet< YaspGrid< dim, Coordinates >, true > LeafIndexSetType
Definition: yaspgrid.hh:727
A pointer to a YaspGrid::Entity.
Definition: yaspgrid.hh:59
reverse communication direction
Definition: gridenums.hh:170
Id Set Interface.
Definition: common/grid.hh:362
static const bool v
Definition: common/capabilities.hh:50
void boundarysegmentssize()
Definition: yaspgrid.hh:689
void makelevel(const Coordinates &coords, std::bitset< dim > periodic, iTupel o_interior, int overlap)
Make a new YGridLevel structure.
Definition: yaspgrid.hh:331
Traits::template Codim< Seed::codimension >::Entity entity(const Seed &seed) const
Definition: yaspgrid.hh:1370
Overlap overlap
PartitionSet for the overlap partition.
Definition: partitionset.hh:232
Definition: ygrid.hh:71
static const bool v
Definition: common/capabilities.hh:118
Specialize with 'true' if implementation guarantees conforming level grids. (default=false) ...
Definition: common/capabilities.hh:98
Intersection of a mesh entities of codimension 0 ("elements") with a "neighboring" element or with th...
Definition: albertagrid/dgfparser.hh:26
const iTupel & dims() const
return dimensions of torus
Definition: torus.hh:114
STL namespace.
the YaspEntity class and its specializations
Index Set Interface base class.
Definition: common/grid.hh:361
const int yaspgrid_dim_bits
Definition: yaspgrid.hh:49
void refineOptions(bool keepPhysicalOverlap)
set options for refinement
Definition: yaspgrid.hh:1234
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:16
Implement the default load balance strategy of yaspgrid.
Definition: partitioning.hh:34
YaspIndexSet< YaspGrid< dim, Coordinates >, false > LevelIndexSetType
Definition: yaspgrid.hh:726
The YaspGeometry class and its specializations.
const Traits::GlobalIdSet & globalIdSet() const
Definition: yaspgrid.hh:1718
Specialize with 'true' for if the codimension 0 entity of the grid has only one possible geometry typ...
Definition: common/capabilities.hh:26
Specialize with 'true' if implementation provides backup and restore facilities. (default=false) ...
Definition: common/capabilities.hh:116
Specialize with 'true' if implementation guarantees a conforming leaf grid. (default=false) ...
Definition: common/capabilities.hh:107
const Traits::LevelIndexSet & levelIndexSet(int level) const
Definition: yaspgrid.hh:1728
Container for equidistant coordinates in a YaspGrid.
Definition: coordinates.hh:26
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lend(int level) const
version without second template parameter for convenience
Definition: yaspgrid.hh:1314
The YaspIntersectionIterator class.
ProcListIterator recvend() const
last process in receive list
Definition: torus.hh:357
static const bool v
Definition: common/capabilities.hh:91
Provides base classes for index and id sets.
YaspIntersectionIterator enables iteration over intersections with neighboring codim 0 entities...
Definition: yaspgrid.hh:62
static const bool v
Definition: common/capabilities.hh:28
interior and border entities
Definition: gridenums.hh:136
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafend() const
return LeafIterator which points behind the last entity in maxLevel
Definition: yaspgrid.hh:1328
static const unsigned int topologyId
Definition: common/capabilities.hh:31
int size(int level, int codim) const
number of entities per level and codim in this process
Definition: yaspgrid.hh:1409
GridTraits< dim, dim, Dune::YaspGrid< dim, Coordinates >, YaspGeometry, YaspEntity, YaspEntityPointer, 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 > Traits
Definition: yaspgrid.hh:118
YGridLevelIterator begin() const
return iterator pointing to coarsest level
Definition: yaspgrid.hh:297
static const bool v
Definition: common/capabilities.hh:59
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:1300
Implementation of Level- and LeafIndexSets for YaspGrid.
Definition: yaspgrid.hh:65
Definition: ygrid.hh:843
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
This provides a YGrid, the elemental component of the yaspgrid implementation.
void communicate(CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir, int level) const
Definition: yaspgrid.hh:1456
Types for GridView.
Definition: common/grid.hh:420
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator lbegin(int level) const
one past the end on this level
Definition: yaspgrid.hh:1293
CommunicationDirection
Define a type for communication direction parameter.
Definition: gridenums.hh:168
Specialize with 'true' if implementation supports parallelism. (default=false)
Definition: common/capabilities.hh:68
This provides container classes for the coordinates to be used in YaspGrid Upon implementation of the...
YaspGrid(Dune::FieldVector< ctype, dim > L, std::array< int, dim > s, std::bitset< dim > periodic=std::bitset< dim >(0ULL), int overlap=1, CollectiveCommunicationType comm=CollectiveCommunicationType(), const YLoadBalance< dim > *lb=defaultLoadbalancer())
Definition: yaspgrid.hh:738
iTupel coord() const
return own coordinates
Definition: torus.hh:102
Iterator end() const
return iterator pointing to the end of the container
Definition: ygrid.hh:927
only interior entities
Definition: gridenums.hh:135
YaspGrid(Dune::FieldVector< ctype, dim > lowerleft, Dune::FieldVector< ctype, dim > upperright, std::array< int, dim > s, std::bitset< dim > periodic=std::bitset< dim >(0ULL), int overlap=1, CollectiveCommunicationType comm=CollectiveCommunicationType(), const YLoadBalance< dim > *lb=defaultLoadbalancer())
Definition: yaspgrid.hh:797
bigunsignedint< dim *yaspgrid_dim_bits+yaspgrid_level_bits+dim > PersistentIndexType
Definition: yaspgrid.hh:718
Definition: yaspgrid.hh:89
const CollectiveCommunicationType & comm() const
return a collective communication object
Definition: yaspgrid.hh:1741
iTupel levelSize(int l) const
return size vector of the grid (in cells) on level l
Definition: yaspgrid.hh:274
A traits struct that collects all associated types of one grid model.
Definition: common/grid.hh:1343
iTupel globalSize() const
return number of cells on finest level on all processors
Definition: yaspgrid.hh:262
The YaspEntitySeed class.
YGridLevelIterator begin(int i) const
return iterator pointing to given level
Definition: yaspgrid.hh:303
specialize with 'true' for all codims that a grid can communicate data on (default=false) ...
Definition: common/capabilities.hh:89
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:363
A set of traits classes to store static information about grid implementation.
Wrapper class for entities.
Definition: common/entity.hh:61
YaspGlobalIdSet< YaspGrid< dim, Coordinates > > GlobalIdSetType
Definition: yaspgrid.hh:728
implements a collection of YGridComponents which form a codimension Entities of given codimension c n...
Definition: ygrid.hh:547