Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
quadrature_point_data.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2016 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_quadrature_point_data_h
17#define dealii_quadrature_point_data_h
18
19#include <deal.II/base/config.h>
20
24
26
27#include <deal.II/fe/fe.h>
28#include <deal.II/fe/fe_tools.h>
29
30#include <deal.II/grid/tria.h>
32
33#include <deal.II/lac/vector.h>
34
35#include <map>
36#include <type_traits>
37#include <vector>
38
40
43
61template <typename CellIteratorType, typename DataType>
63{
64public:
68 CellDataStorage() = default;
69
73 ~CellDataStorage() override = default;
74
100 template <typename T = DataType>
101 void
102 initialize(const CellIteratorType &cell,
103 const unsigned int number_of_data_points_per_cell);
104
110 template <typename T = DataType>
111 void
112 initialize(const CellIteratorType &cell_start,
113 const CellIteratorType &cell_end,
114 const unsigned int number_of_data_points_per_cell);
115
125 bool
126 erase(const CellIteratorType &cell);
127
131 void
133
148 template <typename T = DataType>
149 std::vector<std::shared_ptr<T>>
150 get_data(const CellIteratorType &cell);
151
166 template <typename T = DataType>
167 std::vector<std::shared_ptr<const T>>
168 get_data(const CellIteratorType &cell) const;
169
186 template <typename T = DataType>
187 std_cxx17::optional<std::vector<std::shared_ptr<T>>>
188 try_get_data(const CellIteratorType &cell);
189
206 template <typename T = DataType>
207 std_cxx17::optional<std::vector<std::shared_ptr<const T>>>
208 try_get_data(const CellIteratorType &cell) const;
209
210private:
214 static constexpr unsigned int dimension =
215 CellIteratorType::AccessorType::dimension;
216
220 static constexpr unsigned int space_dimension =
221 CellIteratorType::AccessorType::space_dimension;
222
231
237 std::map<CellId, std::vector<std::shared_ptr<DataType>>> map;
238
244 "Cell data is being retrieved with a type which is different than the type used to initialize it");
245
251 "The provided cell iterator does not belong to the triangulation that corresponds to the CellDataStorage object.");
252};
253
254
276{
277public:
282
287
293 virtual unsigned int
294 number_of_values() const = 0;
295
306 virtual void
307 pack_values(std::vector<double> &values) const = 0;
308
317 virtual void
318 unpack_values(const std::vector<double> &values) = 0;
319};
320
321
322#ifdef DEAL_II_WITH_P4EST
323namespace parallel
324{
325 namespace distributed
326 {
436 template <int dim, typename DataType>
438 {
439 public:
440 static_assert(
441 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
442 "User's DataType class should be derived from TransferableQuadraturePointData");
443
449
461 const Quadrature<dim> &mass_quadrature,
462 const Quadrature<dim> &data_quadrature);
463
474 void
478
489 void
491
492 private:
498 std::vector<char>
501 &cell,
503 status);
504
510 void
513 &cell,
515 status,
516 const boost::iterator_range<std::vector<char>::const_iterator>
517 &data_range);
518
522 const std::unique_ptr<const FiniteElement<dim>> projection_fe;
523
529
533 const unsigned int n_q_points;
534
540
546
554
559
565
570 unsigned int handle;
571
576
582 };
583
584 } // namespace distributed
585
586} // namespace parallel
587
588#endif
589
592#ifndef DOXYGEN
593
594// ------------------- inline and template functions ----------------
595
596//--------------------------------------------------------------------
597// CellDataStorage
598//--------------------------------------------------------------------
599
600template <typename CellIteratorType, typename DataType>
601template <typename T>
602inline void
604 const CellIteratorType &cell,
605 const unsigned int n_q_points)
606{
607 static_assert(std::is_base_of<DataType, T>::value,
608 "User's T class should be derived from user's DataType class");
609 // The first time this method is called, it has to initialize the reference
610 // to the triangulation object
611 if (!tria)
612 tria = &cell->get_triangulation();
613 Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
614
615 const auto key = cell->id();
616 if (map.find(key) == map.end())
617 {
618 map[key] = std::vector<std::shared_ptr<DataType>>(n_q_points);
619 // we need to initialize one-by-one as the std::vector<>(q, T())
620 // will end with a single same T object stored in each element of the
621 // vector:
622 const auto it = map.find(key);
623 for (unsigned int q = 0; q < n_q_points; ++q)
624 it->second[q] = std::make_shared<T>();
625 }
626}
627
628
629
630template <typename CellIteratorType, typename DataType>
631template <typename T>
632inline void
634 const CellIteratorType &cell_start,
635 const CellIteratorType &cell_end,
636 const unsigned int number)
637{
638 for (CellIteratorType it = cell_start; it != cell_end; ++it)
639 if (it->is_locally_owned())
640 initialize<T>(it, number);
641}
642
643
644
645template <typename CellIteratorType, typename DataType>
646inline bool
647CellDataStorage<CellIteratorType, DataType>::erase(const CellIteratorType &cell)
648{
649 const auto key = cell->id();
650 const auto it = map.find(key);
651 if (it == map.end())
652 return false;
653 Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
654 for (unsigned int i = 0; i < it->second.size(); ++i)
655 {
656 Assert(
657 it->second[i].unique(),
659 "Can not erase the cell data multiple objects reference its data."));
660 }
661
662 return (map.erase(key) == 1);
663}
664
665
666
667template <typename CellIteratorType, typename DataType>
668inline void
670{
671 // Do not call
672 // map.clear();
673 // as we want to be sure no one uses the stored objects. Loop manually:
674 auto it = map.begin();
675 while (it != map.end())
676 {
677 // loop over all objects and see if no one is using them
678 for (unsigned int i = 0; i < it->second.size(); ++i)
679 {
680 Assert(
681 it->second[i].unique(),
683 "Can not erase the cell data, multiple objects reference it."));
684 }
685 it = map.erase(it);
686 }
687}
688
689
690
691template <typename CellIteratorType, typename DataType>
692template <typename T>
693inline std::vector<std::shared_ptr<T>>
695 const CellIteratorType &cell)
696{
697 static_assert(std::is_base_of<DataType, T>::value,
698 "User's T class should be derived from user's DataType class");
699 Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
700
701 const auto it = map.find(cell->id());
702 Assert(it != map.end(), ExcMessage("Could not find data for the cell"));
703
704 // It would be nice to have a specialized version of this function for
705 // T==DataType. However explicit (i.e full) specialization of a member
706 // template is only allowed when the enclosing class is also explicitly (i.e
707 // fully) specialized. Thus, stick with copying of shared pointers even when
708 // the T==DataType:
709 std::vector<std::shared_ptr<T>> res(it->second.size());
710 for (unsigned int q = 0; q < res.size(); ++q)
711 {
712 res[q] = std::dynamic_pointer_cast<T>(it->second[q]);
713 Assert(res[q], ExcCellDataTypeMismatch());
714 }
715 return res;
716}
717
718
719
720template <typename CellIteratorType, typename DataType>
721template <typename T>
722inline std::vector<std::shared_ptr<const T>>
724 const CellIteratorType &cell) const
725{
726 static_assert(std::is_base_of<DataType, T>::value,
727 "User's T class should be derived from user's DataType class");
728 Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
729
730 const auto it = map.find(cell->id());
731 Assert(it != map.end(), ExcMessage("Could not find QP data for the cell"));
732
733 // Cast base class to the desired class. This has to be done irrespectively of
734 // T==DataType as we need to return shared_ptr<const T> to make sure the user
735 // does not modify the content of QP objects
736 std::vector<std::shared_ptr<const T>> res(it->second.size());
737 for (unsigned int q = 0; q < res.size(); ++q)
738 {
739 res[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
740 Assert(res[q], ExcCellDataTypeMismatch());
741 }
742 return res;
743}
744
745template <typename CellIteratorType, typename DataType>
746template <typename T>
747inline std_cxx17::optional<std::vector<std::shared_ptr<T>>>
749 const CellIteratorType &cell)
750{
751 static_assert(std::is_base_of<DataType, T>::value,
752 "User's T class should be derived from user's DataType class");
753 Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
754
755 const auto it = map.find(cell->id());
756 if (it != map.end())
757 {
758 // Cast base class to the desired class. This has to be done
759 // irrespectively of T==DataType as we need to return
760 // shared_ptr<const T> to make sure the user
761 // does not modify the content of QP objects
762 std::vector<std::shared_ptr<T>> result(it->second.size());
763 for (unsigned int q = 0; q < result.size(); ++q)
764 {
765 result[q] = std::dynamic_pointer_cast<T>(it->second[q]);
766 Assert(result[q], ExcCellDataTypeMismatch());
767 }
768 return {result};
769 }
770 else
771 {
772 return {};
773 }
774}
775
776template <typename CellIteratorType, typename DataType>
777template <typename T>
778inline std_cxx17::optional<std::vector<std::shared_ptr<const T>>>
780 const CellIteratorType &cell) const
781{
782 static_assert(std::is_base_of<DataType, T>::value,
783 "User's T class should be derived from user's DataType class");
784 Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
785
786 const auto it = map.find(cell->id());
787 if (it != map.end())
788 {
789 // Cast base class to the desired class. This has to be done
790 // irrespectively of T==DataType as we need to return
791 // shared_ptr<const T> to make sure the user
792 // does not modify the content of QP objects
793 std::vector<std::shared_ptr<const T>> result(it->second.size());
794 for (unsigned int q = 0; q < result.size(); ++q)
795 {
796 result[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
797 Assert(result[q], ExcCellDataTypeMismatch());
798 }
799 return {result};
800 }
801 else
802 {
803 return {};
804 }
805}
806
807//--------------------------------------------------------------------
808// ContinuousQuadratureDataTransfer
809//--------------------------------------------------------------------
810
811
812/*
813 * Pack cell data of type @p DataType stored using @p data_storage in @p cell
814 * at each quadrature point to @p matrix_data. Here @p matrix_data is a matrix
815 * whose first index corresponds to different quadrature points on the cell
816 * whereas the second index represents different values stored at each
817 * quadrature point in the DataType class.
818 */
819template <typename CellIteratorType, typename DataType>
820inline void
821pack_cell_data(const CellIteratorType & cell,
823 FullMatrix<double> & matrix_data)
824{
825 static_assert(
826 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
827 "User's DataType class should be derived from QPData");
828
829 if (const auto qpd = data_storage->try_get_data(cell))
830 {
831 const unsigned int m = qpd->size();
832 Assert(m > 0, ExcInternalError());
833 const unsigned int n = (*qpd)[0]->number_of_values();
834 matrix_data.reinit(m, n);
835
836 std::vector<double> single_qp_data(n);
837 for (unsigned int q = 0; q < m; ++q)
838 {
839 (*qpd)[q]->pack_values(single_qp_data);
840 AssertDimension(single_qp_data.size(), n);
841
842 for (unsigned int i = 0; i < n; ++i)
843 matrix_data(q, i) = single_qp_data[i];
844 }
845 }
846 else
847 {
848 matrix_data.reinit({0, 0});
849 }
850}
851
852
853
854/*
855 * the opposite of the pack function above.
856 */
857template <typename CellIteratorType, typename DataType>
858inline void
859unpack_to_cell_data(const CellIteratorType & cell,
860 const FullMatrix<double> & values_at_qp,
862{
863 static_assert(
864 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
865 "User's DataType class should be derived from QPData");
866
867 if (const auto qpd = data_storage->try_get_data(cell))
868 {
869 const unsigned int n = values_at_qp.n();
870 AssertDimension((*qpd)[0]->number_of_values(), n);
871
872 std::vector<double> single_qp_data(n);
873 AssertDimension(qpd->size(), values_at_qp.m());
874
875 for (unsigned int q = 0; q < qpd->size(); ++q)
876 {
877 for (unsigned int i = 0; i < n; ++i)
878 single_qp_data[i] = values_at_qp(q, i);
879 (*qpd)[q]->unpack_values(single_qp_data);
880 }
881 }
882}
883
884
885# ifdef DEAL_II_WITH_P4EST
886
887namespace parallel
888{
889 namespace distributed
890 {
891 template <int dim, typename DataType>
894 const Quadrature<dim> & lhs_quadrature,
895 const Quadrature<dim> & rhs_quadrature)
896 : projection_fe(
897 std::unique_ptr<const FiniteElement<dim>>(projection_fe_.clone()))
898 , data_size_in_bytes(0)
899 , n_q_points(rhs_quadrature.size())
900 , project_to_fe_matrix(projection_fe->n_dofs_per_cell(), n_q_points)
901 , project_to_qp_matrix(n_q_points, projection_fe->n_dofs_per_cell())
902 , handle(numbers::invalid_unsigned_int)
903 , data_storage(nullptr)
904 , triangulation(nullptr)
905 {
906 Assert(
907 projection_fe->n_components() == 1,
909 "ContinuousQuadratureDataTransfer requires scalar FiniteElement"));
910
912 *projection_fe.get(),
913 lhs_quadrature,
914 rhs_quadrature,
916
918 *projection_fe.get(), rhs_quadrature, project_to_qp_matrix);
919 }
920
921
922
923 template <int dim, typename DataType>
924 inline void
929 {
930 Assert(data_storage == nullptr,
931 ExcMessage("This function can be called only once"));
932 triangulation = &tr_;
933 data_storage = &data_storage_;
934
935 handle = triangulation->register_data_attach(
936 [this](
938 dim>::cell_iterator &cell,
940 status) { return this->pack_function(cell, status); },
941 /*returns_variable_size_data=*/true);
942 }
943
944
945
946 template <int dim, typename DataType>
947 inline void
949 {
950 triangulation->notify_ready_to_unpack(
951 handle,
952 [this](
954 dim>::cell_iterator &cell,
956 status,
957 const boost::iterator_range<std::vector<char>::const_iterator>
958 &data_range) { this->unpack_function(cell, status, data_range); });
959
960 // invalidate the pointers
961 data_storage = nullptr;
962 triangulation = nullptr;
963 }
964
965
966
967 template <int dim, typename DataType>
968 inline std::vector<char>
971 &cell,
973 dim>::CellStatus /*status*/)
974 {
975 pack_cell_data(cell, data_storage, matrix_quadrature);
976
977 // project to FE
978 const unsigned int number_of_values = matrix_quadrature.n();
979 matrix_dofs.reinit(project_to_fe_matrix.m(), number_of_values);
980 if (number_of_values > 0)
981 project_to_fe_matrix.mmult(matrix_dofs, matrix_quadrature);
982
983 return Utilities::pack(matrix_dofs, /*allow_compression=*/false);
984 }
985
986
987
988 template <int dim, typename DataType>
989 inline void
992 &cell,
994 status,
995 const boost::iterator_range<std::vector<char>::const_iterator>
996 &data_range)
997 {
998 Assert((status !=
1001 (void)status;
1002
1003 matrix_dofs =
1004 Utilities::unpack<FullMatrix<double>>(data_range.begin(),
1005 data_range.end(),
1006 /*allow_compression=*/false);
1007 const unsigned int number_of_values = matrix_dofs.n();
1008 if (number_of_values == 0)
1009 return;
1010
1011 matrix_quadrature.reinit(n_q_points, number_of_values);
1012
1013 if (cell->has_children())
1014 {
1015 // we need to first use prolongation matrix to get dofvalues on child
1016 // cells based on dofvalues stored in the parent's data_store
1017 matrix_dofs_child.reinit(projection_fe->n_dofs_per_cell(),
1018 number_of_values);
1019 for (unsigned int child = 0; child < cell->n_children(); ++child)
1020 if (cell->child(child)->is_locally_owned())
1021 {
1022 projection_fe
1023 ->get_prolongation_matrix(child, cell->refinement_case())
1024 .mmult(matrix_dofs_child, matrix_dofs);
1025
1026 // now we do the usual business of evaluating FE on quadrature
1027 // points:
1028 project_to_qp_matrix.mmult(matrix_quadrature,
1029 matrix_dofs_child);
1030
1031 // finally, put back into the map:
1032 unpack_to_cell_data(cell->child(child),
1033 matrix_quadrature,
1034 data_storage);
1035 }
1036 }
1037 else
1038 {
1039 // if there are no children, evaluate FE field at
1040 // rhs_quadrature points.
1041 project_to_qp_matrix.mmult(matrix_quadrature, matrix_dofs);
1042
1043 // finally, put back into the map:
1044 unpack_to_cell_data(cell, matrix_quadrature, data_storage);
1045 }
1046 }
1047
1048 } // namespace distributed
1049
1050} // namespace parallel
1051
1052# endif // DEAL_II_WITH_P4EST
1053
1054#endif // DOXYGEN
1056
1057#endif
std_cxx17::optional< std::vector< std::shared_ptr< const T > > > try_get_data(const CellIteratorType &cell) const
CellDataStorage()=default
void initialize(const CellIteratorType &cell, const unsigned int number_of_data_points_per_cell)
bool erase(const CellIteratorType &cell)
static constexpr unsigned int space_dimension
std_cxx17::optional< std::vector< std::shared_ptr< T > > > try_get_data(const CellIteratorType &cell)
SmartPointer< const Triangulation< dimension, space_dimension >, CellDataStorage< CellIteratorType, DataType > > tria
std::vector< std::shared_ptr< const T > > get_data(const CellIteratorType &cell) const
void initialize(const CellIteratorType &cell_start, const CellIteratorType &cell_end, const unsigned int number_of_data_points_per_cell)
~CellDataStorage() override=default
std::map< CellId, std::vector< std::shared_ptr< DataType > > > map
std::vector< std::shared_ptr< T > > get_data(const CellIteratorType &cell)
static constexpr unsigned int dimension
size_type n() const
size_type m() const
virtual ~TransferableQuadraturePointData()=default
virtual unsigned int number_of_values() const =0
virtual void unpack_values(const std::vector< double > &values)=0
virtual void pack_values(std::vector< double > &values) const =0
Triangulation< dim, spacedim > & get_triangulation()
const std::unique_ptr< const FiniteElement< dim > > projection_fe
ContinuousQuadratureDataTransfer(const FiniteElement< dim > &projection_fe, const Quadrature< dim > &mass_quadrature, const Quadrature< dim > &data_quadrature)
parallel::distributed::Triangulation< dim > * triangulation
typename parallel::distributed::Triangulation< dim >::cell_iterator CellIteratorType
void unpack_function(const typename parallel::distributed::Triangulation< dim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim >::CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range)
void prepare_for_coarsening_and_refinement(parallel::distributed::Triangulation< dim > &tria, CellDataStorage< CellIteratorType, DataType > &data_storage)
std::vector< char > pack_function(const typename parallel::distributed::Triangulation< dim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim >::CellStatus status)
CellDataStorage< CellIteratorType, DataType > * data_storage
typename ::Triangulation< dim, spacedim >::CellStatus CellStatus
Definition: tria.h:294
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & ExcCellDataTypeMismatch()
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcTriangulationMismatch()
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
Definition: tria.h:270
void compute_interpolation_to_quadrature_points_matrix(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, FullMatrix< double > &I_q)
void compute_projection_from_quadrature_points_matrix(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &lhs_quadrature, const Quadrature< dim > &rhs_quadrature, FullMatrix< double > &X)
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
Definition: utilities.h:1483
STL namespace.
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
const ::Triangulation< dim, spacedim > & tria