Dune Core Modules (2.9.1)

yaspgrid.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_GRID_YASPGRID_HH
6#define DUNE_GRID_YASPGRID_HH
7
8#include <cstdint>
9#include <iostream>
10#include <iterator>
11#include <vector>
12#include <algorithm>
13#include <stack>
14#include <type_traits>
15
16#include <dune/grid/common/backuprestore.hh>
17#include <dune/grid/common/grid.hh> // the grid base classes
18#include <dune/grid/common/capabilities.hh> // the capabilities
19#include <dune/common/hybridutilities.hh>
21#include <dune/common/math.hh>
27#include <dune/geometry/type.hh>
30
31
32#if HAVE_MPI
34#endif
35
43namespace Dune {
44
45 /* some sizes for building global ids
46 */
47 const int yaspgrid_dim_bits = 24; // bits for encoding each dimension
48 const int yaspgrid_level_bits = 5; // bits for encoding level number
49
50
51 //************************************************************************
52 // forward declaration of templates
53
54 template<int dim, class Coordinates> class YaspGrid;
55 template<int mydim, int cdim, class GridImp> class YaspGeometry;
56 template<int codim, int dim, class GridImp> class YaspEntity;
57 template<int codim, class GridImp> class YaspEntitySeed;
58 template<int codim, PartitionIteratorType pitype, class GridImp> class YaspLevelIterator;
59 template<class GridImp> class YaspIntersectionIterator;
60 template<class GridImp> class YaspIntersection;
61 template<class GridImp> class YaspHierarchicIterator;
62 template<class GridImp, bool isLeafIndexSet> class YaspIndexSet;
63 template<class GridImp> class YaspGlobalIdSet;
64 template<class GridImp> class YaspPersistentContainerIndex;
65
66} // namespace Dune
67
79#include <dune/grid/yaspgrid/yaspgrididset.hh>
81
82namespace Dune {
83
84#if HAVE_MPI
85 using YaspCommunication = Communication<MPI_Comm>;
86#else
87 using YaspCommunication = Communication<No_Comm>;
88#endif
89
90 template<int dim, class Coordinates>
91 struct YaspGridFamily
92 {
93 typedef YaspCommunication CCType;
94
95 typedef GridTraits<dim, // dimension of the grid
96 dim, // dimension of the world space
98 YaspGeometry,YaspEntity,
99 YaspLevelIterator, // type used for the level iterator
100 YaspIntersection, // leaf intersection
101 YaspIntersection, // level intersection
102 YaspIntersectionIterator, // leaf intersection iter
103 YaspIntersectionIterator, // level intersection iter
104 YaspHierarchicIterator,
105 YaspLevelIterator, // type used for the leaf(!) iterator
106 YaspIndexSet< const YaspGrid< dim, Coordinates >, false >, // level index set
107 YaspIndexSet< const YaspGrid< dim, Coordinates >, true >, // leaf index set
108 YaspGlobalIdSet<const YaspGrid<dim, Coordinates> >,
109 bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim>,
110 YaspGlobalIdSet<const YaspGrid<dim, Coordinates> >,
111 bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim>,
112 CCType,
113 DefaultLevelGridViewTraits, DefaultLeafGridViewTraits,
114 YaspEntitySeed>
115 Traits;
116 };
117
118#ifndef DOXYGEN
119 template<int dim, int codim>
120 struct YaspCommunicateMeta {
121 template<class G, class DataHandle>
122 static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
123 {
124 if (data.contains(dim,codim))
125 {
126 g.template communicateCodim<DataHandle,codim>(data,iftype,dir,level);
127 }
128 YaspCommunicateMeta<dim,codim-1>::comm(g,data,iftype,dir,level);
129 }
130 };
131
132 template<int dim>
133 struct YaspCommunicateMeta<dim,0> {
134 template<class G, class DataHandle>
135 static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
136 {
137 if (data.contains(dim,0))
138 g.template communicateCodim<DataHandle,0>(data,iftype,dir,level);
139 }
140 };
141#endif
142
143 //************************************************************************
160 template<int dim, class Coordinates = EquidistantCoordinates<double, dim> >
162 : public GridDefaultImplementation<dim,dim,typename Coordinates::ctype,YaspGridFamily<dim, Coordinates> >
163 {
164
165 template<int, PartitionIteratorType, typename>
166 friend class YaspLevelIterator;
167
168 template<typename>
169 friend class YaspHierarchicIterator;
170
172 protected:
173
174 public:
176 typedef typename Coordinates::ctype ctype;
177
178 using Communication = typename Base::Communication;
179 using CollectiveCommunicationType [[deprecated("Use Communication. Will be removed after release 2.9")]] = typename Base::Communication;
180
181#ifndef DOXYGEN
182 typedef typename Dune::YGrid<Coordinates> YGrid;
184
187 struct YGridLevel {
188
190 int level() const
191 {
192 return level_;
193 }
194
195 Coordinates coords;
196
197 std::array<YGrid, dim+1> overlapfront;
198 std::array<YGridComponent<Coordinates>, Dune::power(2,dim)> overlapfront_data;
199 std::array<YGrid, dim+1> overlap;
200 std::array<YGridComponent<Coordinates>, Dune::power(2,dim)> overlap_data;
201 std::array<YGrid, dim+1> interiorborder;
202 std::array<YGridComponent<Coordinates>, Dune::power(2,dim)> interiorborder_data;
203 std::array<YGrid, dim+1> interior;
204 std::array<YGridComponent<Coordinates>, Dune::power(2,dim)> interior_data;
205
206 std::array<YGridList<Coordinates>,dim+1> send_overlapfront_overlapfront;
207 std::array<std::deque<Intersection>, Dune::power(2,dim)> send_overlapfront_overlapfront_data;
208 std::array<YGridList<Coordinates>,dim+1> recv_overlapfront_overlapfront;
209 std::array<std::deque<Intersection>, Dune::power(2,dim)> recv_overlapfront_overlapfront_data;
210
211 std::array<YGridList<Coordinates>,dim+1> send_overlap_overlapfront;
212 std::array<std::deque<Intersection>, Dune::power(2,dim)> send_overlap_overlapfront_data;
213 std::array<YGridList<Coordinates>,dim+1> recv_overlapfront_overlap;
214 std::array<std::deque<Intersection>, Dune::power(2,dim)> recv_overlapfront_overlap_data;
215
216 std::array<YGridList<Coordinates>,dim+1> send_interiorborder_interiorborder;
217 std::array<std::deque<Intersection>, Dune::power(2,dim)> send_interiorborder_interiorborder_data;
218 std::array<YGridList<Coordinates>,dim+1> recv_interiorborder_interiorborder;
219 std::array<std::deque<Intersection>, Dune::power(2,dim)> recv_interiorborder_interiorborder_data;
220
221 std::array<YGridList<Coordinates>,dim+1> send_interiorborder_overlapfront;
222 std::array<std::deque<Intersection>, Dune::power(2,dim)> send_interiorborder_overlapfront_data;
223 std::array<YGridList<Coordinates>,dim+1> recv_overlapfront_interiorborder;
224 std::array<std::deque<Intersection>, Dune::power(2,dim)> recv_overlapfront_interiorborder_data;
225
226 // general
227 YaspGrid<dim,Coordinates>* mg; // each grid level knows its multigrid
228 int overlapSize; // in mesh cells on this level
229 bool keepOverlap;
230
232 int level_;
233 };
234
236 typedef std::array<int, dim> iTupel;
237 typedef FieldVector<ctype, dim> fTupel;
238
239 // communication tag used by multigrid
240 enum { tag = 17 };
241#endif
242
245 {
246 return _torus;
247 }
248
250 int globalSize(int i) const
251 {
252 return levelSize(maxLevel(),i);
253 }
254
256 iTupel globalSize() const
257 {
258 return levelSize(maxLevel());
259 }
260
262 int levelSize(int l, int i) const
263 {
264 return _coarseSize[i] * (1 << l);
265 }
266
268 iTupel levelSize(int l) const
269 {
270 iTupel s;
271 for (int i=0; i<dim; ++i)
272 s[i] = levelSize(l,i);
273 return s;
274 }
275
277 bool isPeriodic(int i) const
278 {
279 return _periodic[i];
280 }
281
282 bool getRefineOption() const
283 {
284 return keep_ovlp;
285 }
286
289
292 {
293 return _levels.begin();
294 }
295
298 {
299 if (i<0 || i>maxLevel())
300 DUNE_THROW(GridError, "level not existing");
301 return std::next(_levels.begin(), i);
302 }
303
306 {
307 return _levels.end();
308 }
309
310 // static method to create the default load balance strategy
311 [[deprecated("use defaultPartitioner")]]
312 static const YLoadBalanceDefault<dim>* defaultLoadbalancer()
313 {
314 static YLoadBalanceDefault<dim> lb;
315 return & lb;
316 }
317
318 // static method to create the default partitioning strategy
319 static const Yasp::Partitioning<dim>* defaultPartitioner()
320 {
321 static Yasp::DefaultPartitioning<dim> lb;
322 return & lb;
323 }
324
325 protected:
333 void makelevel (const Coordinates& coords, std::bitset<dim> periodic, iTupel o_interior, int overlap)
334 {
335 YGridLevel& g = _levels.back();
336 g.overlapSize = overlap;
337 g.mg = this;
338 g.level_ = maxLevel();
339 g.coords = coords;
340 g.keepOverlap = keep_ovlp;
341
342 using Dune::power;
343
344 // set the inserting positions in the corresponding arrays of YGridLevelStructure
345 typename std::array<YGridComponent<Coordinates>, power(2,dim)>::iterator overlapfront_it = g.overlapfront_data.begin();
346 typename std::array<YGridComponent<Coordinates>, power(2,dim)>::iterator overlap_it = g.overlap_data.begin();
347 typename std::array<YGridComponent<Coordinates>, power(2,dim)>::iterator interiorborder_it = g.interiorborder_data.begin();
348 typename std::array<YGridComponent<Coordinates>, power(2,dim)>::iterator interior_it = g.interior_data.begin();
349
350 typename std::array<std::deque<Intersection>, power(2,dim)>::iterator
351 send_overlapfront_overlapfront_it = g.send_overlapfront_overlapfront_data.begin();
352 typename std::array<std::deque<Intersection>, power(2,dim)>::iterator
353 recv_overlapfront_overlapfront_it = g.recv_overlapfront_overlapfront_data.begin();
354
355 typename std::array<std::deque<Intersection>, power(2,dim)>::iterator
356 send_overlap_overlapfront_it = g.send_overlap_overlapfront_data.begin();
357 typename std::array<std::deque<Intersection>, power(2,dim)>::iterator
358 recv_overlapfront_overlap_it = g.recv_overlapfront_overlap_data.begin();
359
360 typename std::array<std::deque<Intersection>, power(2,dim)>::iterator
361 send_interiorborder_interiorborder_it = g.send_interiorborder_interiorborder_data.begin();
362 typename std::array<std::deque<Intersection>, power(2,dim)>::iterator
363 recv_interiorborder_interiorborder_it = g.recv_interiorborder_interiorborder_data.begin();
364
365 typename std::array<std::deque<Intersection>, power(2,dim)>::iterator
366 send_interiorborder_overlapfront_it = g.send_interiorborder_overlapfront_data.begin();
367 typename std::array<std::deque<Intersection>, power(2,dim)>::iterator
368 recv_overlapfront_interiorborder_it = g.recv_overlapfront_interiorborder_data.begin();
369
370 // have a null array for constructor calls around
371 std::array<int,dim> n;
372 std::fill(n.begin(), n.end(), 0);
373
374 // determine origin of the grid with overlap and store whether an overlap area exists in direction i.
375 std::bitset<dim> ovlp_low(0ULL);
376 std::bitset<dim> ovlp_up(0ULL);
377
378 iTupel o_overlap;
379 iTupel s_overlap;
380
381 // determine at where we have overlap and how big the size of the overlap partition is
382 for (int i=0; i<dim; i++)
383 {
384 // the coordinate container has been contructed to hold the entire grid on
385 // this processor, including overlap. this is the element size.
386 s_overlap[i] = g.coords.size(i);
387
388 //in the periodic case there is always overlap
389 if (periodic[i])
390 {
391 o_overlap[i] = o_interior[i]-overlap;
392 ovlp_low[i] = true;
393 ovlp_up[i] = true;
394 }
395 else
396 {
397 //check lower boundary
398 if (o_interior[i] - overlap < 0)
399 o_overlap[i] = 0;
400 else
401 {
402 o_overlap[i] = o_interior[i] - overlap;
403 ovlp_low[i] = true;
404 }
405
406 //check upper boundary
407 if (o_overlap[i] + g.coords.size(i) < globalSize(i))
408 ovlp_up[i] = true;
409 }
410 }
411
412 for (unsigned int codim = 0; codim < dim + 1; codim++)
413 {
414 // set the begin iterator for the corresponding ygrids
415 g.overlapfront[codim].setBegin(overlapfront_it);
416 g.overlap[codim].setBegin(overlap_it);
417 g.interiorborder[codim].setBegin(interiorborder_it);
418 g.interior[codim].setBegin(interior_it);
419 g.send_overlapfront_overlapfront[codim].setBegin(send_overlapfront_overlapfront_it);
420 g.recv_overlapfront_overlapfront[codim].setBegin(recv_overlapfront_overlapfront_it);
421 g.send_overlap_overlapfront[codim].setBegin(send_overlap_overlapfront_it);
422 g.recv_overlapfront_overlap[codim].setBegin(recv_overlapfront_overlap_it);
423 g.send_interiorborder_interiorborder[codim].setBegin(send_interiorborder_interiorborder_it);
424 g.recv_interiorborder_interiorborder[codim].setBegin(recv_interiorborder_interiorborder_it);
425 g.send_interiorborder_overlapfront[codim].setBegin(send_interiorborder_overlapfront_it);
426 g.recv_overlapfront_interiorborder[codim].setBegin(recv_overlapfront_interiorborder_it);
427
428 // find all combinations of unit vectors that span entities of the given codimension
429 for (unsigned int index = 0; index < (1<<dim); index++)
430 {
431 // check whether the given shift is of our codimension
432 std::bitset<dim> r(index);
433 if (r.count() != dim-codim)
434 continue;
435
436 // get an origin and a size array for subsequent modification
437 std::array<int,dim> origin(o_overlap);
438 std::array<int,dim> size(s_overlap);
439
440 // build overlapfront
441 // we have to extend the element size by one in all directions without shift.
442 for (int i=0; i<dim; i++)
443 if (!r[i])
444 size[i]++;
445 *overlapfront_it = YGridComponent<Coordinates>(origin, r, &g.coords, size, n, size);
446
447 // build overlap
448 for (int i=0; i<dim; i++)
449 {
450 if (!r[i])
451 {
452 if (ovlp_low[i])
453 {
454 origin[i]++;
455 size[i]--;
456 }
457 if (ovlp_up[i])
458 size[i]--;
459 }
460 }
461 *overlap_it = YGridComponent<Coordinates>(origin,size,*overlapfront_it);
462
463 // build interiorborder
464 for (int i=0; i<dim; i++)
465 {
466 if (ovlp_low[i])
467 {
468 origin[i] += overlap;
469 size[i] -= overlap;
470 if (!r[i])
471 {
472 origin[i]--;
473 size[i]++;
474 }
475 }
476 if (ovlp_up[i])
477 {
478 size[i] -= overlap;
479 if (!r[i])
480 size[i]++;
481 }
482 }
483 *interiorborder_it = YGridComponent<Coordinates>(origin,size,*overlapfront_it);
484
485 // build interior
486 for (int i=0; i<dim; i++)
487 {
488 if (!r[i])
489 {
490 if (ovlp_low[i])
491 {
492 origin[i]++;
493 size[i]--;
494 }
495 if (ovlp_up[i])
496 size[i]--;
497 }
498 }
499 *interior_it = YGridComponent<Coordinates>(origin, size, *overlapfront_it);
500
501 intersections(*overlapfront_it,*overlapfront_it,*send_overlapfront_overlapfront_it, *recv_overlapfront_overlapfront_it);
502 intersections(*overlap_it,*overlapfront_it,*send_overlap_overlapfront_it, *recv_overlapfront_overlap_it);
503 intersections(*interiorborder_it,*interiorborder_it,*send_interiorborder_interiorborder_it,*recv_interiorborder_interiorborder_it);
504 intersections(*interiorborder_it,*overlapfront_it,*send_interiorborder_overlapfront_it,*recv_overlapfront_interiorborder_it);
505
506 // advance all iterators pointing to the next insertion point
507 ++overlapfront_it;
508 ++overlap_it;
509 ++interiorborder_it;
510 ++interior_it;
511 ++send_overlapfront_overlapfront_it;
512 ++recv_overlapfront_overlapfront_it;
513 ++send_overlap_overlapfront_it;
514 ++recv_overlapfront_overlap_it;
515 ++send_interiorborder_interiorborder_it;
516 ++recv_interiorborder_interiorborder_it;
517 ++send_interiorborder_overlapfront_it;
518 ++recv_overlapfront_interiorborder_it;
519 }
520
521 // set end iterators in the corresonding ygrids
522 g.overlapfront[codim].finalize(overlapfront_it);
523 g.overlap[codim].finalize(overlap_it);
524 g.interiorborder[codim].finalize(interiorborder_it);
525 g.interior[codim].finalize(interior_it);
526 g.send_overlapfront_overlapfront[codim].finalize(send_overlapfront_overlapfront_it,g.overlapfront[codim]);
527 g.recv_overlapfront_overlapfront[codim].finalize(recv_overlapfront_overlapfront_it,g.overlapfront[codim]);
528 g.send_overlap_overlapfront[codim].finalize(send_overlap_overlapfront_it,g.overlapfront[codim]);
529 g.recv_overlapfront_overlap[codim].finalize(recv_overlapfront_overlap_it,g.overlapfront[codim]);
530 g.send_interiorborder_interiorborder[codim].finalize(send_interiorborder_interiorborder_it,g.overlapfront[codim]);
531 g.recv_interiorborder_interiorborder[codim].finalize(recv_interiorborder_interiorborder_it,g.overlapfront[codim]);
532 g.send_interiorborder_overlapfront[codim].finalize(send_interiorborder_overlapfront_it,g.overlapfront[codim]);
533 g.recv_overlapfront_interiorborder[codim].finalize(recv_overlapfront_interiorborder_it,g.overlapfront[codim]);
534 }
535 }
536
537#ifndef DOXYGEN
546 struct mpifriendly_ygrid {
547 mpifriendly_ygrid ()
548 {
549 std::fill(origin.begin(), origin.end(), 0);
550 std::fill(size.begin(), size.end(), 0);
551 }
552 mpifriendly_ygrid (const YGridComponent<Coordinates>& grid)
553 : origin(grid.origin()), size(grid.size())
554 {}
555 iTupel origin;
556 iTupel size;
557 };
558#endif
559
569 std::deque<Intersection>& sendlist, std::deque<Intersection>& recvlist)
570 {
571 iTupel size = globalSize();
572
573 // the exchange buffers
574 std::vector<YGridComponent<Coordinates> > send_recvgrid(_torus.neighbors());
575 std::vector<YGridComponent<Coordinates> > recv_recvgrid(_torus.neighbors());
576 std::vector<YGridComponent<Coordinates> > send_sendgrid(_torus.neighbors());
577 std::vector<YGridComponent<Coordinates> > recv_sendgrid(_torus.neighbors());
578
579 // new exchange buffers to send simple struct without virtual functions
580 std::vector<mpifriendly_ygrid> mpifriendly_send_recvgrid(_torus.neighbors());
581 std::vector<mpifriendly_ygrid> mpifriendly_recv_recvgrid(_torus.neighbors());
582 std::vector<mpifriendly_ygrid> mpifriendly_send_sendgrid(_torus.neighbors());
583 std::vector<mpifriendly_ygrid> mpifriendly_recv_sendgrid(_torus.neighbors());
584
585 // fill send buffers; iterate over all neighboring processes
586 // non-periodic case is handled automatically because intersection will be zero
587 for (typename Torus<Communication,dim>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
588 {
589 // determine if we communicate with this neighbor (and what)
590 bool skip = false;
591 iTupel coord = _torus.coord(); // my coordinates
592 iTupel delta = i.delta(); // delta to neighbor
593 iTupel nb = coord; // the neighbor
594 for (int k=0; k<dim; k++) nb[k] += delta[k];
595 iTupel v; // grid movement
596 std::fill(v.begin(), v.end(), 0);
597
598 for (int k=0; k<dim; k++)
599 {
600 if (nb[k]<0)
601 {
602 if (_periodic[k])
603 v[k] += size[k];
604 else
605 skip = true;
606 }
607 if (nb[k]>=_torus.dims(k))
608 {
609 if (_periodic[k])
610 v[k] -= size[k];
611 else
612 skip = true;
613 }
614 // neither might be true, then v=0
615 }
616
617 // store moved grids in send buffers
618 if (!skip)
619 {
620 send_sendgrid[i.index()] = sendgrid.move(v);
621 send_recvgrid[i.index()] = recvgrid.move(v);
622 }
623 else
624 {
625 send_sendgrid[i.index()] = YGridComponent<Coordinates>();
626 send_recvgrid[i.index()] = YGridComponent<Coordinates>();
627 }
628 }
629
630 // issue send requests for sendgrid being sent to all neighbors
631 for (typename Torus<Communication,dim>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
632 {
633 mpifriendly_send_sendgrid[i.index()] = mpifriendly_ygrid(send_sendgrid[i.index()]);
634 _torus.send(i.rank(), &mpifriendly_send_sendgrid[i.index()], sizeof(mpifriendly_ygrid));
635 }
636
637 // issue recv requests for sendgrids of neighbors
638 for (typename Torus<Communication,dim>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
639 _torus.recv(i.rank(), &mpifriendly_recv_sendgrid[i.index()], sizeof(mpifriendly_ygrid));
640
641 // exchange the sendgrids
642 _torus.exchange();
643
644 // issue send requests for recvgrid being sent to all neighbors
645 for (typename Torus<Communication,dim>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
646 {
647 mpifriendly_send_recvgrid[i.index()] = mpifriendly_ygrid(send_recvgrid[i.index()]);
648 _torus.send(i.rank(), &mpifriendly_send_recvgrid[i.index()], sizeof(mpifriendly_ygrid));
649 }
650
651 // issue recv requests for recvgrid of neighbors
652 for (typename Torus<Communication,dim>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
653 _torus.recv(i.rank(), &mpifriendly_recv_recvgrid[i.index()], sizeof(mpifriendly_ygrid));
654
655 // exchange the recvgrid
656 _torus.exchange();
657
658 // process receive buffers and compute intersections
659 for (typename Torus<Communication,dim>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
660 {
661 // what must be sent to this neighbor
662 Intersection send_intersection;
663 mpifriendly_ygrid yg = mpifriendly_recv_recvgrid[i.index()];
664 recv_recvgrid[i.index()] = YGridComponent<Coordinates>(yg.origin,yg.size);
665 send_intersection.grid = sendgrid.intersection(recv_recvgrid[i.index()]);
666 send_intersection.rank = i.rank();
667 send_intersection.distance = i.distance();
668 if (!send_intersection.grid.empty()) sendlist.push_front(send_intersection);
669
670 Intersection recv_intersection;
671 yg = mpifriendly_recv_sendgrid[i.index()];
672 recv_sendgrid[i.index()] = YGridComponent<Coordinates>(yg.origin,yg.size);
673 recv_intersection.grid = recvgrid.intersection(recv_sendgrid[i.index()]);
674 recv_intersection.rank = i.rank();
675 recv_intersection.distance = i.distance();
676 if(!recv_intersection.grid.empty()) recvlist.push_back(recv_intersection);
677 }
678 }
679
680 protected:
681
682 typedef const YaspGrid<dim,Coordinates> GridImp;
683
684 void init()
685 {
686 indexsets.push_back( std::make_shared< YaspIndexSet<const YaspGrid<dim, Coordinates>, false > >(*this,0) );
687 boundarysegmentssize();
688 }
689
690 void boundarysegmentssize()
691 {
692 // sizes of local macro grid
693 std::array<int, dim> sides;
694 {
695 for (int i=0; i<dim; i++)
696 {
697 sides[i] =
698 ((begin()->overlap[0].dataBegin()->origin(i) == 0)+
699 (begin()->overlap[0].dataBegin()->origin(i) + begin()->overlap[0].dataBegin()->size(i)
700 == levelSize(0,i)));
701 }
702 }
703 nBSegments = 0;
704 for (int k=0; k<dim; k++)
705 {
706 int offset = 1;
707 for (int l=0; l<dim; l++)
708 {
709 if (l==k) continue;
710 offset *= begin()->overlap[0].dataBegin()->size(l);
711 }
712 nBSegments += sides[k]*offset;
713 }
714 }
715
716 public:
717
718 // define the persistent index type
719 typedef bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+dim> PersistentIndexType;
720
722 typedef YaspGridFamily<dim, Coordinates> GridFamily;
723 // the Traits
725
726 // need for friend declarations in entity
730
738 YaspGrid (const Coordinates& coordinates,
739 std::bitset<dim> periodic = std::bitset<dim>(0ULL),
740 int overlap = 1,
741 Communication comm = Communication(),
742 const Yasp::Partitioning<dim>* partitioner = defaultPartitioner())
743 : ccobj(comm)
744 , leafIndexSet_(*this)
745 , _periodic(periodic)
746 , _overlap(overlap)
747 , keep_ovlp(true)
748 , adaptRefCount(0)
749 , adaptActive(false)
750 {
751 _levels.resize(1);
752
753 // Number of elements per coordinate direction on the coarsest level
754 for (std::size_t i=0; i<dim; i++)
755 _coarseSize[i] = coordinates.size(i);
756
757 // Construct the communication torus
758 _torus = decltype(_torus)(comm,tag,_coarseSize,overlap,partitioner);
759
760 iTupel o;
761 std::fill(o.begin(), o.end(), 0);
762 iTupel o_interior(o);
763 iTupel s_interior(_coarseSize);
764
765 _torus.partition(_torus.rank(),o,_coarseSize,o_interior,s_interior);
766
767 // Set domain size
768 if (std::is_same<Coordinates,EquidistantCoordinates<ctype,dim> >::value
769 || std::is_same<Coordinates,EquidistantOffsetCoordinates<ctype,dim> >::value)
770 {
771 for (std::size_t i=0; i<dim; i++)
772 _L[i] = coordinates.size(i) * coordinates.meshsize(i,0);
773 }
774 if (std::is_same<Coordinates,TensorProductCoordinates<ctype,dim> >::value)
775 {
776 //determine sizes of vector to correctly construct torus structure and store for later size requests
777 for (int i=0; i<dim; i++)
778 _L[i] = coordinates.coordinate(i,_coarseSize[i]) - coordinates.coordinate(i,0);
779 }
780
781#if HAVE_MPI
782 // TODO: Settle on a single value for all coordinate types
783 int mysteryFactor = (std::is_same<Coordinates,EquidistantCoordinates<ctype,dim> >::value) ? 1 : 2;
784
785 // check whether the grid is large enough to be overlapping
786 for (int i=0; i<dim; i++)
787 {
788 // find out whether the grid is too small to
789 int toosmall = (s_interior[i] <= mysteryFactor * overlap) && // interior is very small
790 (periodic[i] || (s_interior[i] != _coarseSize[i])); // there is an overlap in that direction
791 // communicate the result to all those processes to have all processors error out if one process failed.
792 int global = 0;
793 MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR, comm);
794 if (global)
795 DUNE_THROW(Dune::GridError,"YaspGrid does not support degrees of freedom shared by more than immediately neighboring subdomains."
796 " Note that this also holds for DOFs on subdomain boundaries."
797 " Increase grid elements or decrease overlap accordingly.");
798 }
799#endif // #if HAVE_MPI
800
801 if (std::is_same<Coordinates,EquidistantCoordinates<ctype,dim> >::value
802 || std::is_same<Coordinates,EquidistantOffsetCoordinates<ctype,dim> >::value)
803 {
804 iTupel s_overlap(s_interior);
805 for (int i=0; i<dim; i++)
806 {
807 if ((o_interior[i] - overlap > 0) || (periodic[i]))
808 s_overlap[i] += overlap;
809 if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
810 s_overlap[i] += overlap;
811 }
812
813 FieldVector<ctype,dim> upperRightWithOverlap;
814 for (int i=0; i<dim; i++)
815 upperRightWithOverlap[i] = coordinates.coordinate(i,0) + coordinates.meshsize(i,0) * s_overlap[i];
816
817 if constexpr (std::is_same_v<Coordinates,EquidistantCoordinates<ctype,dim>>)
818 {
819 // New coordinate object that additionally contains the overlap elements
820 EquidistantCoordinates<ctype,dim> coordinatesWithOverlap(upperRightWithOverlap,s_overlap);
821
822 // add level (the this-> is needed to make g++-6 happy)
823 this->makelevel(coordinatesWithOverlap,periodic,o_interior,overlap);
824 }
825
826 if constexpr (std::is_same_v<Coordinates,EquidistantOffsetCoordinates<ctype,dim>>)
827 {
829 for (int i=0; i<dim; i++)
830 lowerleft[i] = coordinates.origin(i);
831
832 // New coordinate object that additionally contains the overlap elements
833 EquidistantOffsetCoordinates<ctype,dim> coordinatesWithOverlap(lowerleft,upperRightWithOverlap,s_overlap);
834
835 // add level (the this-> is needed to make g++-6 happy)
836 this->makelevel(coordinatesWithOverlap,periodic,o_interior,overlap);
837 }
838 }
839
840 if constexpr (std::is_same_v<Coordinates,TensorProductCoordinates<ctype,dim>>)
841 {
842 std::array<std::vector<ctype>,dim> newCoords;
843 std::array<int, dim> offset(o_interior);
844
845 // find the relevant part of the coords vector for this processor and copy it to newCoords
846 for (int i=0; i<dim; ++i)
847 {
848 //define the coordinate range to be used
849 std::size_t begin = o_interior[i];
850 std::size_t end = begin + s_interior[i] + 1;
851
852 // check whether we are not at the physical boundary. In that case overlap is a simple
853 // extension of the coordinate range to be used
854 if (o_interior[i] - overlap > 0)
855 {
856 begin = begin - overlap;
857 offset[i] -= overlap;
858 }
859 if (o_interior[i] + s_interior[i] + overlap < _coarseSize[i])
860 end = end + overlap;
861
862 //copy the selected part in the new coord vector
863 newCoords[i].resize(end-begin);
864 auto newCoordsIt = newCoords[i].begin();
865 for (std::size_t j=begin; j<end; j++)
866 {
867 *newCoordsIt = coordinates.coordinate(i, j);
868 newCoordsIt++;
869 }
870
871 // Check whether we are at the physical boundary and have a periodic grid.
872 // In this case the coordinate vector has to be tweaked manually.
873 if ((periodic[i]) && (o_interior[i] + s_interior[i] + overlap >= _coarseSize[i]))
874 {
875 // we need to add the first <overlap> cells to the end of newcoords
876 for (int j=0; j<overlap; ++j)
877 newCoords[i].push_back(newCoords[i].back() - coordinates.coordinate(i,j) + coordinates.coordinate(i,j+1));
878 }
879
880 if ((periodic[i]) && (o_interior[i] - overlap <= 0))
881 {
882 offset[i] -= overlap;
883
884 // we need to add the last <overlap> cells to the begin of newcoords
885 std::size_t reverseCounter = coordinates.size(i);
886 for (int j=0; j<overlap; ++j, --reverseCounter)
887 newCoords[i].insert(newCoords[i].begin(), newCoords[i].front()
888 - coordinates.coordinate(i,reverseCounter) + coordinates.coordinate(i,reverseCounter-1));
889 }
890 }
891
892 TensorProductCoordinates<ctype,dim> coordinatesWithOverlap(newCoords, offset);
893
894 // add level (the this-> is needed to make g++-6 happy)
895 this->makelevel(coordinatesWithOverlap,periodic,o_interior,overlap);
896 }
897
898 init();
899 }
900
909 template<class C = Coordinates,
910 typename std::enable_if_t< std::is_same_v<C, EquidistantCoordinates<ctype,dim> >, int> = 0>
912 std::array<int, std::size_t{dim}> s,
913 std::bitset<std::size_t{dim}> periodic = std::bitset<std::size_t{dim}>{0ULL},
914 int overlap = 1,
915 Communication comm = Communication(),
916 const Yasp::Partitioning<dim>* partitioner = defaultPartitioner())
917 : ccobj(comm), _torus(comm,tag,s,overlap,partitioner), leafIndexSet_(*this),
918 _L(L), _periodic(periodic), _coarseSize(s), _overlap(overlap),
919 keep_ovlp(true), adaptRefCount(0), adaptActive(false)
920 {
921 _levels.resize(1);
922
923 iTupel o;
924 std::fill(o.begin(), o.end(), 0);
925 iTupel o_interior(o);
926 iTupel s_interior(s);
927
928 _torus.partition(_torus.rank(),o,s,o_interior,s_interior);
929
930#if HAVE_MPI
931 // check whether the grid is large enough to be overlapping
932 for (int i=0; i<dim; i++)
933 {
934 // find out whether the grid is too small to
935 int toosmall = (s_interior[i] < 2*overlap) && // interior is very small
936 (periodic[i] || (s_interior[i] != s[i])); // there is an overlap in that direction
937 // communicate the result to all those processes to have all processors error out if one process failed.
938 int global = 0;
939 MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR, comm);
940 if (global)
941 DUNE_THROW(Dune::GridError,"YaspGrid does not support degrees of freedom shared by more than immediately neighboring subdomains."
942 " Note that this also holds for DOFs on subdomain boundaries."
943 " Increase grid elements or decrease overlap accordingly.");
944 }
945#endif // #if HAVE_MPI
946
947 iTupel s_overlap(s_interior);
948 for (int i=0; i<dim; i++)
949 {
950 if ((o_interior[i] - overlap > 0) || (periodic[i]))
951 s_overlap[i] += overlap;
952 if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
953 s_overlap[i] += overlap;
954 }
955
956 FieldVector<ctype,dim> upperRightWithOverlap;
957
958 for (int i=0; i<dim; i++)
959 upperRightWithOverlap[i] = (L[i] / s[i]) * s_overlap[i];
960
961 // New coordinate object that additionally contains the overlap elements
962 EquidistantCoordinates<ctype,dim> cc(upperRightWithOverlap,s_overlap);
963
964 // add level
965 makelevel(cc,periodic,o_interior,overlap);
966
967 init();
968 }
969
979 template<class C = Coordinates,
980 typename std::enable_if_t< std::is_same_v<C, EquidistantOffsetCoordinates<ctype,dim> >, int> = 0>
983 std::array<int, std::size_t{dim}> s,
984 std::bitset<std::size_t{dim}> periodic = std::bitset<std::size_t{dim}>(0ULL),
985 int overlap = 1,
986 Communication comm = Communication(),
987 const Yasp::Partitioning<dim>* partitioner = defaultPartitioner())
988 : ccobj(comm), _torus(comm,tag,s,overlap,partitioner), leafIndexSet_(*this),
989 _L(upperright - lowerleft),
990 _periodic(periodic), _coarseSize(s), _overlap(overlap),
991 keep_ovlp(true), adaptRefCount(0), adaptActive(false)
992 {
993 _levels.resize(1);
994
995 iTupel o;
996 std::fill(o.begin(), o.end(), 0);
997 iTupel o_interior(o);
998 iTupel s_interior(s);
999
1000 _torus.partition(_torus.rank(),o,s,o_interior,s_interior);
1001
1002#if HAVE_MPI
1003 // check whether the grid is large enough to be overlapping
1004 for (int i=0; i<dim; i++)
1005 {
1006 // find out whether the grid is too small to
1007 int toosmall = (s_interior[i] < 2*overlap) && // interior is very small
1008 (periodic[i] || (s_interior[i] != s[i])); // there is an overlap in that direction
1009 // communicate the result to all those processes to have all processors error out if one process failed.
1010 int global = 0;
1011 MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR, comm);
1012 if (global)
1013 DUNE_THROW(Dune::GridError,"YaspGrid does not support degrees of freedom shared by more than immediately neighboring subdomains."
1014 " Note that this also holds for DOFs on subdomain boundaries."
1015 " Increase grid elements or decrease overlap accordingly.");
1016 }
1017#endif // #if HAVE_MPI
1018
1019 iTupel s_overlap(s_interior);
1020 for (int i=0; i<dim; i++)
1021 {
1022 if ((o_interior[i] - overlap > 0) || (periodic[i]))
1023 s_overlap[i] += overlap;
1024 if ((o_interior[i] + s_interior[i] + overlap <= _coarseSize[i]) || (periodic[i]))
1025 s_overlap[i] += overlap;
1026 }
1027
1028 FieldVector<ctype,dim> upperRightWithOverlap;
1029 for (int i=0; i<dim; i++)
1030 upperRightWithOverlap[i] = lowerleft[i]
1031 + s_overlap[i] * (upperright[i]-lowerleft[i]) / s[i];
1032
1033 EquidistantOffsetCoordinates<ctype,dim> cc(lowerleft,upperRightWithOverlap,s_overlap);
1034
1035 // add level
1036 makelevel(cc,periodic,o_interior,overlap);
1037
1038 init();
1039 }
1040
1048 template<class C = Coordinates,
1049 typename std::enable_if_t< std::is_same_v<C, TensorProductCoordinates<ctype,dim> >, int> = 0>
1050 YaspGrid (std::array<std::vector<ctype>, std::size_t{dim}> coords,
1051 std::bitset<std::size_t{dim}> periodic = std::bitset<std::size_t{dim}>(0ULL),
1052 int overlap = 1,
1053 Communication comm = Communication(),
1054 const Yasp::Partitioning<dim>* partitioner = defaultPartitioner())
1055 : ccobj(comm), _torus(comm,tag,Dune::Yasp::sizeArray<dim>(coords),overlap,partitioner),
1056 leafIndexSet_(*this), _periodic(periodic), _overlap(overlap),
1057 keep_ovlp(true), adaptRefCount(0), adaptActive(false)
1058 {
1059 if (!Dune::Yasp::checkIfMonotonous(coords))
1060 DUNE_THROW(Dune::GridError,"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
1061
1062 _levels.resize(1);
1063
1064 //determine sizes of vector to correctly construct torus structure and store for later size requests
1065 for (int i=0; i<dim; i++) {
1066 _coarseSize[i] = coords[i].size() - 1;
1067 _L[i] = coords[i][_coarseSize[i]] - coords[i][0];
1068 }
1069
1070 iTupel o;
1071 std::fill(o.begin(), o.end(), 0);
1072 iTupel o_interior(o);
1073 iTupel s_interior(_coarseSize);
1074
1075 _torus.partition(_torus.rank(),o,_coarseSize,o_interior,s_interior);
1076
1077#if HAVE_MPI
1078 // check whether the grid is large enough to be overlapping
1079 for (int i=0; i<dim; i++)
1080 {
1081 // find out whether the grid is too small to
1082 int toosmall = (s_interior[i] < 2*overlap) && // interior is very small
1083 (periodic[i] || (s_interior[i] != _coarseSize[i])); // there is an overlap in that direction
1084 // communicate the result to all those processes to have all processors error out if one process failed.
1085 int global = 0;
1086 MPI_Allreduce(&toosmall, &global, 1, MPI_INT, MPI_LOR, comm);
1087 if (global)
1088 DUNE_THROW(Dune::GridError,"YaspGrid does not support degrees of freedom shared by more than immediately neighboring subdomains."
1089 " Note that this also holds for DOFs on subdomain boundaries."
1090 " Increase grid elements or decrease overlap accordingly.");
1091 }
1092#endif // #if HAVE_MPI
1093
1094
1095 std::array<std::vector<ctype>,dim> newcoords;
1096 std::array<int, dim> offset(o_interior);
1097
1098 // find the relevant part of the coords vector for this processor and copy it to newcoords
1099 for (int i=0; i<dim; ++i)
1100 {
1101 //define iterators on coords that specify the coordinate range to be used
1102 typename std::vector<ctype>::iterator begin = coords[i].begin() + o_interior[i];
1103 typename std::vector<ctype>::iterator end = begin + s_interior[i] + 1;
1104
1105 // check whether we are not at the physical boundary. In that case overlap is a simple
1106 // extension of the coordinate range to be used
1107 if (o_interior[i] - overlap > 0)
1108 {
1109 begin = begin - overlap;
1110 offset[i] -= overlap;
1111 }
1112 if (o_interior[i] + s_interior[i] + overlap < _coarseSize[i])
1113 end = end + overlap;
1114
1115 //copy the selected part in the new coord vector
1116 newcoords[i].resize(end-begin);
1117 std::copy(begin, end, newcoords[i].begin());
1118
1119 // check whether we are at the physical boundary and a have a periodic grid.
1120 // In this case the coordinate vector has to be tweaked manually.
1121 if ((periodic[i]) && (o_interior[i] + s_interior[i] + overlap >= _coarseSize[i]))
1122 {
1123 // we need to add the first <overlap> cells to the end of newcoords
1124 typename std::vector<ctype>::iterator it = coords[i].begin();
1125 for (int j=0; j<overlap; ++j)
1126 newcoords[i].push_back(newcoords[i].back() - *it + *(++it));
1127 }
1128
1129 if ((periodic[i]) && (o_interior[i] - overlap <= 0))
1130 {
1131 offset[i] -= overlap;
1132
1133 // we need to add the last <overlap> cells to the begin of newcoords
1134 typename std::vector<ctype>::iterator it = coords[i].end() - 1;
1135 for (int j=0; j<overlap; ++j)
1136 newcoords[i].insert(newcoords[i].begin(), newcoords[i].front() - *it + *(--it));
1137 }
1138 }
1139
1140 TensorProductCoordinates<ctype,dim> cc(newcoords, offset);
1141
1142 // add level
1143 makelevel(cc,periodic,o_interior,overlap);
1144 init();
1145 }
1146
1147 private:
1148
1163 YaspGrid (std::array<std::vector<ctype>, std::size_t{dim}> coords,
1164 std::bitset<std::size_t{dim}> periodic,
1165 int overlap,
1166 Communication comm,
1167 std::array<int,dim> coarseSize,
1168 const Yasp::Partitioning<dim>* partitioner = defaultPartitioner())
1169 : ccobj(comm), _torus(comm,tag,coarseSize,overlap,partitioner), leafIndexSet_(*this),
1170 _periodic(periodic), _coarseSize(coarseSize), _overlap(overlap),
1171 keep_ovlp(true), adaptRefCount(0), adaptActive(false)
1172 {
1173 // check whether YaspGrid has been given the correct template parameter
1174 static_assert(std::is_same<Coordinates,TensorProductCoordinates<ctype,dim> >::value,
1175 "YaspGrid coordinate container template parameter and given constructor values do not match!");
1176
1177 if (!Dune::Yasp::checkIfMonotonous(coords))
1178 DUNE_THROW(Dune::GridError,"Setup of a tensorproduct grid requires monotonous sequences of coordinates.");
1179
1180 for (int i=0; i<dim; i++)
1181 _L[i] = coords[i][coords[i].size() - 1] - coords[i][0];
1182
1183 _levels.resize(1);
1184
1185 std::array<int,dim> o;
1186 std::fill(o.begin(), o.end(), 0);
1187 std::array<int,dim> o_interior(o);
1188 std::array<int,dim> s_interior(coarseSize);
1189
1190 _torus.partition(_torus.rank(),o,coarseSize,o_interior,s_interior);
1191
1192 // get offset by modifying o_interior according to overlap
1193 std::array<int,dim> offset(o_interior);
1194 for (int i=0; i<dim; i++)
1195 if ((periodic[i]) || (o_interior[i] > 0))
1196 offset[i] -= overlap;
1197
1198 TensorProductCoordinates<ctype,dim> cc(coords, offset);
1199
1200 // add level
1201 makelevel(cc,periodic,o_interior,overlap);
1202
1203 init();
1204 }
1205
1206 // the backup restore facility needs to be able to use above constructor
1207 friend struct BackupRestoreFacility<YaspGrid<dim,Coordinates> >;
1208
1209 // do not copy this class
1210 YaspGrid(const YaspGrid&);
1211
1212 public:
1213
1217 int maxLevel() const
1218 {
1219 return _levels.size()-1;
1220 }
1221
1223 void globalRefine (int refCount)
1224 {
1225 if (refCount < -maxLevel())
1226 DUNE_THROW(GridError, "Only " << maxLevel() << " levels left. " <<
1227 "Coarsening " << -refCount << " levels requested!");
1228
1229 // If refCount is negative then coarsen the grid
1230 for (int k=refCount; k<0; k++)
1231 {
1232 // create an empty grid level
1233 YGridLevel empty;
1234 _levels.back() = empty;
1235 // reduce maxlevel
1236 _levels.pop_back();
1237
1238 indexsets.pop_back();
1239 }
1240
1241 // If refCount is positive refine the grid
1242 for (int k=0; k<refCount; k++)
1243 {
1244 // access to coarser grid level
1245 YGridLevel& cg = _levels[maxLevel()];
1246
1247 std::bitset<dim> ovlp_low(0ULL), ovlp_up(0ULL);
1248 for (int i=0; i<dim; i++)
1249 {
1250 if (cg.overlap[0].dataBegin()->origin(i) > 0 || _periodic[i])
1251 ovlp_low[i] = true;
1252 if (cg.overlap[0].dataBegin()->max(i) + 1 < globalSize(i) || _periodic[i])
1253 ovlp_up[i] = true;
1254 }
1255
1256 Coordinates newcont(cg.coords.refine(ovlp_low, ovlp_up, cg.overlapSize, keep_ovlp));
1257
1258 int overlap = (keep_ovlp) ? 2*cg.overlapSize : cg.overlapSize;
1259
1260 //determine new origin
1261 iTupel o_interior;
1262 for (int i=0; i<dim; i++)
1263 o_interior[i] = 2*cg.interior[0].dataBegin()->origin(i);
1264
1265 // add level
1266 _levels.resize(_levels.size() + 1);
1267 makelevel(newcont,_periodic,o_interior,overlap);
1268
1269 indexsets.push_back( std::make_shared<YaspIndexSet<const YaspGrid<dim,Coordinates>, false > >(*this,maxLevel()) );
1270 }
1271 }
1272
1277 void refineOptions (bool keepPhysicalOverlap)
1278 {
1279 keep_ovlp = keepPhysicalOverlap;
1280 }
1281
1293 bool mark( int refCount, const typename Traits::template Codim<0>::Entity & e )
1294 {
1295 assert(adaptActive == false);
1296 if (e.level() != maxLevel()) return false;
1297 adaptRefCount = std::max(adaptRefCount, refCount);
1298 return true;
1299 }
1300
1307 int getMark ( const typename Traits::template Codim<0>::Entity &e ) const
1308 {
1309 return ( e.level() == maxLevel() ) ? adaptRefCount : 0;
1310 }
1311
1313 bool adapt ()
1314 {
1315 globalRefine(adaptRefCount);
1316 return (adaptRefCount > 0);
1317 }
1318
1320 bool preAdapt ()
1321 {
1322 adaptActive = true;
1323 adaptRefCount = comm().max(adaptRefCount);
1324 return (adaptRefCount < 0);
1325 }
1326
1329 {
1330 adaptActive = false;
1331 adaptRefCount = 0;
1332 }
1333
1335 template<int cd, PartitionIteratorType pitype>
1336 typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lbegin (int level) const
1337 {
1338 return levelbegin<cd,pitype>(level);
1339 }
1340
1342 template<int cd, PartitionIteratorType pitype>
1343 typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lend (int level) const
1344 {
1345 return levelend<cd,pitype>(level);
1346 }
1347
1349 template<int cd>
1350 typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lbegin (int level) const
1351 {
1352 return levelbegin<cd,All_Partition>(level);
1353 }
1354
1356 template<int cd>
1357 typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lend (int level) const
1358 {
1359 return levelend<cd,All_Partition>(level);
1360 }
1361
1363 template<int cd, PartitionIteratorType pitype>
1364 typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafbegin () const
1365 {
1366 return levelbegin<cd,pitype>(maxLevel());
1367 }
1368
1370 template<int cd, PartitionIteratorType pitype>
1371 typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafend () const
1372 {
1373 return levelend<cd,pitype>(maxLevel());
1374 }
1375
1377 template<int cd>
1378 typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafbegin () const
1379 {
1380 return levelbegin<cd,All_Partition>(maxLevel());
1381 }
1382
1384 template<int cd>
1385 typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafend () const
1386 {
1387 return levelend<cd,All_Partition>(maxLevel());
1388 }
1389
1390 // \brief obtain Entity from EntitySeed. */
1391 template <typename Seed>
1392 typename Traits::template Codim<Seed::codimension>::Entity
1393 entity(const Seed& seed) const
1394 {
1395 const int codim = Seed::codimension;
1396 YGridLevelIterator g = begin(seed.impl().level());
1397
1398 typedef typename Traits::template Codim<Seed::codimension>::Entity Entity;
1399 typedef YaspEntity<codim,dim,const YaspGrid> EntityImp;
1400 typedef typename YGrid::Iterator YIterator;
1401
1402 return Entity(EntityImp(g,YIterator(g->overlapfront[codim],seed.impl().coord(),seed.impl().offset())));
1403 }
1404
1406 int overlapSize (int level, [[maybe_unused]] int codim) const
1407 {
1408 YGridLevelIterator g = begin(level);
1409 return g->overlapSize;
1410 }
1411
1413 int overlapSize ([[maybe_unused]] int odim) const
1414 {
1416 return g->overlapSize;
1417 }
1418
1420 int ghostSize ([[maybe_unused]] int level, [[maybe_unused]] int codim) const
1421 {
1422 return 0;
1423 }
1424
1426 int ghostSize ([[maybe_unused]] int codim) const
1427 {
1428 return 0;
1429 }
1430
1432 int size (int level, int codim) const
1433 {
1434 YGridLevelIterator g = begin(level);
1435
1436 // sum over all components of the codimension
1437 int count = 0;
1438 typedef typename std::array<YGridComponent<Coordinates>, Dune::power(2,dim)>::iterator DAI;
1439 for (DAI it = g->overlapfront[codim].dataBegin(); it != g->overlapfront[codim].dataEnd(); ++it)
1440 count += it->totalsize();
1441
1442 return count;
1443 }
1444
1446 int size (int codim) const
1447 {
1448 return size(maxLevel(),codim);
1449 }
1450
1452 int size (int level, GeometryType type) const
1453 {
1454 return (type.isCube()) ? size(level,dim-type.dim()) : 0;
1455 }
1456
1458 int size (GeometryType type) const
1459 {
1460 return size(maxLevel(),type);
1461 }
1462
1464 size_t numBoundarySegments () const
1465 {
1466 return nBSegments;
1467 }
1468
1471 return _L;
1472 }
1473
1478 template<class DataHandleImp, class DataType>
1480 {
1481 YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,level);
1482 }
1483
1488 template<class DataHandleImp, class DataType>
1490 {
1491 YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,this->maxLevel());
1492 }
1493
1498 template<class DataHandle, int codim>
1499 void communicateCodim (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level) const
1500 {
1501 // check input
1502 if (!data.contains(dim,codim)) return; // should have been checked outside
1503
1504 // data types
1505 typedef typename DataHandle::DataType DataType;
1506
1507 // access to grid level
1508 YGridLevelIterator g = begin(level);
1509
1510 // find send/recv lists or throw error
1511 const YGridList<Coordinates>* sendlist = 0;
1512 const YGridList<Coordinates>* recvlist = 0;
1513
1515 {
1516 sendlist = &g->send_interiorborder_interiorborder[codim];
1517 recvlist = &g->recv_interiorborder_interiorborder[codim];
1518 }
1519 if (iftype==InteriorBorder_All_Interface)
1520 {
1521 sendlist = &g->send_interiorborder_overlapfront[codim];
1522 recvlist = &g->recv_overlapfront_interiorborder[codim];
1523 }
1525 {
1526 sendlist = &g->send_overlap_overlapfront[codim];
1527 recvlist = &g->recv_overlapfront_overlap[codim];
1528 }
1529 if (iftype==All_All_Interface)
1530 {
1531 sendlist = &g->send_overlapfront_overlapfront[codim];
1532 recvlist = &g->recv_overlapfront_overlapfront[codim];
1533 }
1534
1535 // change communication direction?
1536 if (dir==BackwardCommunication)
1537 std::swap(sendlist,recvlist);
1538
1539 int cnt;
1540
1541 // Size computation (requires communication if variable size)
1542 std::vector<int> send_size(sendlist->size(),-1); // map rank to total number of objects (of type DataType) to be sent
1543 std::vector<int> recv_size(recvlist->size(),-1); // map rank to total number of objects (of type DataType) to be recvd
1544 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
1545 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
1546
1547 // define type to iterate over send and recv lists
1548 typedef typename YGridList<Coordinates>::Iterator ListIt;
1549
1550 if (data.fixedSize(dim,codim))
1551 {
1552 // fixed size: just take a dummy entity, size can be computed without communication
1553 cnt=0;
1554 for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1555 {
1556 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1558 send_size[cnt] = is->grid.totalsize() * data.size(*it);
1559 cnt++;
1560 }
1561 cnt=0;
1562 for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1563 {
1564 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1566 recv_size[cnt] = is->grid.totalsize() * data.size(*it);
1567 cnt++;
1568 }
1569 }
1570 else
1571 {
1572 // variable size case: sender side determines the size
1573 cnt=0;
1574 for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1575 {
1576 // allocate send buffer for sizes per entitiy
1577 size_t *buf = new size_t[is->grid.totalsize()];
1578 send_sizes[cnt] = buf;
1579
1580 // loop over entities and ask for size
1581 int i=0; size_t n=0;
1582 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1584 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1585 itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1586 for ( ; it!=itend; ++it)
1587 {
1588 buf[i] = data.size(*it);
1589 n += buf[i];
1590 i++;
1591 }
1592
1593 // now we know the size for this rank
1594 send_size[cnt] = n;
1595
1596 // hand over send request to torus class
1597 torus().send(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
1598 cnt++;
1599 }
1600
1601 // allocate recv buffers for sizes and store receive request
1602 cnt=0;
1603 for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1604 {
1605 // allocate recv buffer
1606 size_t *buf = new size_t[is->grid.totalsize()];
1607 recv_sizes[cnt] = buf;
1608
1609 // hand over recv request to torus class
1610 torus().recv(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
1611 cnt++;
1612 }
1613
1614 // exchange all size buffers now
1615 torus().exchange();
1616
1617 // release send size buffers
1618 cnt=0;
1619 for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1620 {
1621 delete[] send_sizes[cnt];
1622 send_sizes[cnt] = 0;
1623 cnt++;
1624 }
1625
1626 // process receive size buffers
1627 cnt=0;
1628 for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1629 {
1630 // get recv buffer
1631 size_t *buf = recv_sizes[cnt];
1632
1633 // compute total size
1634 size_t n=0;
1635 for (int i=0; i<is->grid.totalsize(); ++i)
1636 n += buf[i];
1637
1638 // ... and store it
1639 recv_size[cnt] = n;
1640 ++cnt;
1641 }
1642 }
1643
1644
1645 // allocate & fill the send buffers & store send request
1646 std::vector<DataType*> sends(sendlist->size(), static_cast<DataType*>(0)); // store pointers to send buffers
1647 cnt=0;
1648 for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1649 {
1650 // allocate send buffer
1651 DataType *buf = new DataType[send_size[cnt]];
1652
1653 // remember send buffer
1654 sends[cnt] = buf;
1655
1656 // make a message buffer
1657 MessageBuffer<DataType> mb(buf);
1658
1659 // fill send buffer; iterate over cells in intersection
1660 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1662 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1663 itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1664 for ( ; it!=itend; ++it)
1665 data.gather(mb,*it);
1666
1667 // hand over send request to torus class
1668 torus().send(is->rank,buf,send_size[cnt]*sizeof(DataType));
1669 cnt++;
1670 }
1671
1672 // allocate recv buffers and store receive request
1673 std::vector<DataType*> recvs(recvlist->size(),static_cast<DataType*>(0)); // store pointers to send buffers
1674 cnt=0;
1675 for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1676 {
1677 // allocate recv buffer
1678 DataType *buf = new DataType[recv_size[cnt]];
1679
1680 // remember recv buffer
1681 recvs[cnt] = buf;
1682
1683 // hand over recv request to torus class
1684 torus().recv(is->rank,buf,recv_size[cnt]*sizeof(DataType));
1685 cnt++;
1686 }
1687
1688 // exchange all buffers now
1689 torus().exchange();
1690
1691 // release send buffers
1692 cnt=0;
1693 for (ListIt is=sendlist->begin(); is!=sendlist->end(); ++is)
1694 {
1695 delete[] sends[cnt];
1696 sends[cnt] = 0;
1697 cnt++;
1698 }
1699
1700 // process receive buffers and delete them
1701 cnt=0;
1702 for (ListIt is=recvlist->begin(); is!=recvlist->end(); ++is)
1703 {
1704 // get recv buffer
1705 DataType *buf = recvs[cnt];
1706
1707 // make a message buffer
1708 MessageBuffer<DataType> mb(buf);
1709
1710 // copy data from receive buffer; iterate over cells in intersection
1711 if (data.fixedSize(dim,codim))
1712 {
1713 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1715 size_t n=data.size(*it);
1716 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1717 itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1718 for ( ; it!=itend; ++it)
1719 data.scatter(mb,*it,n);
1720 }
1721 else
1722 {
1723 int i=0;
1724 size_t *sbuf = recv_sizes[cnt];
1725 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1727 typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
1728 itend(YaspLevelIterator<codim,All_Partition,GridImp>(g, typename YGrid::Iterator(is->yg,true)));
1729 for ( ; it!=itend; ++it)
1730 data.scatter(mb,*it,sbuf[i++]);
1731 delete[] sbuf;
1732 }
1733
1734 // delete buffer
1735 delete[] buf; // hier krachts !
1736 cnt++;
1737 }
1738 }
1739
1740 // The new index sets from DDM 11.07.2005
1741 const typename Traits::GlobalIdSet& globalIdSet() const
1742 {
1743 return theglobalidset;
1744 }
1745
1746 const typename Traits::LocalIdSet& localIdSet() const
1747 {
1748 return theglobalidset;
1749 }
1750
1751 const typename Traits::LevelIndexSet& levelIndexSet(int level) const
1752 {
1753 if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1754 return *(indexsets[level]);
1755 }
1756
1757 const typename Traits::LeafIndexSet& leafIndexSet() const
1758 {
1759 return leafIndexSet_;
1760 }
1761
1764 const Communication& comm () const
1765 {
1766 return ccobj;
1767 }
1768
1769 private:
1770
1771 // number of boundary segments of the level 0 grid
1772 int nBSegments;
1773
1774 // Index classes need access to the real entity
1775 friend class Dune::YaspIndexSet<const Dune::YaspGrid<dim, Coordinates>, true >;
1776 friend class Dune::YaspIndexSet<const Dune::YaspGrid<dim, Coordinates>, false >;
1777 friend class Dune::YaspGlobalIdSet<const Dune::YaspGrid<dim, Coordinates> >;
1778 friend class Dune::YaspPersistentContainerIndex<const Dune::YaspGrid<dim, Coordinates> >;
1779
1780 friend class Dune::YaspIntersectionIterator<const Dune::YaspGrid<dim, Coordinates> >;
1781 friend class Dune::YaspIntersection<const Dune::YaspGrid<dim, Coordinates> >;
1782 friend class Dune::YaspEntity<0, dim, const Dune::YaspGrid<dim, Coordinates> >;
1783
1784 template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
1785 friend class Entity;
1786
1787 template<class DT>
1788 class MessageBuffer {
1789 public:
1790 // Constructor
1791 MessageBuffer (DT *p)
1792 {
1793 a=p;
1794 i=0;
1795 j=0;
1796 }
1797
1798 // write data to message buffer, acts like a stream !
1799 template<class Y>
1800 void write (const Y& data)
1801 {
1802 static_assert(( std::is_same<DT,Y>::value ), "DataType mismatch");
1803 a[i++] = data;
1804 }
1805
1806 // read data from message buffer, acts like a stream !
1807 template<class Y>
1808 void read (Y& data) const
1809 {
1810 static_assert(( std::is_same<DT,Y>::value ), "DataType mismatch");
1811 data = a[j++];
1812 }
1813
1814 private:
1815 DT *a;
1816 int i;
1817 mutable int j;
1818 };
1819
1821 template<int cd, PartitionIteratorType pitype>
1822 YaspLevelIterator<cd,pitype,GridImp> levelbegin (int level) const
1823 {
1824 YGridLevelIterator g = begin(level);
1825 if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1826
1827 if (pitype==Interior_Partition)
1828 return YaspLevelIterator<cd,pitype,GridImp>(g,g->interior[cd].begin());
1829 if (pitype==InteriorBorder_Partition)
1830 return YaspLevelIterator<cd,pitype,GridImp>(g,g->interiorborder[cd].begin());
1831 if (pitype==Overlap_Partition)
1832 return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlap[cd].begin());
1833 if (pitype<=All_Partition)
1834 return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlapfront[cd].begin());
1835 if (pitype==Ghost_Partition)
1836 return levelend <cd, pitype> (level);
1837
1838 DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
1839 }
1840
1842 template<int cd, PartitionIteratorType pitype>
1843 YaspLevelIterator<cd,pitype,GridImp> levelend (int level) const
1844 {
1845 YGridLevelIterator g = begin(level);
1846 if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
1847
1848 if (pitype==Interior_Partition)
1849 return YaspLevelIterator<cd,pitype,GridImp>(g,g->interior[cd].end());
1850 if (pitype==InteriorBorder_Partition)
1851 return YaspLevelIterator<cd,pitype,GridImp>(g,g->interiorborder[cd].end());
1852 if (pitype==Overlap_Partition)
1853 return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlap[cd].end());
1854 if (pitype<=All_Partition || pitype == Ghost_Partition)
1855 return YaspLevelIterator<cd,pitype,GridImp>(g,g->overlapfront[cd].end());
1856
1857 DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
1858 }
1859
1860 Communication ccobj;
1861
1862 Torus<Communication,dim> _torus;
1863
1864 std::vector< std::shared_ptr< YaspIndexSet<const YaspGrid<dim,Coordinates>, false > > > indexsets;
1865 YaspIndexSet<const YaspGrid<dim,Coordinates>, true> leafIndexSet_;
1866 YaspGlobalIdSet<const YaspGrid<dim,Coordinates> > theglobalidset;
1867
1869 iTupel _s;
1870 std::bitset<dim> _periodic;
1871 iTupel _coarseSize;
1872 ReservedVector<YGridLevel,32> _levels;
1873 int _overlap;
1874 bool keep_ovlp;
1875 int adaptRefCount;
1876 bool adaptActive;
1877 };
1878
1879#if __cpp_deduction_guides >= 201611
1880 // Class template deduction guides
1881 template<typename ctype, int dim>
1882 YaspGrid(FieldVector<ctype, dim>,
1883 std::array<int, std::size_t{dim}>,
1884 std::bitset<std::size_t{dim}> = std::bitset<std::size_t{dim}>{0ULL},
1885 int = 1,
1886 YaspCommunication = YaspCommunication(),
1887 const YLoadBalance<dim>* = YaspGrid< dim, EquidistantCoordinates<ctype, dim> >::defaultLoadbalancer())
1888 -> YaspGrid< dim, EquidistantCoordinates<ctype, dim> >;
1889
1890 template<typename ctype, int dim>
1891 YaspGrid(FieldVector<ctype, dim>,
1892 FieldVector<ctype, dim>,
1893 std::array<int, std::size_t{dim}>,
1894 std::bitset<std::size_t{dim}> = std::bitset<std::size_t{dim}>{0ULL},
1895 int = 1,
1896 YaspCommunication = YaspCommunication(),
1897 const YLoadBalance<dim>* = YaspGrid< dim, EquidistantOffsetCoordinates<ctype, dim> >::defaultLoadbalancer())
1898 -> YaspGrid< dim, EquidistantOffsetCoordinates<ctype, dim> >;
1899
1900 template<typename ctype, std::size_t dim>
1901 YaspGrid(std::array<std::vector<ctype>, dim>,
1902 std::bitset<dim> = std::bitset<dim>{0ULL},
1903 int = 1,
1904 YaspCommunication = YaspCommunication(),
1905 const YLoadBalance<int{dim}>* = YaspGrid< int{dim}, TensorProductCoordinates<ctype, int{dim}> >::defaultLoadbalancer())
1906 -> YaspGrid< int{dim}, TensorProductCoordinates<ctype, int{dim}> >;
1907#endif
1908
1910 template <int d, class CC>
1911 std::ostream& operator<< (std::ostream& s, const YaspGrid<d,CC>& grid)
1912 {
1913 int rank = grid.torus().rank();
1914
1915 s << "[" << rank << "]:" << " YaspGrid maxlevel=" << grid.maxLevel() << std::endl;
1916
1917 s << "Printing the torus: " <<std::endl;
1918 s << grid.torus() << std::endl;
1919
1920 for (typename YaspGrid<d,CC>::YGridLevelIterator g=grid.begin(); g!=grid.end(); ++g)
1921 {
1922 s << "[" << rank << "]: " << std::endl;
1923 s << "[" << rank << "]: " << "==========================================" << std::endl;
1924 s << "[" << rank << "]: " << "level=" << g->level() << std::endl;
1925
1926 for (int codim = 0; codim < d + 1; ++codim)
1927 {
1928 s << "[" << rank << "]: " << "overlapfront[" << codim << "]: " << g->overlapfront[codim] << std::endl;
1929 s << "[" << rank << "]: " << "overlap[" << codim << "]: " << g->overlap[codim] << std::endl;
1930 s << "[" << rank << "]: " << "interiorborder[" << codim << "]: " << g->interiorborder[codim] << std::endl;
1931 s << "[" << rank << "]: " << "interior[" << codim << "]: " << g->interior[codim] << std::endl;
1932
1933 typedef typename YGridList<CC>::Iterator I;
1934 for (I i=g->send_overlapfront_overlapfront[codim].begin();
1935 i!=g->send_overlapfront_overlapfront[codim].end(); ++i)
1936 s << "[" << rank << "]: " << " s_of_of[" << codim << "] to rank "
1937 << i->rank << " " << i->grid << std::endl;
1938
1939 for (I i=g->recv_overlapfront_overlapfront[codim].begin();
1940 i!=g->recv_overlapfront_overlapfront[codim].end(); ++i)
1941 s << "[" << rank << "]: " << " r_of_of[" << codim << "] to rank "
1942 << i->rank << " " << i->grid << std::endl;
1943
1944 for (I i=g->send_overlap_overlapfront[codim].begin();
1945 i!=g->send_overlap_overlapfront[codim].end(); ++i)
1946 s << "[" << rank << "]: " << " s_o_of[" << codim << "] to rank "
1947 << i->rank << " " << i->grid << std::endl;
1948
1949 for (I i=g->recv_overlapfront_overlap[codim].begin();
1950 i!=g->recv_overlapfront_overlap[codim].end(); ++i)
1951 s << "[" << rank << "]: " << " r_of_o[" << codim << "] to rank "
1952 << i->rank << " " << i->grid << std::endl;
1953
1954 for (I i=g->send_interiorborder_interiorborder[codim].begin();
1955 i!=g->send_interiorborder_interiorborder[codim].end(); ++i)
1956 s << "[" << rank << "]: " << " s_ib_ib[" << codim << "] to rank "
1957 << i->rank << " " << i->grid << std::endl;
1958
1959 for (I i=g->recv_interiorborder_interiorborder[codim].begin();
1960 i!=g->recv_interiorborder_interiorborder[codim].end(); ++i)
1961 s << "[" << rank << "]: " << " r_ib_ib[" << codim << "] to rank "
1962 << i->rank << " " << i->grid << std::endl;
1963
1964 for (I i=g->send_interiorborder_overlapfront[codim].begin();
1965 i!=g->send_interiorborder_overlapfront[codim].end(); ++i)
1966 s << "[" << rank << "]: " << " s_ib_of[" << codim << "] to rank "
1967 << i->rank << " " << i->grid << std::endl;
1968
1969 for (I i=g->recv_overlapfront_interiorborder[codim].begin();
1970 i!=g->recv_overlapfront_interiorborder[codim].end(); ++i)
1971 s << "[" << rank << "]: " << " r_of_ib[" << codim << "] to rank "
1972 << i->rank << " " << i->grid << std::endl;
1973 }
1974 }
1975
1976 s << std::endl;
1977
1978 return s;
1979 }
1980
1981 namespace Capabilities
1982 {
1983
1991 template<int dim, class Coordinates>
1992 struct hasBackupRestoreFacilities< YaspGrid<dim, Coordinates> >
1993 {
1994 static const bool v = true;
1995 };
1996
2000 template<int dim, class Coordinates>
2001 struct hasSingleGeometryType< YaspGrid<dim, Coordinates> >
2002 {
2003 static const bool v = true;
2004 static const unsigned int topologyId = GeometryTypes::cube(dim).id();
2005 };
2006
2010 template<int dim, class Coordinates>
2011 struct isCartesian< YaspGrid<dim, Coordinates> >
2012 {
2013 static const bool v = true;
2014 };
2015
2019 template<int dim, class Coordinates, int codim>
2020 struct hasEntity< YaspGrid<dim, Coordinates>, codim>
2021 {
2022 static const bool v = true;
2023 };
2024
2029 template<int dim, class Coordinates, int codim>
2030 struct hasEntityIterator<YaspGrid<dim, Coordinates>, codim>
2031 {
2032 static const bool v = true;
2033 };
2034
2038 template<int dim, int codim, class Coordinates>
2039 struct canCommunicate< YaspGrid< dim, Coordinates>, codim >
2040 {
2041 static const bool v = true;
2042 };
2043
2047 template<int dim, class Coordinates>
2048 struct isLevelwiseConforming< YaspGrid<dim, Coordinates> >
2049 {
2050 static const bool v = true;
2051 };
2052
2056 template<int dim, class Coordinates>
2057 struct isLeafwiseConforming< YaspGrid<dim, Coordinates> >
2058 {
2059 static const bool v = true;
2060 };
2061
2065 template<int dim, class Coordinates>
2066 struct viewThreadSafe< YaspGrid<dim, Coordinates> >
2067 {
2068 static const bool v = true;
2069 };
2070
2071 }
2072
2073} // end namespace
2074
2075// Include the specialization of the StructuredGridFactory class for YaspGrid
2077// Include the specialization of the BackupRestoreFacility class for YaspGrid
2078#include <dune/grid/yaspgrid/backuprestore.hh>
2079
2080#endif
A geometry implementation for axis-aligned hypercubes.
Portable very large unsigned integers.
CommDataHandleIF describes the features of a data handle for communication in parallel runs using the...
Definition: datahandleif.hh:78
Wrapper class for entities.
Definition: entity.hh:66
Container for equidistant coordinates in a YaspGrid.
Definition: coordinates.hh:29
Container for equidistant coordinates in a YaspGrid with non-trivial origin.
Definition: coordinates.hh:131
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:126
constexpr unsigned int dim() const
Return dimension of the type.
Definition: type.hh:372
constexpr bool isCube() const
Return true if entity is a cube of any dimension.
Definition: type.hh:336
constexpr unsigned int id() const
Return the topology id of the type.
Definition: type.hh:377
Definition: grid.hh:862
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:20
detected_or_fallback_t< DeprecatedCollectiveCommunication_t, Communication_t, typename GridFamily::Traits > Communication
A type that is a model of Dune::Communication. It provides a portable way for communication on the se...
Definition: grid.hh:525
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition: intersection.hh:164
A Vector class with statically reserved memory.
Definition: reservedvector.hh:47
constexpr iterator end() noexcept
Returns an iterator pointing to the end of the vector.
Definition: reservedvector.hh:273
constexpr size_type size() const noexcept
Returns number of elements in the vector.
Definition: reservedvector.hh:387
constexpr iterator begin() noexcept
Returns a iterator pointing to the beginning of the vector.
Definition: reservedvector.hh:237
constexpr void pop_back() noexcept
Erases the last element of the vector, O(1) time.
Definition: reservedvector.hh:227
constexpr reference back() noexcept
Returns reference to last element of vector.
Definition: reservedvector.hh:357
constexpr void resize(size_type s) noexcept
Specifies a new size for the vector.
Definition: reservedvector.hh:187
Coordinate container for a tensor product YaspGrid.
Definition: coordinates.hh:245
double partition(int rank, iTupel origin_in, iTupel size_in, iTupel &origin_out, iTupel &size_out) const
partition the given grid onto the torus and return the piece of the process with given rank; returns ...
Definition: torus.hh:239
iTupel coord() const
return own coordinates
Definition: torus.hh:100
int rank() const
return own rank
Definition: torus.hh:94
const iTupel & dims() const
return dimensions of torus
Definition: torus.hh:112
void recv(int rank, void *buffer, int size) const
store a receive request; buffers are received in order; handles also local requests with memcpy
Definition: torus.hh:374
void send(int rank, void *buffer, int size) const
store a send request; buffers are sent in order; handles also local requests with memcpy
Definition: torus.hh:361
int neighbors() const
return the number of neighbors, which is
Definition: torus.hh:203
ProcListIterator recvend() const
last process in receive list
Definition: torus.hh:355
ProcListIterator sendend() const
end of send list
Definition: torus.hh:343
ProcListIterator sendbegin() const
first process in send list
Definition: torus.hh:337
void exchange() const
exchange messages stored in request buffers; clear request buffers afterwards
Definition: torus.hh:387
ProcListIterator recvbegin() const
first process in receive list
Definition: torus.hh:349
Definition: ygrid.hh:75
YGridComponent< Coordinates > move(iTupel v) const
return grid moved by the vector v
Definition: ygrid.hh:263
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:271
implements a collection of multiple std::deque<Intersection> Intersections with neighboring processor...
Definition: ygrid.hh:823
int size() const
return the size of the container, this is the sum of the sizes of all deques
Definition: ygrid.hh:953
Iterator end() const
return iterator pointing to the end of the container
Definition: ygrid.hh:929
Iterator begin() const
return iterator pointing to the begin of the container
Definition: ygrid.hh:923
Iterator over a collection o YGrids A YGrid::Iterator is the heart of an entity in YaspGrid.
Definition: ygrid.hh:594
implements a collection of YGridComponents which form a codimension Entities of given codimension c n...
Definition: ygrid.hh:551
Implement the default load balance strategy of yaspgrid.
Definition: partitioning.hh:206
persistent, globally unique Ids
Definition: yaspgrididset.hh:25
[ provides Dune::Grid ]
Definition: yaspgrid.hh:163
YaspGridFamily< dim, Coordinates > GridFamily
the GridFamily of this grid
Definition: yaspgrid.hh:722
int ghostSize(int level, int codim) const
return size (= distance in graph) of ghost region
Definition: yaspgrid.hh:1420
const Torus< Communication, dim > & torus() const
return reference to torus
Definition: yaspgrid.hh:244
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lbegin(int level) const
version without second template parameter for convenience
Definition: yaspgrid.hh:1350
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: yaspgrid.hh:1452
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity in maxLevel
Definition: yaspgrid.hh:1378
void globalRefine(int refCount)
refine the grid refCount times.
Definition: yaspgrid.hh:1223
int getMark(const typename Traits::template Codim< 0 >::Entity &e) const
returns adaptation mark for given entity
Definition: yaspgrid.hh:1307
YaspGrid(std::array< std::vector< ctype >, std::size_t{dim}> coords, std::bitset< std::size_t{dim}> periodic=std::bitset< std::size_t{dim}>(0ULL), int overlap=1, Communication comm=Communication(), const Yasp::Partitioning< dim > *partitioner=defaultPartitioner())
Standard constructor for a tensorproduct YaspGrid.
Definition: yaspgrid.hh:1050
int ghostSize(int codim) const
return size (= distance in graph) of ghost region
Definition: yaspgrid.hh:1426
int globalSize(int i) const
return number of cells on finest level in given direction on all processors
Definition: yaspgrid.hh:250
YaspGrid(Dune::FieldVector< ctype, dim > lowerleft, Dune::FieldVector< ctype, dim > upperright, std::array< int, std::size_t{dim}> s, std::bitset< std::size_t{dim}> periodic=std::bitset< std::size_t{dim}>(0ULL), int overlap=1, Communication comm=Communication(), const Yasp::Partitioning< dim > *partitioner=defaultPartitioner())
Definition: yaspgrid.hh:981
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafend() const
return LeafIterator which points behind the last entity in maxLevel
Definition: yaspgrid.hh:1371
void postAdapt()
clean up some markers
Definition: yaspgrid.hh:1328
const Dune::FieldVector< ctype, dim > & domainSize() const
returns the size of the physical domain
Definition: yaspgrid.hh:1470
YGridLevelIterator end() const
return iterator pointing to one past the finest level
Definition: yaspgrid.hh:305
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: yaspgrid.hh:1458
int maxLevel() const
Definition: yaspgrid.hh:1217
YaspGrid(Dune::FieldVector< ctype, dim > L, std::array< int, std::size_t{dim}> s, std::bitset< std::size_t{dim}> periodic=std::bitset< std::size_t{dim}>{0ULL}, int overlap=1, Communication comm=Communication(), const Yasp::Partitioning< dim > *partitioner=defaultPartitioner())
Definition: yaspgrid.hh:911
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafend() const
return LeafIterator which points behind the last entity in maxLevel
Definition: yaspgrid.hh:1385
void communicate(CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir, int level) const
Definition: yaspgrid.hh:1479
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:568
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity in maxLevel
Definition: yaspgrid.hh:1364
bool preAdapt()
returns true, if the grid will be coarsened
Definition: yaspgrid.hh:1320
iTupel levelSize(int l) const
return size vector of the grid (in cells) on level l
Definition: yaspgrid.hh:268
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:1293
int overlapSize(int level, int codim) const
return size (= distance in graph) of overlap region
Definition: yaspgrid.hh:1406
void communicate(CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir) const
Definition: yaspgrid.hh:1489
int size(int codim) const
number of leaf entities per codim in this process
Definition: yaspgrid.hh:1446
bool isPeriodic(int i) const
return whether the grid is periodic in direction i
Definition: yaspgrid.hh:277
void refineOptions(bool keepPhysicalOverlap)
set options for refinement
Definition: yaspgrid.hh:1277
int size(int level, int codim) const
number of entities per level and codim in this process
Definition: yaspgrid.hh:1432
YaspGrid(const Coordinates &coordinates, std::bitset< dim > periodic=std::bitset< dim >(0ULL), int overlap=1, Communication comm=Communication(), const Yasp::Partitioning< dim > *partitioner=defaultPartitioner())
Definition: yaspgrid.hh:738
bool adapt()
map adapt to global refine
Definition: yaspgrid.hh:1313
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lend(int level) const
version without second template parameter for convenience
Definition: yaspgrid.hh:1357
ReservedVector< YGridLevel, 32 >::const_iterator YGridLevelIterator
Iterator over the grid levels.
Definition: yaspgrid.hh:288
size_t numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition: yaspgrid.hh:1464
void communicateCodim(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int level) const
Definition: yaspgrid.hh:1499
void makelevel(const Coordinates &coords, std::bitset< dim > periodic, iTupel o_interior, int overlap)
Make a new YGridLevel structure.
Definition: yaspgrid.hh:333
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:1343
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator lbegin(int level) const
one past the end on this level
Definition: yaspgrid.hh:1336
iTupel globalSize() const
return number of cells on finest level on all processors
Definition: yaspgrid.hh:256
const Communication & comm() const
return a communication object
Definition: yaspgrid.hh:1764
int overlapSize(int odim) const
return size (= distance in graph) of overlap region
Definition: yaspgrid.hh:1413
YGridLevelIterator begin(int i) const
return iterator pointing to given level
Definition: yaspgrid.hh:297
Coordinates::ctype ctype
Type used for coordinates.
Definition: yaspgrid.hh:176
int levelSize(int l, int i) const
return size of the grid (in cells) on level l in direction i
Definition: yaspgrid.hh:262
YGridLevelIterator begin() const
return iterator pointing to coarsest level
Definition: yaspgrid.hh:291
YaspHierarchicIterator enables iteration over son entities of codim 0.
Definition: yaspgridhierarchiciterator.hh:20
Implementation of Level- and LeafIndexSets for YaspGrid.
Definition: yaspgridindexsets.hh:25
YaspIntersectionIterator enables iteration over intersections with neighboring codim 0 entities.
Definition: yaspgridintersectioniterator.hh:22
YaspIntersection provides data about intersection with neighboring codim 0 entities.
Definition: yaspgridintersection.hh:22
Iterates over entities of one grid level.
Definition: yaspgridleveliterator.hh:19
implement a consecutive index for all entities of given codim of a YaspGrid
Definition: yaspgridpersistentcontainer.hh:35
a base class for the yaspgrid partitioning strategy
Definition: partitioning.hh:39
A set of traits classes to store static information about grid implementation.
Different resources needed by all grid implementations.
Implements an utility class that provides collective communication methods for sequential programs.
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.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
CommunicationDirection
Define a type for communication direction parameter.
Definition: gridenums.hh:170
InterfaceType
Parameter to be used for the communication functions.
Definition: gridenums.hh:86
@ All_Partition
all entities
Definition: gridenums.hh:141
@ Interior_Partition
only interior entities
Definition: gridenums.hh:137
@ InteriorBorder_Partition
interior and border entities
Definition: gridenums.hh:138
@ Overlap_Partition
interior, border, and overlap entities
Definition: gridenums.hh:139
@ Ghost_Partition
only ghost entities
Definition: gridenums.hh:142
@ BackwardCommunication
reverse communication direction
Definition: gridenums.hh:172
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition: gridenums.hh:88
@ All_All_Interface
send all and receive all entities
Definition: gridenums.hh:91
@ Overlap_All_Interface
send overlap, receive all entities
Definition: gridenums.hh:90
@ Overlap_OverlapFront_Interface
send overlap, receive overlap and front entities
Definition: gridenums.hh:89
@ InteriorBorder_InteriorBorder_Interface
send/receive interior and border entities
Definition: gridenums.hh:87
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:473
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
Provides base classes for index and id sets.
Some useful basic math stuff.
Implements an utility class that provides MPI's collective communication methods.
Helpers for dealing with MPI.
constexpr Front front
PartitionSet for the front partition.
Definition: partitionset.hh:281
Dune namespace.
Definition: alignedallocator.hh:13
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
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:97
specialize with 'true' for all codims that a grid provides an iterator for (default=hasEntity<codim>:...
Definition: capabilities.hh:74
Specialize with 'true' for all codims that a grid implements entities for. (default=false)
Definition: capabilities.hh:58
Specialize with 'true' for if the codimension 0 entity of the grid has only one possible geometry typ...
Definition: capabilities.hh:27
Specialize with 'true' if the grid is a Cartesian grid. Cartesian grids satisfy the following propert...
Definition: capabilities.hh:48
Specialize with 'true' if implementation guarantees a conforming leaf grid. (default=false)
Definition: capabilities.hh:115
Specialize with 'true' if implementation guarantees conforming level grids. (default=false)
Definition: capabilities.hh:106
Specialize with 'true' if the grid implementation is thread safe, while it is not modified....
Definition: capabilities.hh:169
Static tag representing a codimension.
Definition: dimension.hh:24
A traits struct that collects all associated types of one grid model.
Definition: grid.hh:995
IdSet< const GridImp, LocalIdSetImp, LIDType > LocalIdSet
The type of the local id set.
Definition: grid.hh:1068
IdSet< const GridImp, GlobalIdSetImp, GIDType > GlobalIdSet
The type of the global id set.
Definition: grid.hh:1066
IndexSet< const GridImp, LeafIndexSetImp > LeafIndexSet
The type of the leaf index set.
Definition: grid.hh:1064
IndexSet< const GridImp, LevelIndexSetImp > LevelIndexSet
The type of the level index set.
Definition: grid.hh:1062
type describing an intersection with a neighboring processor
Definition: ygrid.hh:829
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.111.3 (Nov 24, 23:30, 2024)