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