DUNE PDELab (2.7)

genericdatahandle.hh
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_PDELAB_GRIDFUNCTIONSPACE_GENERICDATAHANDLE_HH
4#define DUNE_PDELAB_GRIDFUNCTIONSPACE_GENERICDATAHANDLE_HH
5
6#include <vector>
7#include <set>
8#include <limits>
9
12
14#include <dune/grid/common/gridenums.hh>
15
16#include <dune/pdelab/common/polymorphicbufferwrapper.hh>
17#include <dune/pdelab/gridfunctionspace/entityindexcache.hh>
18
19namespace Dune {
20 namespace PDELab {
21
23 template<typename E, bool transmit_rank = false>
25 {
26
27 typedef char DataType;
28
30 typedef std::size_t size_type;
31
32 // Wrap the grid's communication buffer to enable sending leaf ordering sizes along with the data
33 static const bool wrap_buffer = true;
34
35 // Determine whether the data handle should send the MPI rank to the receiver
36 static constexpr bool transmitRank()
37 {
38 return transmit_rank;
39 }
40
41 // export original data type to fix up size information forwarded to standard gather / scatter functors
42 typedef E OriginalDataType;
43
44 template<typename GFS>
45 bool contains(const GFS& gfs, int dim, int codim) const
46 {
47 return gfs.dataHandleContains(codim);
48 }
49
50 template<typename GFS>
51 bool fixedSize(const GFS& gfs, int dim, int codim) const
52 {
53 return gfs.dataHandleFixedSize(codim);
54 }
55
56 template<typename GFS, typename Entity>
57 std::size_t size(const GFS& gfs, const Entity& e) const
58 {
59 // include size of leaf ordering offsets and of the rank if necessary
60 return gfs.dataHandleSize(e) * sizeof(E)
61 + (gfs.sendLeafSizes() ? TypeTree::TreeInfo<typename GFS::Ordering>::leafCount * sizeof(size_type) : 0)
62 + (transmitRank() ? sizeof(int) : 0);
63 }
64
65 };
66
68 template<typename E, bool transmit_rank = false>
70 {
71
73 typedef std::size_t size_type;
74
75 // Data is per entity, so we don' need to send leaf ordering size and thus can avoid wrapping the
76 // grid's communication buffer unless we have to send our rank to the receiver
77 static const bool wrap_buffer = transmit_rank;
78
79 // Determine whether the data handle should send the MPI rank to the receiver
80 static constexpr bool transmitRank()
81 {
82 return transmit_rank;
83 }
84
85 // export original data type to fix up size information forwarded to standard gather / scatter functors
86 typedef E OriginalDataType;
87
88 using DataType = std::conditional_t<wrap_buffer,char,E>;
89
90 template<typename GFS>
91 bool contains(const GFS& gfs, int dim, int codim) const
92 {
93 return gfs.dataHandleContains(codim);
94 }
95
96 template<typename GFS>
97 bool fixedSize(const GFS& gfs, int dim, int codim) const
98 {
99 return true;
100 }
101
102 template<typename GFS, typename Entity>
103 std::size_t size(const GFS& gfs, const Entity& e) const
104 {
105 if (wrap_buffer)
106 return (gfs.dataHandleContains(Entity::codimension) and gfs.entitySet().contains(e) ? _count * sizeof(E) : 0)
107 + (transmitRank() ? sizeof(int) : 0);
108 else
109 return gfs.dataHandleContains(Entity::codimension) && gfs.entitySet().contains(e) ? _count : 0;
110 }
111
112 explicit EntityDataCommunicationDescriptor(std::size_t count = 1)
113 : _count(count)
114 {}
115
116 private:
117
118 const std::size_t _count;
119
120 };
121
123
129 template<typename GFS, typename V, typename GatherScatter, typename CommunicationDescriptor = DOFDataCommunicationDescriptor<typename V::ElementType> >
131 : public Dune::CommDataHandleIF<GFSDataHandle<GFS,V,GatherScatter,CommunicationDescriptor>,typename CommunicationDescriptor::DataType>
132 {
133
134 public:
135
136 typedef typename CommunicationDescriptor::DataType DataType;
137 typedef typename GFS::Traits::SizeType size_type;
138
139 static const size_type leaf_count = TypeTree::TreeInfo<typename GFS::Ordering>::leafCount;
140
141 GFSDataHandle(const GFS& gfs, V& v, GatherScatter gather_scatter = GatherScatter(), CommunicationDescriptor communication_descriptor = CommunicationDescriptor())
142 : _gfs(gfs)
143 , _index_cache(gfs)
144 , _local_view(v)
145 , _gather_scatter(gather_scatter)
146 , _communication_descriptor(communication_descriptor)
147 , _rank(gfs.gridView().comm().rank())
148 {}
149
151 bool contains(int dim, int codim) const
152 {
153 return _communication_descriptor.contains(_gfs,dim,codim);
154 }
155
157 bool fixedSize(int dim, int codim) const
158 {
159 return _communication_descriptor.fixedSize(_gfs,dim,codim);
160 }
161
166 template<typename Entity>
167 size_type size(const Entity& e) const
168 {
169 return _communication_descriptor.size(_gfs,e);
170 }
171
173 template<typename MessageBuffer, typename Entity>
174 typename std::enable_if<
175 CommunicationDescriptor::wrap_buffer && AlwaysTrue<Entity>::value // we can only support this if the buffer is wrapped
176 >::type
177 gather(MessageBuffer& buff, const Entity& e) const
178 {
180 buff,
182 _rank,
183 _communication_descriptor.transmitRank()
184 );
185 _index_cache.update(e);
186 _local_view.bind(_index_cache);
187 if (_gfs.sendLeafSizes())
188 {
189 // send the leaf ordering offsets as exported by the EntityIndexCache
190 for (auto it = _index_cache.offsets().begin() + 1,
191 end_it = _index_cache.offsets().end();
192 it != end_it;
193 ++it)
194 {
195 buf_wrapper.write(static_cast<typename CommunicationDescriptor::size_type>(*it));
196 }
197 }
198 // send the normal data
199 if (_gather_scatter.gather(buf_wrapper,e,_local_view))
200 _local_view.commit();
201 _local_view.unbind();
202 }
203
205 template<typename MessageBuffer, typename Entity>
206 typename std::enable_if<
207 !CommunicationDescriptor::wrap_buffer && AlwaysTrue<Entity>::value
208 >::type
209 gather(MessageBuffer& buff, const Entity& e) const
210 {
211 _index_cache.update(e);
212 _local_view.bind(_index_cache);
213 if (_gather_scatter.gather(buff,e,_local_view))
214 _local_view.commit();
215 _local_view.unbind();
216 }
217
224 template<typename MessageBuffer, typename Entity>
225 typename std::enable_if<
226 CommunicationDescriptor::wrap_buffer && AlwaysTrue<Entity>::value // we require the buffer to be wrapped
227 >::type
228 scatter(MessageBuffer& buff, const Entity& e, size_type n)
229 {
231 buff,
233 _rank,
234 _communication_descriptor.transmitRank()
235 );
236 _index_cache.update(e);
237 _local_view.bind(_index_cache);
238 bool needs_commit = false;
239 if (_gfs.sendLeafSizes())
240 {
241 // receive leaf ordering offsets and store in local array
242 typename IndexCache::Offsets remote_offsets = {{0}};
243 for (auto it = remote_offsets.begin() + 1,
244 end_it = remote_offsets.end();
245 it != end_it;
246 ++it)
247 {
248 typename CommunicationDescriptor::size_type data = 0;
249 buf_wrapper.read(data);
250 *it = data;
251 }
252 // call special version of scatter() that can handle empty leafs in the ordering tree
253 needs_commit = _gather_scatter.scatter(buf_wrapper,remote_offsets,_index_cache.offsets(),e,_local_view);
254 }
255 else
256 {
257 // call standard version of scatter - make sure to fix the reported communication size
258 size_type size = n;
259 if (_communication_descriptor.transmitRank())
260 {
261 // this takes care of the rank transmission done inside the buffer
262 size -= sizeof(int);
263 }
264 size /= sizeof(typename CommunicationDescriptor::OriginalDataType);
265 needs_commit = _gather_scatter.scatter(buf_wrapper,size,e,_local_view);
266 }
267
268 if (needs_commit)
269 _local_view.commit();
270
271 _local_view.unbind();
272 }
273
280 template<typename MessageBuffer, typename Entity>
281 typename std::enable_if<
282 !CommunicationDescriptor::wrap_buffer && AlwaysTrue<Entity>::value
283 >::type
284 scatter(MessageBuffer& buff, const Entity& e, size_type n)
285 {
286 _index_cache.update(e);
287 _local_view.bind(_index_cache);
288
289 if (_gather_scatter.scatter(buff,n,e,_local_view))
290 _local_view.commit();
291
292 _local_view.unbind();
293 }
294
295
296 private:
297
298 typedef EntityIndexCache<GFS> IndexCache;
299 typedef typename V::template LocalView<IndexCache> LocalView;
300
301 const GFS& _gfs;
302 mutable IndexCache _index_cache;
303 mutable LocalView _local_view;
304 mutable GatherScatter _gather_scatter;
305 CommunicationDescriptor _communication_descriptor;
306 int _rank;
307
308 };
309
310
311 template<typename GatherScatter>
312 class DataGatherScatter
313 {
314
315 public:
316
317 typedef std::size_t size_type;
318
319 template<typename MessageBuffer, typename Entity, typename LocalView>
320 bool gather(MessageBuffer& buff, const Entity& e, const LocalView& local_view) const
321 {
322 for (std::size_t i = 0; i < local_view.size(); ++i)
323 _gather_scatter.gather(buff,local_view[i]);
324 return false;
325 }
326
327 // default scatter - requires function space structure to be identical on sender and receiver side
328 template<typename MessageBuffer, typename Entity, typename LocalView>
329 bool scatter(MessageBuffer& buff, size_type n, const Entity& e, LocalView& local_view) const
330 {
331 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
332 {
333 if (local_view.size() != n)
334 DUNE_THROW(Exception,"size mismatch in GridFunctionSpace data handle, have " << local_view.size() << "DOFs, but received " << n);
335
336 for (std::size_t i = 0; i < local_view.size(); ++i)
337 _gather_scatter.scatter(buff,local_view[i]);
338 return true;
339 }
340 else
341 {
342 if (local_view.size() != 0)
343 DUNE_THROW(Exception,"expected no DOFs in partition '" << e.partitionType() << "', but have " << local_view.size());
344
345 for (std::size_t i = 0; i < local_view.size(); ++i)
346 {
347 typename LocalView::ElementType dummy;
348 buff.read(dummy);
349 }
350 return false;
351 }
352 }
353
354 // enhanced scatter with support for function spaces with different structure on sender and receiver side
355 template<typename MessageBuffer, typename Offsets, typename Entity, typename LocalView>
356 bool scatter(MessageBuffer& buff, const Offsets& remote_offsets, const Offsets& local_offsets, const Entity& e, LocalView& local_view) const
357 {
358 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
359 {
360 // the idea here is this:
361 // the compile time structure of the overall function space (and its ordering) will be identical on both sides
362 // of the communication, but one side may be missing information on some leaf orderings, e.g. because the DOF
363 // belongs to a MultiDomain subdomain that only has an active grid cell on one side of the communication.
364 // So we step through the leaves and simply ignore any block where one of the two sides is of size 0.
365 // Otherwise, it's business as usual: we make sure that the sizes from both sides match and then process all
366 // data with the DOF-local gather / scatter functor.
367 size_type remote_i = 0;
368 size_type local_i = 0;
369 bool needs_commit = false;
370 for (size_type block = 1; block < local_offsets.size(); ++block)
371 {
372 // we didn't get any data - just ignore
373 if (remote_offsets[block] == remote_i)
374 {
375 local_i = local_offsets[block];
376 continue;
377 }
378
379 // we got data for DOFs we don't have - drain buffer
380 if (local_offsets[block] == local_i)
381 {
382 for (; remote_i < remote_offsets[block]; ++remote_i)
383 {
384 typename LocalView::ElementType dummy;
385 buff.read(dummy);
386 }
387 continue;
388 }
389
390 if (remote_offsets[block] - remote_i != local_offsets[block] - local_i)
391 DUNE_THROW(Exception,"size mismatch in GridFunctionSpace data handle block " << block << ", have " << local_offsets[block] - local_i << "DOFs, but received " << remote_offsets[block] - remote_i);
392
393 for (; local_i < local_offsets[block]; ++local_i)
394 _gather_scatter.scatter(buff,local_view[local_i]);
395
396 remote_i = remote_offsets[block];
397 needs_commit = true;
398 }
399 return needs_commit;
400 }
401 else
402 {
403 if (local_view.size() != 0)
404 DUNE_THROW(Exception,"expected no DOFs in partition '" << e.partitionType() << "', but have " << local_view.size());
405
406 for (std::size_t i = 0; i < remote_offsets.back(); ++i)
407 {
408 typename LocalView::ElementType dummy;
409 buff.read(dummy);
410 }
411 return false;
412 }
413 }
414
415 DataGatherScatter(GatherScatter gather_scatter = GatherScatter())
416 : _gather_scatter(gather_scatter)
417 {}
418
419 private:
420
421 GatherScatter _gather_scatter;
422
423 };
424
425
426 template<typename GatherScatter>
427 class DataEntityGatherScatter
428 {
429
430 public:
431
432 typedef std::size_t size_type;
433
434 template<typename MessageBuffer, typename Entity, typename LocalView>
435 bool gather(MessageBuffer& buff, const Entity& e, const LocalView& local_view) const
436 {
437 for (std::size_t i = 0; i < local_view.size(); ++i)
438 _gather_scatter.gather(buff,e,local_view[i]);
439 return false;
440 }
441
442 // see documentation in DataGatherScatter for further info on the scatter() implementations
443 template<typename MessageBuffer, typename Entity, typename LocalView>
444 bool scatter(MessageBuffer& buff, size_type n, const Entity& e, LocalView& local_view) const
445 {
446 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
447 {
448 if (local_view.size() != n)
449 DUNE_THROW(Exception,"size mismatch in GridFunctionSpace data handle, have " << local_view.size() << "DOFs, but received " << n);
450
451 for (std::size_t i = 0; i < local_view.size(); ++i)
452 _gather_scatter.scatter(buff,e,local_view[i]);
453 return true;
454 }
455 else
456 {
457 if (local_view.size() != 0)
458 DUNE_THROW(Exception,"expected no DOFs in partition '" << e.partitionType() << "', but have " << local_view.size());
459
460 for (std::size_t i = 0; i < local_view.size(); ++i)
461 {
462 typename LocalView::ElementType dummy;
463 buff.read(dummy);
464 }
465 return false;
466 }
467 }
468
469 // see documentation in DataGatherScatter for further info on the scatter() implementations
470 template<typename MessageBuffer, typename Offsets, typename Entity, typename LocalView>
471 bool scatter(MessageBuffer& buff, const Offsets& remote_offsets, const Offsets& local_offsets, const Entity& e, LocalView& local_view) const
472 {
473 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
474 {
475 size_type remote_i = 0;
476 size_type local_i = 0;
477 bool needs_commit = false;
478 for (size_type block = 1; block < local_offsets.size(); ++block)
479 {
480
481 // we didn't get any data - just ignore
482 if (remote_offsets[block] == remote_i)
483 {
484 local_i = local_offsets[block];
485 continue;
486 }
487
488 // we got data for DOFs we don't have - drain buffer
489 if (local_offsets[block] == local_i)
490 {
491 for (; remote_i < remote_offsets[block]; ++remote_i)
492 {
493 typename LocalView::ElementType dummy;
494 buff.read(dummy);
495 }
496 continue;
497 }
498
499 if (remote_offsets[block] - remote_i != local_offsets[block] - local_i)
500 DUNE_THROW(Exception,"size mismatch in GridFunctionSpace data handle block " << block << ", have " << local_offsets[block] - local_i << "DOFs, but received " << remote_offsets[block] - remote_i);
501
502 for (; local_i < local_offsets[block]; ++local_i)
503 _gather_scatter.scatter(buff,e,local_view[local_i]);
504
505 remote_i = remote_offsets[block];
506 needs_commit = true;
507 }
508 return needs_commit;
509 }
510 else
511 {
512 if (local_view.size() != 0)
513 DUNE_THROW(Exception,"expected no DOFs in partition '" << e.partitionType() << "', but have " << local_view.size());
514
515 for (std::size_t i = 0; i < remote_offsets.back(); ++i)
516 {
517 typename LocalView::ElementType dummy;
518 buff.read(dummy);
519 }
520 return false;
521 }
522 }
523
524 DataEntityGatherScatter(GatherScatter gather_scatter = GatherScatter())
525 : _gather_scatter(gather_scatter)
526 {}
527
528 private:
529
530 GatherScatter _gather_scatter;
531
532 };
533
534
535 template<typename GatherScatter>
536 class DataContainerIndexGatherScatter
537 {
538
539 public:
540
541 typedef std::size_t size_type;
542
543 template<typename MessageBuffer, typename Entity, typename LocalView>
544 bool gather(MessageBuffer& buff, const Entity& e, const LocalView& local_view) const
545 {
546 for (std::size_t i = 0; i < local_view.size(); ++i)
547 _gather_scatter.gather(buff,local_view.cache().containerIndex(i),local_view[i]);
548 return false;
549 }
550
551 // see documentation in DataGatherScatter for further info on the scatter() implementations
552 template<typename MessageBuffer, typename Entity, typename LocalView>
553 bool scatter(MessageBuffer& buff, size_type n, const Entity& e, LocalView& local_view) const
554 {
555 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
556 {
557 if (local_view.size() != n)
558 DUNE_THROW(Exception,"size mismatch in GridFunctionSpace data handle, have " << local_view.size() << "DOFs, but received " << n);
559
560 for (std::size_t i = 0; i < local_view.size(); ++i)
561 _gather_scatter.scatter(buff,local_view.cache().containerIndex(i),local_view[i]);
562
563 return true;
564 }
565 else
566 {
567 if (local_view.size() != 0)
568 DUNE_THROW(Exception,"expected no DOFs in partition '" << e.partitionType() << "', but have " << local_view.size());
569
570 for (std::size_t i = 0; i < local_view.size(); ++i)
571 {
572 typename LocalView::ElementType dummy;
573 buff.read(dummy);
574 }
575 return false;
576 }
577 }
578
579 // see documentation in DataGatherScatter for further info on the scatter() implementations
580 template<typename MessageBuffer, typename Offsets, typename Entity, typename LocalView>
581 bool scatter(MessageBuffer& buff, const Offsets& remote_offsets, const Offsets& local_offsets, const Entity& e, LocalView& local_view) const
582 {
583 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
584 {
585 size_type remote_i = 0;
586 size_type local_i = 0;
587 bool needs_commit = false;
588 for (size_type block = 1; block < local_offsets.size(); ++block)
589 {
590
591 // we didn't get any data - just ignore
592 if (remote_offsets[block] == remote_i)
593 {
594 local_i = local_offsets[block];
595 continue;
596 }
597
598 // we got data for DOFs we don't have - drain buffer
599 if (local_offsets[block] == local_i)
600 {
601 for (; remote_i < remote_offsets[block]; ++remote_i)
602 {
603 typename LocalView::ElementType dummy;
604 buff.read(dummy);
605 }
606 continue;
607 }
608
609 if (remote_offsets[block] - remote_i != local_offsets[block] - local_i)
610 DUNE_THROW(Exception,"size mismatch in GridFunctionSpace data handle block " << block << ", have " << local_offsets[block] - local_i << "DOFs, but received " << remote_offsets[block] - remote_i);
611
612 for (; local_i < local_offsets[block]; ++local_i)
613 _gather_scatter.scatter(buff,local_view.cache().containerIndex(local_i),local_view[local_i]);
614
615 remote_i = remote_offsets[block];
616 needs_commit = true;
617 }
618 return needs_commit;
619 }
620 else
621 {
622 if (local_view.size() != 0)
623 DUNE_THROW(Exception,"expected no DOFs in partition '" << e.partitionType() << "', but have " << local_view.size());
624
625 for (std::size_t i = 0; i < remote_offsets.back(); ++i)
626 {
627 typename LocalView::ElementType dummy;
628 buff.read(dummy);
629 }
630 return false;
631 }
632 }
633
634
635 DataContainerIndexGatherScatter(GatherScatter gather_scatter = GatherScatter())
636 : _gather_scatter(gather_scatter)
637 {}
638
639 private:
640
641 GatherScatter _gather_scatter;
642
643 };
644
645
646 class AddGatherScatter
647 {
648 public:
649 template<class MessageBuffer, class DataType>
650 void gather (MessageBuffer& buff, DataType& data) const
651 {
652 buff.write(data);
653 }
654
655 template<class MessageBuffer, class DataType>
656 void scatter (MessageBuffer& buff, DataType& data) const
657 {
658 DataType x;
659 buff.read(x);
660 data += x;
661 }
662 };
663
664 template<class GFS, class V>
665 class AddDataHandle
666 : public GFSDataHandle<GFS,V,DataGatherScatter<AddGatherScatter> >
667 {
668 typedef GFSDataHandle<GFS,V,DataGatherScatter<AddGatherScatter> > BaseT;
669
670 public:
671
672 AddDataHandle (const GFS& gfs_, V& v_)
673 : BaseT(gfs_,v_)
674 {}
675 };
676
677 class AddClearGatherScatter
678 {
679 public:
680 template<class MessageBuffer, class DataType>
681 void gather (MessageBuffer& buff, DataType& data) const
682 {
683 buff.write(data);
684 data = (DataType) 0;
685 }
686
687 template<class MessageBuffer, class DataType>
688 void scatter (MessageBuffer& buff, DataType& data) const
689 {
690 DataType x;
691 buff.read(x);
692 data += x;
693 }
694 };
695
696 template<class GFS, class V>
697 class AddClearDataHandle
698 : public GFSDataHandle<GFS,V,DataGatherScatter<AddClearGatherScatter> >
699 {
700 typedef GFSDataHandle<GFS,V,DataGatherScatter<AddClearGatherScatter> > BaseT;
701
702 public:
703
704 AddClearDataHandle (const GFS& gfs_, V& v_)
705 : BaseT(gfs_,v_)
706 {}
707 };
708
709 class CopyGatherScatter
710 {
711 public:
712 template<class MessageBuffer, class DataType>
713 void gather (MessageBuffer& buff, DataType& data) const
714 {
715 buff.write(data);
716 }
717
718 template<class MessageBuffer, class DataType>
719 void scatter (MessageBuffer& buff, DataType& data) const
720 {
721 DataType x;
722 buff.read(x);
723 data = x;
724 }
725 };
726
727 template<class GFS, class V>
728 class CopyDataHandle
729 : public GFSDataHandle<GFS,V,DataGatherScatter<CopyGatherScatter> >
730 {
731 typedef GFSDataHandle<GFS,V,DataGatherScatter<CopyGatherScatter> > BaseT;
732
733 public:
734
735 CopyDataHandle (const GFS& gfs_, V& v_)
736 : BaseT(gfs_,v_)
737 {}
738 };
739
740 class MinGatherScatter
741 {
742 public:
743 template<class MessageBuffer, class DataType>
744 void gather (MessageBuffer& buff, DataType& data) const
745 {
746 buff.write(data);
747 }
748
749 template<class MessageBuffer, class DataType>
750 void scatter (MessageBuffer& buff, DataType& data) const
751 {
752 DataType x;
753 buff.read(x);
754 data = std::min(data,x);
755 }
756 };
757
758 template<class GFS, class V>
759 class MinDataHandle
760 : public GFSDataHandle<GFS,V,DataGatherScatter<MinGatherScatter> >
761 {
762 typedef GFSDataHandle<GFS,V,DataGatherScatter<MinGatherScatter> > BaseT;
763
764 public:
765
766 MinDataHandle (const GFS& gfs_, V& v_)
767 : BaseT(gfs_,v_)
768 {}
769 };
770
771 class MaxGatherScatter
772 {
773 public:
774 template<class MessageBuffer, class DataType>
775 void gather (MessageBuffer& buff, DataType& data) const
776 {
777 buff.write(data);
778 }
779
780 template<class MessageBuffer, class DataType>
781 void scatter (MessageBuffer& buff, DataType& data) const
782 {
783 DataType x;
784 buff.read(x);
785 data = std::max(data,x);
786 }
787 };
788
789 template<class GFS, class V>
790 class MaxDataHandle
791 : public GFSDataHandle<GFS,V,DataGatherScatter<MaxGatherScatter> >
792 {
793 typedef GFSDataHandle<GFS,V,DataGatherScatter<MaxGatherScatter> > BaseT;
794
795 public:
796
797 MaxDataHandle (const GFS& gfs_, V& v_)
798 : BaseT(gfs_,v_)
799 {}
800 };
801
802
804
812 {
813 public:
814
815 template<typename MessageBuffer, typename Entity, typename LocalView>
816 bool gather(MessageBuffer& buff, const Entity& e, LocalView& local_view) const
817 {
818 // Figure out where we are...
819 const bool ghost = e.partitionType()!=Dune::InteriorEntity && e.partitionType()!=Dune::BorderEntity;
820
821 // ... and send something (doesn't really matter what, we'll throw it away on the receiving side).
822 buff.write(ghost);
823
824 return false;
825 }
826
827 template<typename MessageBuffer, typename Entity, typename LocalView>
828 bool scatter(MessageBuffer& buff, std::size_t n, const Entity& e, LocalView& local_view) const
829 {
830 // Figure out where we are - we have to do this again on the receiving side due to the asymmetric
831 // communication interface!
832 const bool ghost = e.partitionType()!=Dune::InteriorEntity && e.partitionType()!=Dune::BorderEntity;
833
834 // drain buffer
835 bool dummy;
836 buff.read(dummy);
837
838 for (std::size_t i = 0; i < local_view.size(); ++i)
839 local_view[i] = ghost;
840
841 return true;
842 }
843
844 };
845
847
854 template<class GFS, class V>
856 : public Dune::PDELab::GFSDataHandle<GFS,
857 V,
858 GhostGatherScatter,
859 EntityDataCommunicationDescriptor<bool> >
860 {
862 GFS,
863 V,
866 > BaseT;
867
868 static_assert((std::is_same<typename V::ElementType,bool>::value),
869 "GhostDataHandle expects a vector of bool values");
870
871 public:
872
874
883 GhostDataHandle(const GFS& gfs_, V& v_, bool init_vector = true)
884 : BaseT(gfs_,v_)
885 {
886 if (init_vector)
887 v_ = false;
888 }
889 };
890
891
893
901 template<typename RankIndex>
903 {
904
905 public:
906
907 template<typename MessageBuffer, typename Entity, typename LocalView>
908 bool gather(MessageBuffer& buff, const Entity& e, LocalView& local_view) const
909 {
910 // We only gather from interior and border entities, so we can throw in our ownership
911 // claim without any further checks.
912 buff.write(_rank);
913
914 return true;
915 }
916
917 template<typename MessageBuffer, typename Entity, typename LocalView>
918 bool scatter(MessageBuffer& buff, std::size_t n, const Entity& e, LocalView& local_view) const
919 {
920 // Value used for DOFs with currently unknown rank.
921 const RankIndex unknown_rank = std::numeric_limits<RankIndex>::max();
922
923 // We can only own this DOF if it is either on the interior or border partition.
924 const bool is_interior_or_border = (e.partitionType()==Dune::InteriorEntity || e.partitionType()==Dune::BorderEntity);
925
926 // Receive data.
927 RankIndex received_rank;
928 buff.read(received_rank);
929
930 for (std::size_t i = 0; i < local_view.size(); ++i)
931 {
932 // Get the currently stored owner rank for this DOF.
933 RankIndex current_rank = local_view[i];
934
935 // We only gather from interior and border entities, so we need to make sure
936 // we relinquish any ownership claims on overlap and ghost entities on the
937 // receiving side. We also need to make sure not to overwrite any data already
938 // received, so we only blank the rank value if the currently stored value is
939 // equal to our own rank.
940 if (!is_interior_or_border && current_rank == _rank)
941 current_rank = unknown_rank;
942
943 // Assign DOFs to minimum rank value.
944 local_view[i] = std::min(current_rank,received_rank);
945 }
946 return true;
947 }
948
950
954 : _rank(rank)
955 {}
956
957 private:
958
959 const RankIndex _rank;
960
961 };
962
964
972 template<class GFS, class V>
974 : public Dune::PDELab::GFSDataHandle<GFS,
975 V,
976 DisjointPartitioningGatherScatter<
977 typename V::ElementType
978 >,
979 EntityDataCommunicationDescriptor<
980 typename V::ElementType
981 >
982 >
983 {
985 GFS,
986 V,
988 typename V::ElementType
989 >,
991 typename V::ElementType
992 >
993 > BaseT;
994
995 public:
996
998
1007 DisjointPartitioningDataHandle(const GFS& gfs_, V& v_, bool init_vector = true)
1008 : BaseT(gfs_,v_,DisjointPartitioningGatherScatter<typename V::ElementType>(gfs_.gridView().comm().rank()))
1009 {
1010 if (init_vector)
1011 v_ = gfs_.gridView().comm().rank();
1012 }
1013 };
1014
1015
1017
1024 {
1025
1026 template<typename MessageBuffer, typename Entity, typename LocalView>
1027 bool gather(MessageBuffer& buff, const Entity& e, LocalView& local_view) const
1028 {
1029 buff.write(local_view.size() > 0);
1030 return false;
1031 }
1032
1033 template<typename MessageBuffer, typename Entity, typename LocalView>
1034 bool scatter(MessageBuffer& buff, std::size_t n, const Entity& e, LocalView& local_view) const
1035 {
1036 bool remote_entity_has_dofs;
1037 buff.read(remote_entity_has_dofs);
1038
1039 for (std::size_t i = 0; i < local_view.size(); ++i)
1040 {
1041 local_view[i] |= remote_entity_has_dofs;
1042 }
1043 return true;
1044 }
1045
1046 };
1047
1048
1050
1056 template<class GFS, class V>
1058 : public Dune::PDELab::GFSDataHandle<GFS,
1059 V,
1060 SharedDOFGatherScatter,
1061 EntityDataCommunicationDescriptor<bool> >
1062 {
1064 GFS,
1065 V,
1068 > BaseT;
1069
1070 static_assert((std::is_same<typename V::ElementType,bool>::value),
1071 "SharedDOFDataHandle expects a vector of bool values");
1072
1073 public:
1074
1076
1085 SharedDOFDataHandle(const GFS& gfs_, V& v_, bool init_vector = true)
1086 : BaseT(gfs_,v_)
1087 {
1088 if (init_vector)
1089 v_ = false;
1090 }
1091 };
1092
1093
1095
1102 template<typename GFS, typename RankIndex>
1104 : public Dune::CommDataHandleIF<GFSNeighborDataHandle<GFS,RankIndex>,RankIndex>
1105 {
1106
1107 // We deliberately avoid using the GFSDataHandle here, as we don't want to incur the
1108 // overhead of invoking the whole GFS infrastructure.
1109
1110 public:
1111
1112 typedef RankIndex DataType;
1113 typedef typename GFS::Traits::SizeType size_type;
1114
1115 GFSNeighborDataHandle(const GFS& gfs, RankIndex rank, std::set<RankIndex>& neighbors)
1116 : _gfs(gfs)
1117 , _rank(rank)
1118 , _neighbors(neighbors)
1119 {}
1120
1121 bool contains(int dim, int codim) const
1122 {
1123 // Only create neighbor relations for codims used by the GFS.
1124 return _gfs.dataHandleContains(codim);
1125 }
1126
1127 bool fixedSize(int dim, int codim) const
1128 {
1129 // We always send a single value, the MPI rank.
1130 return true;
1131 }
1132
1133 template<typename Entity>
1134 size_type size(Entity& e) const
1135 {
1136 return 1;
1137 }
1138
1139 template<typename MessageBuffer, typename Entity>
1140 void gather(MessageBuffer& buff, const Entity& e) const
1141 {
1142 buff.write(_rank);
1143 }
1144
1145 template<typename MessageBuffer, typename Entity>
1146 void scatter(MessageBuffer& buff, const Entity& e, size_type n)
1147 {
1148 RankIndex rank;
1149 buff.read(rank);
1150 _neighbors.insert(rank);
1151 }
1152
1153 private:
1154
1155 const GFS& _gfs;
1156 const RankIndex _rank;
1157 std::set<RankIndex>& _neighbors;
1158
1159 };
1160
1161
1162 } // namespace PDELab
1163} // namespace Dune
1164
1165#endif // DUNE_PDELAB_GRIDFUNCTIONSPACE_GENERICDATAHANDLE_HH
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
@ codimension
Know your own codimension.
Definition: entity.hh:105
GatherScatter data handle for creating a disjoint DOF partitioning.
Definition: genericdatahandle.hh:983
DisjointPartitioningDataHandle(const GFS &gfs_, V &v_, bool init_vector=true)
Creates a new DisjointPartitioningDataHandle.
Definition: genericdatahandle.hh:1007
GatherScatter functor for creating a disjoint DOF partitioning.
Definition: genericdatahandle.hh:903
DisjointPartitioningGatherScatter(RankIndex rank)
Create a DisjointPartitioningGatherScatter object.
Definition: genericdatahandle.hh:953
Implement a data handle with a grid function space.
Definition: genericdatahandle.hh:132
bool contains(int dim, int codim) const
returns true if data for this codim should be communicated
Definition: genericdatahandle.hh:151
std::enable_if< CommunicationDescriptor::wrap_buffer &&AlwaysTrue< Entity >::value >::type gather(MessageBuffer &buff, const Entity &e) const
pack data from user to message buffer - version with support for sending leaf ordering sizes
Definition: genericdatahandle.hh:177
std::enable_if<!CommunicationDescriptor::wrap_buffer &&AlwaysTrue< Entity >::value >::type gather(MessageBuffer &buff, const Entity &e) const
pack data from user to message buffer - version without support for sending leaf ordering sizes
Definition: genericdatahandle.hh:209
size_type size(const Entity &e) const
how many objects of type DataType have to be sent for a given entity
Definition: genericdatahandle.hh:167
bool fixedSize(int dim, int codim) const
returns true if size per entity of given dim and codim is a constant
Definition: genericdatahandle.hh:157
std::enable_if<!CommunicationDescriptor::wrap_buffer &&AlwaysTrue< Entity >::value >::type scatter(MessageBuffer &buff, const Entity &e, size_type n)
unpack data from message buffer to user
Definition: genericdatahandle.hh:284
std::enable_if< CommunicationDescriptor::wrap_buffer &&AlwaysTrue< Entity >::value >::type scatter(MessageBuffer &buff, const Entity &e, size_type n)
unpack data from message buffer to user
Definition: genericdatahandle.hh:228
Data handle for collecting set of neighboring MPI ranks.
Definition: genericdatahandle.hh:1105
Data handle for marking ghost DOFs.
Definition: genericdatahandle.hh:860
GhostDataHandle(const GFS &gfs_, V &v_, bool init_vector=true)
Creates a new GhostDataHandle.
Definition: genericdatahandle.hh:883
GatherScatter functor for marking ghost DOFs.
Definition: genericdatahandle.hh:812
Wrapper for message buffers of grid DataHandles that allows for sending different types of data.
Definition: polymorphicbufferwrapper.hh:32
Data handle for marking shared DOFs.
Definition: genericdatahandle.hh:1062
SharedDOFDataHandle(const GFS &gfs_, V &v_, bool init_vector=true)
Creates a new SharedDOFDataHandle.
Definition: genericdatahandle.hh:1085
Describes the parallel communication interface class for MessageBuffers and DataHandles.
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
@ InteriorEntity
all interior entities
Definition: gridenums.hh:29
@ BorderEntity
on boundary between interior and overlap
Definition: gridenums.hh:30
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:87
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
Helpers for dealing with MPI.
Dune namespace.
Definition: alignedallocator.hh:14
template which always yields a true value
Definition: typetraits.hh:139
Communication descriptor for sending one item of type E per DOF.
Definition: genericdatahandle.hh:25
std::size_t size_type
size type to use if communicating leaf ordering sizes
Definition: genericdatahandle.hh:30
Communication descriptor for sending count items of type E per entity with attached DOFs.
Definition: genericdatahandle.hh:70
std::size_t size_type
size type to use if communicating leaf ordering sizes
Definition: genericdatahandle.hh:73
GatherScatter functor for marking shared DOFs.
Definition: genericdatahandle.hh:1024
Struct for obtaining some basic structural information about a TypeTree.
Definition: utility.hh:66
static const std::size_t leafCount
The number of leaf nodes in the TypeTree.
Definition: utility.hh:81
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)