Dune Core Modules (2.6.0)

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 #include <type_traits>
11 
12 // either include stdint.h or provide fallback for uint8_t
13 #if HAVE_STDINT_H
14 #include <stdint.h>
15 #else
16 typedef unsigned char uint8_t;
17 #endif
18 
19 #include <dune/grid/common/backuprestore.hh>
20 #include <dune/grid/common/grid.hh> // the grid base classes
21 #include <dune/grid/common/capabilities.hh> // the capabilities
22 #include <dune/common/hybridutilities.hh>
23 #include <dune/common/power.hh>
31 #include <dune/geometry/type.hh>
34 
35 
36 #if HAVE_MPI
38 #endif
39 
47 namespace Dune {
48 
49  /* some sizes for building global ids
50  */
51  const int yaspgrid_dim_bits = 24; // bits for encoding each dimension
52  const int yaspgrid_level_bits = 5; // bits for encoding level number
53 
54 
55  //************************************************************************
56  // forward declaration of templates
57 
58  template<int dim, class Coordinates> class YaspGrid;
59  template<int mydim, int cdim, class GridImp> class YaspGeometry;
60  template<int codim, int dim, class GridImp> class YaspEntity;
61  template<int codim, class GridImp> class YaspEntitySeed;
62  template<int codim, PartitionIteratorType pitype, class GridImp> class YaspLevelIterator;
63  template<class GridImp> class YaspIntersectionIterator;
64  template<class GridImp> class YaspIntersection;
65  template<class GridImp> class YaspHierarchicIterator;
66  template<class GridImp, bool isLeafIndexSet> class YaspIndexSet;
67  template<class GridImp> class YaspGlobalIdSet;
68  template<class GridImp> class YaspPersistentContainerIndex;
69 
70 } // namespace Dune
71 
83 #include <dune/grid/yaspgrid/yaspgrididset.hh>
85 
86 namespace Dune {
87 
88  template<int dim, class Coordinates>
89  struct YaspGridFamily
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
100  YaspGeometry,YaspEntity,
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
106  YaspHierarchicIterator,
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
110  YaspGlobalIdSet<const YaspGrid<dim, Coordinates> >,
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,
115  DefaultLevelGridViewTraits, DefaultLeafGridViewTraits,
116  YaspEntitySeed>
117  Traits;
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>
171  friend class YaspHierarchicIterator;
172 
173  protected:
174 
176 
177  public:
179  typedef typename Coordinates::ctype ctype;
180 #if HAVE_MPI
182 #else
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 
294 
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 
311  {
312  return YGridLevelIterator(_levels,maxLevel()+1);
313  }
314 
315  // static method to create the default load balance strategy
316  static const YLoadBalanceDefault<dim>* defaultLoadbalancer()
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
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 
677  typedef const YaspGrid<dim,Coordinates> GridImp;
678 
679  void init()
680  {
681  indexsets.push_back( std::make_shared< YaspIndexSet<const YaspGrid<dim, Coordinates>, false > >(*this,0) );
682  boundarysegmentssize();
683  }
684 
685  void boundarysegmentssize()
686  {
687  // sizes of local macro grid
688  std::array<int, dim> sides;
689  {
690  for (int i=0; i<dim; i++)
691  {
692  sides[i] =
693  ((begin()->overlap[0].dataBegin()->origin(i) == 0)+
694  (begin()->overlap[0].dataBegin()->origin(i) + begin()->overlap[0].dataBegin()->size(i)
695  == levelSize(0,i)));
696  }
697  }
698  nBSegments = 0;
699  for (int k=0; k<dim; k++)
700  {
701  int offset = 1;
702  for (int l=0; l<dim; l++)
703  {
704  if (l==k) continue;
705  offset *= begin()->overlap[0].dataBegin()->size(l);
706  }
707  nBSegments += sides[k]*offset;
708  }
709  }
710 
711  public:
712 
713  // define the persistent index type
714  typedef bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim> PersistentIndexType;
715 
717  typedef YaspGridFamily<dim, Coordinates> GridFamily;
718  // the Traits
720 
721  // need for friend declarations in entity
725 
735  std::array<int, dim> s,
736  std::bitset<dim> periodic = std::bitset<dim>(0ULL),
737  int overlap = 1,
739  const YLoadBalance<dim>* lb = defaultLoadbalancer())
740  : ccobj(comm), _torus(comm,tag,s,lb), leafIndexSet_(*this),
741  _L(L), _periodic(periodic), _coarseSize(s), _overlap(overlap),
742  keep_ovlp(true), adaptRefCount(0), adaptActive(false)
743  {
744  // check whether YaspGrid has been given the correct template parameter
745  static_assert(std::is_same<Coordinates,EquidistantCoordinates<ctype,dim> >::value,
746  "YaspGrid coordinate container template parameter and given constructor values do not match!");
747 
748  _levels.resize(1);
749 
750  iTupel o;
751  std::fill(o.begin(), o.end(), 0);
752  iTupel o_interior(o);
753  iTupel s_interior(s);
754 
755  _torus.partition(_torus.rank(),o,s,o_interior,s_interior);
756 
757 #if HAVE_MPI
758  // check whether the grid is large enough to be overlapping
759  for (int i=0; i<dim; i++)
760  {
761  // find out whether the grid is too small to
762  int toosmall = (s_interior[i] / 2 <= overlap) && // interior is very small
763  (periodic[i] || (s_interior[i] != s[i])); // there is an overlap in that direction
764  // communicate the result to all those processes to have all processors error out if one process failed.
765  int global = 0;
766  MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR, comm);
767  if (global)
768  DUNE_THROW(Dune::GridError,"YaspGrid is too small to be overlapping");
769  }
770 #endif // #if HAVE_MPI
771 
772  fTupel h(L);
773  for (int i=0; i<dim; i++)
774  h[i] /= s[i];
775 
776  iTupel s_overlap(s_interior);
777  for (int i=0; i<dim; i++)
778  {
779  if ((o_interior[i] - overlap > 0) || (periodic[i]))
780  s_overlap[i] += overlap;
781  if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
782  s_overlap[i] += overlap;
783  }
784 
785  EquidistantCoordinates<ctype,dim> cc(h,s_overlap);
786 
787  // add level
788  makelevel(cc,periodic,o_interior,overlap);
789 
790  init();
791  }
792 
804  std::array<int, dim> s,
805  std::bitset<dim> periodic = std::bitset<dim>(0ULL),
806  int overlap = 1,
808  const YLoadBalance<dim>* lb = defaultLoadbalancer())
809  : ccobj(comm), _torus(comm,tag,s,lb), leafIndexSet_(*this),
810  _L(upperright - lowerleft),
811  _periodic(periodic), _coarseSize(s), _overlap(overlap),
812  keep_ovlp(true), adaptRefCount(0), adaptActive(false)
813  {
814  // check whether YaspGrid has been given the correct template parameter
815  static_assert(std::is_same<Coordinates,EquidistantOffsetCoordinates<ctype,dim> >::value,
816  "YaspGrid coordinate container template parameter and given constructor values do not match!");
817 
818  _levels.resize(1);
819 
820  iTupel o;
821  std::fill(o.begin(), o.end(), 0);
822  iTupel o_interior(o);
823  iTupel s_interior(s);
824 
825  _torus.partition(_torus.rank(),o,s,o_interior,s_interior);
826 
827 #if HAVE_MPI
828  // check whether the grid is large enough to be overlapping
829  for (int i=0; i<dim; i++)
830  {
831  // find out whether the grid is too small to
832  int toosmall = (s_interior[i] / 2 <= overlap) && // interior is very small
833  (periodic[i] || (s_interior[i] != s[i])); // there is an overlap in that direction
834  // communicate the result to all those processes to have all processors error out if one process failed.
835  int global = 0;
836  MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR, comm);
837  if (global)
838  DUNE_THROW(Dune::GridError,"YaspGrid is too small to be overlapping");
839  }
840 #endif // #if HAVE_MPI
841 
842  Dune::FieldVector<ctype,dim> extension(upperright);
844  for (int i=0; i<dim; i++)
845  {
846  extension[i] -= lowerleft[i];
847  h[i] = extension[i] / s[i];
848  }
849 
850  iTupel s_overlap(s_interior);
851  for (int i=0; i<dim; i++)
852  {
853  if ((o_interior[i] - overlap > 0) || (periodic[i]))
854  s_overlap[i] += overlap;
855  if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
856  s_overlap[i] += overlap;
857  }
858 
859  EquidistantOffsetCoordinates<ctype,dim> cc(lowerleft,h,s_overlap);
860 
861  // add level
862  makelevel(cc,periodic,o_interior,overlap);
863 
864  init();
865  }
866 
874  YaspGrid (std::array<std::vector<ctype>, dim> coords,
875  std::bitset<dim> periodic = std::bitset<dim>(0ULL),
876  int overlap = 1,
878  const YLoadBalance<dim>* lb = defaultLoadbalancer())
879  : ccobj(comm), _torus(comm,tag,Dune::Yasp::sizeArray<dim>(coords),lb),
880  leafIndexSet_(*this), _periodic(periodic), _overlap(overlap),
881  keep_ovlp(true), adaptRefCount(0), adaptActive(false)
882  {
883  if (!Dune::Yasp::checkIfMonotonous(coords))
884  DUNE_THROW(Dune::GridError,"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
885 
886  // check whether YaspGrid has been given the correct template parameter
887  static_assert(std::is_same<Coordinates,TensorProductCoordinates<ctype,dim> >::value,
888  "YaspGrid coordinate container template parameter and given constructor values do not match!");
889 
890  _levels.resize(1);
891 
892  //determine sizes of vector to correctly construct torus structure and store for later size requests
893  for (int i=0; i<dim; i++) {
894  _coarseSize[i] = coords[i].size() - 1;
895  _L[i] = coords[i][_coarseSize[i]] - coords[i][0];
896  }
897 
898  iTupel o;
899  std::fill(o.begin(), o.end(), 0);
900  iTupel o_interior(o);
901  iTupel s_interior(_coarseSize);
902 
903  _torus.partition(_torus.rank(),o,_coarseSize,o_interior,s_interior);
904 
905 #if HAVE_MPI
906  // check whether the grid is large enough to be overlapping
907  for (int i=0; i<dim; i++)
908  {
909  // find out whether the grid is too small to
910  int toosmall = (s_interior[i] / 2 <= overlap) && // interior is very small
911  (periodic[i] || (s_interior[i] != _coarseSize[i])); // there is an overlap in that direction
912  // communicate the result to all those processes to have all processors error out if one process failed.
913  int global = 0;
914  MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR, comm);
915  if (global)
916  DUNE_THROW(Dune::GridError,"YaspGrid is too small to be overlapping");
917  }
918 #endif // #if HAVE_MPI
919 
920 
921  std::array<std::vector<ctype>,dim> newcoords;
922  std::array<int, dim> offset(o_interior);
923 
924  // find the relevant part of the coords vector for this processor and copy it to newcoords
925  for (int i=0; i<dim; ++i)
926  {
927  //define iterators on coords that specify the coordinate range to be used
928  typename std::vector<ctype>::iterator begin = coords[i].begin() + o_interior[i];
929  typename std::vector<ctype>::iterator end = begin + s_interior[i] + 1;
930 
931  // check whether we are not at the physical boundary. In that case overlap is a simple
932  // extension of the coordinate range to be used
933  if (o_interior[i] - overlap > 0)
934  {
935  begin = begin - overlap;
936  offset[i] -= overlap;
937  }
938  if (o_interior[i] + s_interior[i] + overlap < _coarseSize[i])
939  end = end + overlap;
940 
941  //copy the selected part in the new coord vector
942  newcoords[i].resize(end-begin);
943  std::copy(begin, end, newcoords[i].begin());
944 
945  // check whether we are at the physical boundary and a have a periodic grid.
946  // In this case the coordinate vector has to be tweaked manually.
947  if ((periodic[i]) && (o_interior[i] + s_interior[i] + overlap >= _coarseSize[i]))
948  {
949  // we need to add the first <overlap> cells to the end of newcoords
950  typename std::vector<ctype>::iterator it = coords[i].begin();
951  for (int j=0; j<overlap; ++j)
952  newcoords[i].push_back(newcoords[i].back() - *it + *(++it));
953  }
954 
955  if ((periodic[i]) && (o_interior[i] - overlap <= 0))
956  {
957  offset[i] -= overlap;
958 
959  // we need to add the last <overlap> cells to the begin of newcoords
960  typename std::vector<ctype>::iterator it = coords[i].end() - 1;
961  for (int j=0; j<overlap; ++j)
962  newcoords[i].insert(newcoords[i].begin(), newcoords[i].front() - *it + *(--it));
963  }
964  }
965 
966  TensorProductCoordinates<ctype,dim> cc(newcoords, offset);
967 
968  // add level
969  makelevel(cc,periodic,o_interior,overlap);
970  init();
971  }
972 
973  private:
974 
989  YaspGrid (std::array<std::vector<ctype>, dim> coords,
990  std::bitset<dim> periodic,
991  int overlap,
992  CollectiveCommunicationType comm,
993  std::array<int,dim> coarseSize,
994  const YLoadBalance<dim>* lb = defaultLoadbalancer())
995  : ccobj(comm), _torus(comm,tag,coarseSize,lb), leafIndexSet_(*this),
996  _periodic(periodic), _coarseSize(coarseSize), _overlap(overlap),
997  keep_ovlp(true), adaptRefCount(0), adaptActive(false)
998  {
999  // check whether YaspGrid has been given the correct template parameter
1000  static_assert(std::is_same<Coordinates,TensorProductCoordinates<ctype,dim> >::value,
1001  "YaspGrid coordinate container template parameter and given constructor values do not match!");
1002 
1003  if (!Dune::Yasp::checkIfMonotonous(coords))
1004  DUNE_THROW(Dune::GridError,"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
1005 
1006  for (int i=0; i<dim; i++)
1007  _L[i] = coords[i][coords[i].size() - 1] - coords[i][0];
1008 
1009  _levels.resize(1);
1010 
1011  std::array<int,dim> o;
1012  std::fill(o.begin(), o.end(), 0);
1013  std::array<int,dim> o_interior(o);
1014  std::array<int,dim> s_interior(coarseSize);
1015 
1016  _torus.partition(_torus.rank(),o,coarseSize,o_interior,s_interior);
1017 
1018  // get offset by modifying o_interior according to overlap
1019  std::array<int,dim> offset(o_interior);
1020  for (int i=0; i<dim; i++)
1021  if ((periodic[i]) || (o_interior[i] > 0))
1022  offset[i] -= overlap;
1023 
1024  TensorProductCoordinates<ctype,dim> cc(coords, offset);
1025 
1026  // add level
1027  makelevel(cc,periodic,o_interior,overlap);
1028 
1029  init();
1030  }
1031 
1032  // the backup restore facility needs to be able to use above constructor
1033  friend struct BackupRestoreFacility<YaspGrid<dim,Coordinates> >;
1034 
1035  // do not copy this class
1036  YaspGrid(const YaspGrid&);
1037 
1038  public:
1039 
1043  int maxLevel() const
1044  {
1045  return _levels.size()-1;
1046  }
1047 
1049  void globalRefine (int refCount)
1050  {
1051  if (refCount < -maxLevel())
1052  DUNE_THROW(GridError, "Only " << maxLevel() << " levels left. " <<
1053  "Coarsening " << -refCount << " levels requested!");
1054 
1055  // If refCount is negative then coarsen the grid
1056  for (int k=refCount; k<0; k++)
1057  {
1058  // create an empty grid level
1059  YGridLevel empty;
1060  _levels.back() = empty;
1061  // reduce maxlevel
1062  _levels.pop_back();
1063 
1064  indexsets.pop_back();
1065  }
1066 
1067  // If refCount is positive refine the grid
1068  for (int k=0; k<refCount; k++)
1069  {
1070  // access to coarser grid level
1071  YGridLevel& cg = _levels[maxLevel()];
1072 
1073  std::bitset<dim> ovlp_low(0ULL), ovlp_up(0ULL);
1074  for (int i=0; i<dim; i++)
1075  {
1076  if (cg.overlap[0].dataBegin()->origin(i) > 0 || _periodic[i])
1077  ovlp_low[i] = true;
1078  if (cg.overlap[0].dataBegin()->max(i) + 1 < globalSize(i) || _periodic[i])
1079  ovlp_up[i] = true;
1080  }
1081 
1082  Coordinates newcont(cg.coords.refine(ovlp_low, ovlp_up, cg.overlapSize, keep_ovlp));
1083 
1084  int overlap = (keep_ovlp) ? 2*cg.overlapSize : cg.overlapSize;
1085 
1086  //determine new origin
1087  iTupel o_interior;
1088  for (int i=0; i<dim; i++)
1089  o_interior[i] = 2*cg.interior[0].dataBegin()->origin(i);
1090 
1091  // add level
1092  _levels.resize(_levels.size() + 1);
1093  makelevel(newcont,_periodic,o_interior,overlap);
1094 
1095  indexsets.push_back( std::make_shared<YaspIndexSet<const YaspGrid<dim,Coordinates>, false > >(*this,maxLevel()) );
1096  }
1097  }
1098 
1103  void refineOptions (bool keepPhysicalOverlap)
1104  {
1105  keep_ovlp = keepPhysicalOverlap;
1106  }
1107 
1119  bool mark( int refCount, const typename Traits::template Codim<0>::Entity & e )
1120  {
1121  assert(adaptActive == false);
1122  if (e.level() != maxLevel()) return false;
1123  adaptRefCount = std::max(adaptRefCount, refCount);
1124  return true;
1125  }
1126 
1133  int getMark ( const typename Traits::template Codim<0>::Entity &e ) const
1134  {
1135  return ( e.level() == maxLevel() ) ? adaptRefCount : 0;
1136  }
1137 
1139  bool adapt ()
1140  {
1141  globalRefine(adaptRefCount);
1142  return (adaptRefCount > 0);
1143  }
1144 
1146  bool preAdapt ()
1147  {
1148  adaptActive = true;
1149  adaptRefCount = comm().max(adaptRefCount);
1150  return (adaptRefCount < 0);
1151  }
1152 
1154  void postAdapt()
1155  {
1156  adaptActive = false;
1157  adaptRefCount = 0;
1158  }
1159 
1161  template<int cd, PartitionIteratorType pitype>
1162  typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lbegin (int level) const
1163  {
1164  return levelbegin<cd,pitype>(level);
1165  }
1166 
1168  template<int cd, PartitionIteratorType pitype>
1169  typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lend (int level) const
1170  {
1171  return levelend<cd,pitype>(level);
1172  }
1173 
1175  template<int cd>
1176  typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lbegin (int level) const
1177  {
1178  return levelbegin<cd,All_Partition>(level);
1179  }
1180 
1182  template<int cd>
1183  typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lend (int level) const
1184  {
1185  return levelend<cd,All_Partition>(level);
1186  }
1187 
1189  template<int cd, PartitionIteratorType pitype>
1190  typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafbegin () const
1191  {
1192  return levelbegin<cd,pitype>(maxLevel());
1193  }
1194 
1196  template<int cd, PartitionIteratorType pitype>
1197  typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafend () const
1198  {
1199  return levelend<cd,pitype>(maxLevel());
1200  }
1201 
1203  template<int cd>
1204  typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafbegin () const
1205  {
1206  return levelbegin<cd,All_Partition>(maxLevel());
1207  }
1208 
1210  template<int cd>
1211  typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafend () const
1212  {
1213  return levelend<cd,All_Partition>(maxLevel());
1214  }
1215 
1216  // \brief obtain Entity from EntitySeed. */
1217  template <typename Seed>
1218  typename Traits::template Codim<Seed::codimension>::Entity
1219  entity(const Seed& seed) const
1220  {
1221  const int codim = Seed::codimension;
1222  YGridLevelIterator g = begin(this->getRealImplementation(seed).level());
1223 
1224  typedef typename Traits::template Codim<Seed::codimension>::Entity Entity;
1225  typedef YaspEntity<codim,dim,const YaspGrid> EntityImp;
1226  typedef typename YGrid::Iterator YIterator;
1227 
1228  return Entity(EntityImp(g,YIterator(g->overlapfront[codim],this->getRealImplementation(seed).coord(),this->getRealImplementation(seed).offset())));
1229  }
1230 
1232  int overlapSize (int level, int codim) const
1233  {
1234  YGridLevelIterator g = begin(level);
1235  return g->overlapSize;
1236  }
1237 
1239  int overlapSize (int codim) const
1240  {
1242  return g->overlapSize;
1243  }
1244 
1246  int ghostSize (int level, int codim) const
1247  {
1248  return 0;
1249  }
1250 
1252  int ghostSize (int codim) const
1253  {
1254  return 0;
1255  }
1256 
1258  int size (int level, int codim) const
1259  {
1260  YGridLevelIterator g = begin(level);
1261 
1262  // sum over all components of the codimension
1263  int count = 0;
1264  typedef typename std::array<YGridComponent<Coordinates>, StaticPower<2,dim>::power>::iterator DAI;
1265  for (DAI it = g->overlapfront[codim].dataBegin(); it != g->overlapfront[codim].dataEnd(); ++it)
1266  count += it->totalsize();
1267 
1268  return count;
1269  }
1270 
1272  int size (int codim) const
1273  {
1274  return size(maxLevel(),codim);
1275  }
1276 
1278  int size (int level, GeometryType type) const
1279  {
1280  return (type.isCube()) ? size(level,dim-type.dim()) : 0;
1281  }
1282 
1284  int size (GeometryType type) const
1285  {
1286  return size(maxLevel(),type);
1287  }
1288 
1290  size_t numBoundarySegments () const
1291  {
1292  return nBSegments;
1293  }
1294 
1297  return _L;
1298  }
1299 
1304  template<class DataHandleImp, class DataType>
1306  {
1307  YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,level);
1308  }
1309 
1314  template<class DataHandleImp, class DataType>
1316  {
1317  YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,this->maxLevel());
1318  }
1319 
1324  template<class DataHandle, int codim>
1325  void communicateCodim (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level) const
1326  {
1327  // check input
1328  if (!data.contains(dim,codim)) return; // should have been checked outside
1329 
1330  // data types
1331  typedef typename DataHandle::DataType DataType;
1332 
1333  // access to grid level
1334  YGridLevelIterator g = begin(level);
1335 
1336  // find send/recv lists or throw error
1337  const YGridList<Coordinates>* sendlist = 0;
1338  const YGridList<Coordinates>* recvlist = 0;
1339 
1341  {
1342  sendlist = &g->send_interiorborder_interiorborder[codim];
1343  recvlist = &g->recv_interiorborder_interiorborder[codim];
1344  }
1345  if (iftype==InteriorBorder_All_Interface)
1346  {
1347  sendlist = &g->send_interiorborder_overlapfront[codim];
1348  recvlist = &g->recv_overlapfront_interiorborder[codim];
1349  }
1351  {
1352  sendlist = &g->send_overlap_overlapfront[codim];
1353  recvlist = &g->recv_overlapfront_overlap[codim];
1354  }
1355  if (iftype==All_All_Interface)
1356  {
1357  sendlist = &g->send_overlapfront_overlapfront[codim];
1358  recvlist = &g->recv_overlapfront_overlapfront[codim];
1359  }
1360 
1361  // change communication direction?
1362  if (dir==BackwardCommunication)
1363  std::swap(sendlist,recvlist);
1364 
1365  int cnt;
1366 
1367  // Size computation (requires communication if variable size)
1368  std::vector<int> send_size(sendlist->size(),-1); // map rank to total number of objects (of type DataType) to be sent
1369  std::vector<int> recv_size(recvlist->size(),-1); // map rank to total number of objects (of type DataType) to be recvd
1370  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
1371  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
1372 
1373  // define type to iterate over send and recv lists
1374  typedef typename YGridList<Coordinates>::Iterator ListIt;
1375 
1376  if (data.fixedSize(dim,codim))
1377  {
1378  // fixed size: just take a dummy entity, size can be computed without communication
1379  cnt=0;
1380  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1381  {
1382  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1384  send_size[cnt] = is->grid.totalsize() * data.size(*it);
1385  cnt++;
1386  }
1387  cnt=0;
1388  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1389  {
1390  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1392  recv_size[cnt] = is->grid.totalsize() * data.size(*it);
1393  cnt++;
1394  }
1395  }
1396  else
1397  {
1398  // variable size case: sender side determines the size
1399  cnt=0;
1400  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1401  {
1402  // allocate send buffer for sizes per entitiy
1403  size_t *buf = new size_t[is->grid.totalsize()];
1404  send_sizes[cnt] = buf;
1405 
1406  // loop over entities and ask for size
1407  int i=0; size_t n=0;
1408  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1410  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1411  itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1412  for ( ; it!=itend; ++it)
1413  {
1414  buf[i] = data.size(*it);
1415  n += buf[i];
1416  i++;
1417  }
1418 
1419  // now we know the size for this rank
1420  send_size[cnt] = n;
1421 
1422  // hand over send request to torus class
1423  torus().send(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
1424  cnt++;
1425  }
1426 
1427  // allocate recv buffers for sizes and store receive request
1428  cnt=0;
1429  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1430  {
1431  // allocate recv buffer
1432  size_t *buf = new size_t[is->grid.totalsize()];
1433  recv_sizes[cnt] = buf;
1434 
1435  // hand over recv request to torus class
1436  torus().recv(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
1437  cnt++;
1438  }
1439 
1440  // exchange all size buffers now
1441  torus().exchange();
1442 
1443  // release send size buffers
1444  cnt=0;
1445  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1446  {
1447  delete[] send_sizes[cnt];
1448  send_sizes[cnt] = 0;
1449  cnt++;
1450  }
1451 
1452  // process receive size buffers
1453  cnt=0;
1454  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1455  {
1456  // get recv buffer
1457  size_t *buf = recv_sizes[cnt];
1458 
1459  // compute total size
1460  size_t n=0;
1461  for (int i=0; i<is->grid.totalsize(); ++i)
1462  n += buf[i];
1463 
1464  // ... and store it
1465  recv_size[cnt] = n;
1466  ++cnt;
1467  }
1468  }
1469 
1470 
1471  // allocate & fill the send buffers & store send request
1472  std::vector<DataType*> sends(sendlist->size(), static_cast<DataType*>(0)); // store pointers to send buffers
1473  cnt=0;
1474  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1475  {
1476  // allocate send buffer
1477  DataType *buf = new DataType[send_size[cnt]];
1478 
1479  // remember send buffer
1480  sends[cnt] = buf;
1481 
1482  // make a message buffer
1483  MessageBuffer<DataType> mb(buf);
1484 
1485  // fill send buffer; iterate over cells in intersection
1486  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1488  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1489  itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1490  for ( ; it!=itend; ++it)
1491  data.gather(mb,*it);
1492 
1493  // hand over send request to torus class
1494  torus().send(is->rank,buf,send_size[cnt]*sizeof(DataType));
1495  cnt++;
1496  }
1497 
1498  // allocate recv buffers and store receive request
1499  std::vector<DataType*> recvs(recvlist->size(),static_cast<DataType*>(0)); // store pointers to send buffers
1500  cnt=0;
1501  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1502  {
1503  // allocate recv buffer
1504  DataType *buf = new DataType[recv_size[cnt]];
1505 
1506  // remember recv buffer
1507  recvs[cnt] = buf;
1508 
1509  // hand over recv request to torus class
1510  torus().recv(is->rank,buf,recv_size[cnt]*sizeof(DataType));
1511  cnt++;
1512  }
1513 
1514  // exchange all buffers now
1515  torus().exchange();
1516 
1517  // release send buffers
1518  cnt=0;
1519  for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1520  {
1521  delete[] sends[cnt];
1522  sends[cnt] = 0;
1523  cnt++;
1524  }
1525 
1526  // process receive buffers and delete them
1527  cnt=0;
1528  for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1529  {
1530  // get recv buffer
1531  DataType *buf = recvs[cnt];
1532 
1533  // make a message buffer
1534  MessageBuffer<DataType> mb(buf);
1535 
1536  // copy data from receive buffer; iterate over cells in intersection
1537  if (data.fixedSize(dim,codim))
1538  {
1539  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1541  size_t n=data.size(*it);
1542  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1543  itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1544  for ( ; it!=itend; ++it)
1545  data.scatter(mb,*it,n);
1546  }
1547  else
1548  {
1549  int i=0;
1550  size_t *sbuf = recv_sizes[cnt];
1551  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1553  typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1554  itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1555  for ( ; it!=itend; ++it)
1556  data.scatter(mb,*it,sbuf[i++]);
1557  delete[] sbuf;
1558  }
1559 
1560  // delete buffer
1561  delete[] buf; // hier krachts !
1562  cnt++;
1563  }
1564  }
1565 
1566  // The new index sets from DDM 11.07.2005
1567  const typename Traits::GlobalIdSet& globalIdSet() const
1568  {
1569  return theglobalidset;
1570  }
1571 
1572  const typename Traits::LocalIdSet& localIdSet() const
1573  {
1574  return theglobalidset;
1575  }
1576 
1577  const typename Traits::LevelIndexSet& levelIndexSet(int level) const
1578  {
1579  if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1580  return *(indexsets[level]);
1581  }
1582 
1583  const typename Traits::LeafIndexSet& leafIndexSet() const
1584  {
1585  return leafIndexSet_;
1586  }
1587 
1591  {
1592  return ccobj;
1593  }
1594 
1595  private:
1596 
1597  // number of boundary segments of the level 0 grid
1598  int nBSegments;
1599 
1600  // Index classes need access to the real entity
1601  friend class Dune::YaspIndexSet<const Dune::YaspGrid<dim, Coordinates>, true >;
1602  friend class Dune::YaspIndexSet<const Dune::YaspGrid<dim, Coordinates>, false >;
1603  friend class Dune::YaspGlobalIdSet<const Dune::YaspGrid<dim, Coordinates> >;
1604  friend class Dune::YaspPersistentContainerIndex<const Dune::YaspGrid<dim, Coordinates> >;
1605 
1606  friend class Dune::YaspIntersectionIterator<const Dune::YaspGrid<dim, Coordinates> >;
1607  friend class Dune::YaspIntersection<const Dune::YaspGrid<dim, Coordinates> >;
1608  friend class Dune::YaspEntity<0, dim, const Dune::YaspGrid<dim, Coordinates> >;
1609 
1610  template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
1611  friend class Entity;
1612 
1613  template<class DT>
1614  class MessageBuffer {
1615  public:
1616  // Constructor
1617  MessageBuffer (DT *p)
1618  {
1619  a=p;
1620  i=0;
1621  j=0;
1622  }
1623 
1624  // write data to message buffer, acts like a stream !
1625  template<class Y>
1626  void write (const Y& data)
1627  {
1628  static_assert(( std::is_same<DT,Y>::value ), "DataType mismatch");
1629  a[i++] = data;
1630  }
1631 
1632  // read data from message buffer, acts like a stream !
1633  template<class Y>
1634  void read (Y& data) const
1635  {
1636  static_assert(( std::is_same<DT,Y>::value ), "DataType mismatch");
1637  data = a[j++];
1638  }
1639 
1640  private:
1641  DT *a;
1642  int i;
1643  mutable int j;
1644  };
1645 
1647  template<int cd, PartitionIteratorType pitype>
1648  YaspLevelIterator<cd,pitype,GridImp> levelbegin (int level) const
1649  {
1650  YGridLevelIterator g = begin(level);
1651  if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1652 
1653  if (pitype==Interior_Partition)
1654  return YaspLevelIterator<cd,pitype,GridImp>(g,g->interior[cd].begin());
1655  if (pitype==InteriorBorder_Partition)
1656  return YaspLevelIterator<cd,pitype,GridImp>(g,g->interiorborder[cd].begin());
1657  if (pitype==Overlap_Partition)
1658  return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlap[cd].begin());
1659  if (pitype<=All_Partition)
1660  return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlapfront[cd].begin());
1661  if (pitype==Ghost_Partition)
1662  return levelend <cd, pitype> (level);
1663 
1664  DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
1665  }
1666 
1668  template<int cd, PartitionIteratorType pitype>
1669  YaspLevelIterator<cd,pitype,GridImp> levelend (int level) const
1670  {
1671  YGridLevelIterator g = begin(level);
1672  if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1673 
1674  if (pitype==Interior_Partition)
1675  return YaspLevelIterator<cd,pitype,GridImp>(g,g->interior[cd].end());
1676  if (pitype==InteriorBorder_Partition)
1677  return YaspLevelIterator<cd,pitype,GridImp>(g,g->interiorborder[cd].end());
1678  if (pitype==Overlap_Partition)
1679  return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlap[cd].end());
1680  if (pitype<=All_Partition || pitype == Ghost_Partition)
1681  return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlapfront[cd].end());
1682 
1683  DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
1684  }
1685 
1686  CollectiveCommunicationType ccobj;
1687 
1688  Torus<CollectiveCommunicationType,dim> _torus;
1689 
1690  std::vector< std::shared_ptr< YaspIndexSet<const YaspGrid<dim,Coordinates>, false > > > indexsets;
1691  YaspIndexSet<const YaspGrid<dim,Coordinates>, true> leafIndexSet_;
1692  YaspGlobalIdSet<const YaspGrid<dim,Coordinates> > theglobalidset;
1693 
1695  iTupel _s;
1696  std::bitset<dim> _periodic;
1697  iTupel _coarseSize;
1698  ReservedVector<YGridLevel,32> _levels;
1699  int _overlap;
1700  bool keep_ovlp;
1701  int adaptRefCount;
1702  bool adaptActive;
1703  };
1704 
1706 
1707  template <int d, class CC>
1708  std::ostream& operator<< (std::ostream& s, const YaspGrid<d,CC>& grid)
1709  {
1710  int rank = grid.torus().rank();
1711 
1712  s << "[" << rank << "]:" << " YaspGrid maxlevel=" << grid.maxLevel() << std::endl;
1713 
1714  s << "Printing the torus: " <<std::endl;
1715  s << grid.torus() << std::endl;
1716 
1717  for (typename YaspGrid<d,CC>::YGridLevelIterator g=grid.begin(); g!=grid.end(); ++g)
1718  {
1719  s << "[" << rank << "]: " << std::endl;
1720  s << "[" << rank << "]: " << "==========================================" << std::endl;
1721  s << "[" << rank << "]: " << "level=" << g->level() << std::endl;
1722 
1723  for (int codim = 0; codim < d + 1; ++codim)
1724  {
1725  s << "[" << rank << "]: " << "overlapfront[" << codim << "]: " << g->overlapfront[codim] << std::endl;
1726  s << "[" << rank << "]: " << "overlap[" << codim << "]: " << g->overlap[codim] << std::endl;
1727  s << "[" << rank << "]: " << "interiorborder[" << codim << "]: " << g->interiorborder[codim] << std::endl;
1728  s << "[" << rank << "]: " << "interior[" << codim << "]: " << g->interior[codim] << std::endl;
1729 
1730  typedef typename YGridList<CC>::Iterator I;
1731  for (I i=g->send_overlapfront_overlapfront[codim].begin();
1732  i!=g->send_overlapfront_overlapfront[codim].end(); ++i)
1733  s << "[" << rank << "]: " << " s_of_of[" << codim << "] to rank "
1734  << i->rank << " " << i->grid << std::endl;
1735 
1736  for (I i=g->recv_overlapfront_overlapfront[codim].begin();
1737  i!=g->recv_overlapfront_overlapfront[codim].end(); ++i)
1738  s << "[" << rank << "]: " << " r_of_of[" << codim << "] to rank "
1739  << i->rank << " " << i->grid << std::endl;
1740 
1741  for (I i=g->send_overlap_overlapfront[codim].begin();
1742  i!=g->send_overlap_overlapfront[codim].end(); ++i)
1743  s << "[" << rank << "]: " << " s_o_of[" << codim << "] to rank "
1744  << i->rank << " " << i->grid << std::endl;
1745 
1746  for (I i=g->recv_overlapfront_overlap[codim].begin();
1747  i!=g->recv_overlapfront_overlap[codim].end(); ++i)
1748  s << "[" << rank << "]: " << " r_of_o[" << codim << "] to rank "
1749  << i->rank << " " << i->grid << std::endl;
1750 
1751  for (I i=g->send_interiorborder_interiorborder[codim].begin();
1752  i!=g->send_interiorborder_interiorborder[codim].end(); ++i)
1753  s << "[" << rank << "]: " << " s_ib_ib[" << codim << "] to rank "
1754  << i->rank << " " << i->grid << std::endl;
1755 
1756  for (I i=g->recv_interiorborder_interiorborder[codim].begin();
1757  i!=g->recv_interiorborder_interiorborder[codim].end(); ++i)
1758  s << "[" << rank << "]: " << " r_ib_ib[" << codim << "] to rank "
1759  << i->rank << " " << i->grid << std::endl;
1760 
1761  for (I i=g->send_interiorborder_overlapfront[codim].begin();
1762  i!=g->send_interiorborder_overlapfront[codim].end(); ++i)
1763  s << "[" << rank << "]: " << " s_ib_of[" << codim << "] to rank "
1764  << i->rank << " " << i->grid << std::endl;
1765 
1766  for (I i=g->recv_overlapfront_interiorborder[codim].begin();
1767  i!=g->recv_overlapfront_interiorborder[codim].end(); ++i)
1768  s << "[" << rank << "]: " << " r_of_ib[" << codim << "] to rank "
1769  << i->rank << " " << i->grid << std::endl;
1770  }
1771  }
1772 
1773  s << std::endl;
1774 
1775  return s;
1776  }
1777 
1778  namespace Capabilities
1779  {
1780 
1788  template<int dim, class Coordinates>
1789  struct hasBackupRestoreFacilities< YaspGrid<dim, Coordinates> >
1790  {
1791  static const bool v = true;
1792  };
1793 
1797  template<int dim, class Coordinates>
1798  struct hasSingleGeometryType< YaspGrid<dim, Coordinates> >
1799  {
1800  static const bool v = true;
1801  static const unsigned int topologyId = Impl::CubeTopology< dim >::type::id;
1802  };
1803 
1807  template<int dim, class Coordinates>
1808  struct isCartesian< YaspGrid<dim, Coordinates> >
1809  {
1810  static const bool v = true;
1811  };
1812 
1816  template<int dim, class Coordinates, int codim>
1817  struct hasEntity< YaspGrid<dim, Coordinates>, codim>
1818  {
1819  static const bool v = true;
1820  };
1821 
1826  template<int dim, class Coordinates, int codim>
1827  struct hasEntityIterator<YaspGrid<dim, Coordinates>, codim>
1828  {
1829  static const bool v = true;
1830  };
1831 
1835  template<int dim, int codim, class Coordinates>
1836  struct canCommunicate< YaspGrid< dim, Coordinates>, codim >
1837  {
1838  static const bool v = true;
1839  };
1840 
1844  template<int dim, class Coordinates>
1845  struct isLevelwiseConforming< YaspGrid<dim, Coordinates> >
1846  {
1847  static const bool v = true;
1848  };
1849 
1853  template<int dim, class Coordinates>
1854  struct isLeafwiseConforming< YaspGrid<dim, Coordinates> >
1855  {
1856  static const bool v = true;
1857  };
1858 
1859  }
1860 
1861 } // end namespace
1862 
1863 // Include the specialization of the StructuredGridFactory class for YaspGrid
1865 // Include the specialization of the BackupRestoreFacility class for YaspGrid
1866 #include <dune/grid/yaspgrid/backuprestore.hh>
1867 
1868 #endif
A geometry implementation for axis-aligned hypercubes.
Portable very large unsigned integers.
Specialization of CollectiveCommunication for MPI.
Definition: mpicollectivecommunication.hh:143
T max(const T &in) const
Compute the maximum of the argument over all processes and return the result in every process....
Definition: mpicollectivecommunication.hh:225
CommDataHandleIF describes the features of a data handle for communication in parallel runs using the...
Definition: datahandleif.hh:76
Wrapper class for entities.
Definition: entity.hh:64
Container for equidistant coordinates in a YaspGrid.
Definition: coordinates.hh:27
Container for equidistant coordinates in a YaspGrid with non-trivial origin.
Definition: coordinates.hh:125
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:277
constexpr unsigned int dim() const
Return dimension of the type.
Definition: type.hh:572
constexpr bool isCube() const
Return true if entity is a cube of any dimension.
Definition: type.hh:562
Definition: grid.hh:855
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: grid.hh:1026
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:18
GridFamily::Traits::CollectiveCommunication CollectiveCommunication
A type that is a model of Dune::CollectiveCommunication. It provides a portable way for collective co...
Definition: grid.hh:519
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition: intersection.hh:162
A Vector class with statically reserved memory.
Definition: reservedvector.hh:42
size_type size() const
Returns number of elements in the vector.
Definition: reservedvector.hh:183
void resize(size_t s)
Specifies a new size for the vector.
Definition: reservedvector.hh:97
reference back()
Returns reference to last element of vector.
Definition: reservedvector.hh:165
void pop_back()
Erases the last element of the vector, O(1) time.
Definition: reservedvector.hh:111
Coordinate container for a tensor product YaspGrid.
Definition: coordinates.hh:234
Definition: torus.hh:280
Definition: torus.hh:44
Definition: ygrid.hh:73
YGridComponent< Coordinates > move(iTupel v) const
return grid moved by the vector v
Definition: ygrid.hh:261
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:269
implements a collection of multiple std::deque<Intersection> Intersections with neighboring processor...
Definition: ygrid.hh:821
int size() const
return the size of the container, this is the sum of the sizes of all deques
Definition: ygrid.hh:951
Iterator end() const
return iterator pointing to the end of the container
Definition: ygrid.hh:927
Iterator begin() const
return iterator pointing to the begin of the container
Definition: ygrid.hh:921
Iterator over a collection o YGrids A YGrid::Iterator is the heart of an entity in YaspGrid.
Definition: ygrid.hh:592
implements a collection of YGridComponents which form a codimension Entities of given codimension c n...
Definition: ygrid.hh:549
Implement the default load balance strategy of yaspgrid.
Definition: partitioning.hh:35
a base class for the yaspgrid partitioning strategy The name might be irritating. It will probably ch...
Definition: partitioning.hh:24
persistent, globally unique Ids
Definition: yaspgrididset.hh:23
[ provides Dune::Grid ]
Definition: yaspgrid.hh:165
YaspGridFamily< dim, Coordinates > GridFamily
the GridFamily of this grid
Definition: yaspgrid.hh:717
int ghostSize(int level, int codim) const
return size (= distance in graph) of ghost region
Definition: yaspgrid.hh:1246
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:874
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: yaspgrid.hh:1278
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafend() const
return LeafIterator which points behind the last entity in maxLevel
Definition: yaspgrid.hh:1211
void globalRefine(int refCount)
refine the grid refCount times.
Definition: yaspgrid.hh:1049
int getMark(const typename Traits::template Codim< 0 >::Entity &e) const
returns adaptation mark for given entity
Definition: yaspgrid.hh:1133
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lbegin(int level) const
version without second template parameter for convenience
Definition: yaspgrid.hh:1176
int ghostSize(int codim) const
return size (= distance in graph) of ghost region
Definition: yaspgrid.hh:1252
int overlapSize(int codim) const
return size (= distance in graph) of overlap region
Definition: yaspgrid.hh:1239
int globalSize(int i) const
return number of cells on finest level in given direction on all processors
Definition: yaspgrid.hh:255
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity in maxLevel
Definition: yaspgrid.hh:1204
const Dune::FieldVector< ctype, dim > & domainSize() const
returns the size of the physical domain
Definition: yaspgrid.hh:1296
void postAdapt()
clean up some markers
Definition: yaspgrid.hh:1154
YGridLevelIterator end() const
return iterator pointing to one past the finest level
Definition: yaspgrid.hh:310
const CollectiveCommunicationType & comm() const
return a collective communication object
Definition: yaspgrid.hh:1590
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: yaspgrid.hh:1284
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:802
int maxLevel() const
Definition: yaspgrid.hh:1043
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafend() const
return LeafIterator which points behind the last entity in maxLevel
Definition: yaspgrid.hh:1197
void communicate(CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir, int level) const
Definition: yaspgrid.hh:1305
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
bool preAdapt()
returns true, if the grid will be coarsened
Definition: yaspgrid.hh:1146
iTupel levelSize(int l) const
return size vector of the grid (in cells) on level l
Definition: yaspgrid.hh:273
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:1119
int overlapSize(int level, int codim) const
return size (= distance in graph) of overlap region
Definition: yaspgrid.hh:1232
void communicate(CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir) const
Definition: yaspgrid.hh:1315
int size(int codim) const
number of leaf entities per codim in this process
Definition: yaspgrid.hh:1272
bool isPeriodic(int i) const
return whether the grid is periodic in direction i
Definition: yaspgrid.hh:282
void refineOptions(bool keepPhysicalOverlap)
set options for refinement
Definition: yaspgrid.hh:1103
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:734
int size(int level, int codim) const
number of entities per level and codim in this process
Definition: yaspgrid.hh:1258
const Torus< CollectiveCommunicationType, dim > & torus() const
return reference to torus
Definition: yaspgrid.hh:249
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lend(int level) const
version without second template parameter for convenience
Definition: yaspgrid.hh:1183
bool adapt()
map adapt to global refine
Definition: yaspgrid.hh:1139
ReservedVector< YGridLevel, 32 >::const_iterator YGridLevelIterator
Iterator over the grid levels.
Definition: yaspgrid.hh:293
size_t numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition: yaspgrid.hh:1290
void communicateCodim(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int level) const
Definition: yaspgrid.hh:1325
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:1169
void makelevel(const Coordinates &coords, std::bitset< dim > periodic, iTupel o_interior, int overlap)
Make a new YGridLevel structure.
Definition: yaspgrid.hh:330
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity in maxLevel
Definition: yaspgrid.hh:1190
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator lbegin(int level) const
one past the end on this level
Definition: yaspgrid.hh:1162
iTupel globalSize() const
return number of cells on finest level on all processors
Definition: yaspgrid.hh:261
YGridLevelIterator begin(int i) const
return iterator pointing to given level
Definition: yaspgrid.hh:302
Coordinates::ctype ctype
Type used for coordinates.
Definition: yaspgrid.hh:179
int levelSize(int l, int i) const
return size of the grid (in cells) on level l in direction i
Definition: yaspgrid.hh:267
YGridLevelIterator begin() const
return iterator pointing to coarsest level
Definition: yaspgrid.hh:296
YaspHierarchicIterator enables iteration over son entities of codim 0.
Definition: yaspgridhierarchiciterator.hh:18
Implementation of Level- and LeafIndexSets for YaspGrid.
Definition: yaspgridindexsets.hh:23
YaspIntersectionIterator enables iteration over intersections with neighboring codim 0 entities.
Definition: yaspgridintersectioniterator.hh:20
YaspIntersection provides data about intersection with neighboring codim 0 entities.
Definition: yaspgridintersection.hh:20
Iterates over entities of one grid level.
Definition: yaspgridleveliterator.hh:17
implement a consecutive index for all entities of given codim of a YaspGrid
Definition: yaspgridpersistentcontainer.hh:33
Implements an utility class that provides collective communication methods for sequential programs.
A set of traits classes to store static information about grid implementation.
Different resources needed by all grid implementations.
This provides container classes for the coordinates to be used in YaspGrid Upon implementation of the...
Describes the parallel communication interface class for MessageBuffers and DataHandles.
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
Various implementations of the power function for run-time and static arguments.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
CommunicationDirection
Define a type for communication direction parameter.
Definition: gridenums.hh:168
InterfaceType
Parameter to be used for the communication functions.
Definition: gridenums.hh:84
@ All_Partition
all entities
Definition: gridenums.hh:139
@ Interior_Partition
only interior entities
Definition: gridenums.hh:135
@ InteriorBorder_Partition
interior and border entities
Definition: gridenums.hh:136
@ Overlap_Partition
interior, border, and overlap entities
Definition: gridenums.hh:137
@ Ghost_Partition
only ghost entities
Definition: gridenums.hh:140
@ BackwardCommunication
reverse communication direction
Definition: gridenums.hh:170
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition: gridenums.hh:86
@ All_All_Interface
send all and receive all entities
Definition: gridenums.hh:89
@ Overlap_All_Interface
send overlap, receive all entities
Definition: gridenums.hh:88
@ Overlap_OverlapFront_Interface
send overlap, receive overlap and front entities
Definition: gridenums.hh:87
@ InteriorBorder_InteriorBorder_Interface
send/receive interior and border entities
Definition: gridenums.hh:85
Provides base classes for index and id sets.
Implements an utility class that provides MPI's collective communication methods.
Helpers for dealing with MPI.
Dune namespace.
Definition: alignedallocator.hh:10
An stl-compliant random-access container which stores everything on the stack.
specialize with 'true' for all codims that a grid can communicate data on (default=false)
Definition: capabilities.hh:95
specialize with 'true' for all codims that a grid provides an iterator for (default=false)
Definition: capabilities.hh:72
Specialize with 'true' for all codims that a grid implements entities for. (default=false)
Definition: capabilities.hh:56
Specialize with 'true' for if the codimension 0 entity of the grid has only one possible geometry typ...
Definition: capabilities.hh:25
Specialize with 'true' if the grid is a Cartesian grid. Cartesian grids satisfy the following propert...
Definition: capabilities.hh:46
Specialize with 'true' if implementation guarantees a conforming leaf grid. (default=false)
Definition: capabilities.hh:113
Specialize with 'true' if implementation guarantees conforming level grids. (default=false)
Definition: capabilities.hh:104
Static tag representing a codimension.
Definition: dimension.hh:22
A traits struct that collects all associated types of one grid model.
Definition: grid.hh:1064
IdSet< const GridImp, LocalIdSetImp, LIDType > LocalIdSet
The type of the local id set.
Definition: grid.hh:1137
IdSet< const GridImp, GlobalIdSetImp, GIDType > GlobalIdSet
The type of the global id set.
Definition: grid.hh:1135
IndexSet< const GridImp, LeafIndexSetImp > LeafIndexSet
The type of the leaf index set.
Definition: grid.hh:1133
IndexSet< const GridImp, LevelIndexSetImp > LevelIndexSet
The type of the level index set.
Definition: grid.hh:1131
Calculates m^p at compile time.
Definition: power.hh:20
type describing an intersection with a neighboring processor
Definition: ygrid.hh:827
Specialization of the StructuredGridFactory class for YaspGrid.
This file provides the infrastructure for toroidal communication in YaspGrid.
A unique label for each type of element that can occur in a grid.
Traits for type conversions and type information.
the YaspEntity class and its specializations
The YaspEntitySeed class.
The YaspGeometry class and its specializations.
level-wise, non-persistent, consecutive indices for YaspGrid
The YaspIntersection class.
The YaspIntersectionIterator class.
The YaspLevelIterator class.
Specialization of the PersistentContainer for YaspGrid.
This provides a YGrid, the elemental component of the yaspgrid implementation.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 18, 22:30, 2024)