16 #ifndef dealii_quadrature_point_data_h
17 #define dealii_quadrature_point_data_h
36 #include <type_traits>
63 template <
typename CellIteratorType,
typename DataType>
102 template <
typename T = DataType>
105 const unsigned int number_of_data_points_per_cell);
112 template <
typename T = DataType>
114 initialize(
const CellIteratorType &cell_start,
115 const CellIteratorType &cell_end,
116 const unsigned int number_of_data_points_per_cell);
128 erase(
const CellIteratorType &cell);
150 template <
typename T = DataType>
151 std::vector<std::shared_ptr<T>>
152 get_data(
const CellIteratorType &cell);
168 template <
typename T = DataType>
169 std::vector<std::shared_ptr<const T>>
170 get_data(
const CellIteratorType &cell)
const;
188 template <
typename T = DataType>
189 std_cxx17::optional<std::vector<std::shared_ptr<T>>>
208 template <
typename T = DataType>
209 std_cxx17::optional<std::vector<std::shared_ptr<const T>>>
217 CellIteratorType::AccessorType::dimension;
223 CellIteratorType::AccessorType::space_dimension;
239 std::map<CellId, std::vector<std::shared_ptr<DataType>>>
map;
246 "Cell data is being retrieved with a type which is different than the type used to initialize it");
253 "The provided cell iterator does not belong to the triangulation that corresponds to the CellDataStorage object.");
311 pack_values(std::vector<double> &values)
const = 0;
326 #ifdef DEAL_II_WITH_P4EST
329 namespace distributed
442 template <
int dim,
typename DataType>
448 "User's DataType class should be derived from TransferableQuadraturePointData");
522 const boost::iterator_range<std::vector<char>::const_iterator>
606 template <
typename CellIteratorType,
typename DataType>
607 template <
typename T>
610 const CellIteratorType &cell,
611 const unsigned int n_q_points)
614 "User's T class should be derived from user's DataType class");
618 tria = &cell->get_triangulation();
621 const auto key = cell->id();
622 if (map.find(key) == map.end())
624 map[key] = std::vector<std::shared_ptr<DataType>>(n_q_points);
628 const auto it = map.find(key);
629 for (
unsigned int q = 0; q < n_q_points; q++)
630 it->second[q] = std::make_shared<T>();
636 template <
typename CellIteratorType,
typename DataType>
637 template <
typename T>
640 const CellIteratorType &cell_start,
641 const CellIteratorType &cell_end,
642 const unsigned int number)
644 for (CellIteratorType it = cell_start; it != cell_end; it++)
645 if (it->is_locally_owned())
646 initialize<T>(it, number);
651 template <
typename CellIteratorType,
typename DataType>
655 const auto key = cell->id();
656 const auto it = map.find(key);
660 for (
unsigned int i = 0; i < it->second.size(); i++)
663 it->second[i].unique(),
665 "Can not erase the cell data multiple objects reference its data."));
668 return (map.erase(key) == 1);
673 template <
typename CellIteratorType,
typename DataType>
680 auto it = map.begin();
681 while (it != map.end())
684 for (
unsigned int i = 0; i < it->second.size(); i++)
687 it->second[i].unique(),
689 "Can not erase the cell data, multiple objects reference it."));
697 template <
typename CellIteratorType,
typename DataType>
698 template <
typename T>
699 inline std::vector<std::shared_ptr<T>>
701 const CellIteratorType &cell)
704 "User's T class should be derived from user's DataType class");
707 const auto it = map.find(cell->id());
715 std::vector<std::shared_ptr<T>> res(it->second.size());
716 for (
unsigned int q = 0; q < res.size(); q++)
718 res[q] = std::dynamic_pointer_cast<T>(it->second[q]);
719 Assert(res[q], ExcCellDataTypeMismatch());
726 template <
typename CellIteratorType,
typename DataType>
727 template <
typename T>
728 inline std::vector<std::shared_ptr<const T>>
730 const CellIteratorType &cell)
const
733 "User's T class should be derived from user's DataType class");
736 const auto it = map.find(cell->id());
742 std::vector<std::shared_ptr<const T>> res(it->second.size());
743 for (
unsigned int q = 0; q < res.size(); q++)
745 res[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
746 Assert(res[q], ExcCellDataTypeMismatch());
751 template <
typename CellIteratorType,
typename DataType>
752 template <
typename T>
753 inline std_cxx17::optional<std::vector<std::shared_ptr<T>>>
755 const CellIteratorType &cell)
758 "User's T class should be derived from user's DataType class");
761 const auto it = map.find(cell->id());
768 std::vector<std::shared_ptr<T>> result(it->second.size());
769 for (
unsigned int q = 0; q < result.size(); q++)
771 result[q] = std::dynamic_pointer_cast<T>(it->second[q]);
772 Assert(result[q], ExcCellDataTypeMismatch());
782 template <
typename CellIteratorType,
typename DataType>
783 template <
typename T>
784 inline std_cxx17::optional<std::vector<std::shared_ptr<const T>>>
786 const CellIteratorType &cell)
const
789 "User's T class should be derived from user's DataType class");
792 const auto it = map.find(cell->id());
799 std::vector<std::shared_ptr<const T>> result(it->second.size());
800 for (
unsigned int q = 0; q < result.size(); q++)
802 result[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
803 Assert(result[q], ExcCellDataTypeMismatch());
825 template <
typename CellIteratorType,
typename DataType>
827 pack_cell_data(
const CellIteratorType & cell,
833 "User's DataType class should be derived from QPData");
837 const unsigned int m = qpd->size();
839 const unsigned int n = (*qpd)[0]->number_of_values();
842 std::vector<double> single_qp_data(n);
843 for (
unsigned int q = 0; q < m; ++q)
845 (*qpd)[q]->pack_values(single_qp_data);
848 for (
unsigned int i = 0; i < n; ++i)
849 matrix_data(q, i) = single_qp_data[i];
854 matrix_data.
reinit({0, 0});
863 template <
typename CellIteratorType,
typename DataType>
865 unpack_to_cell_data(
const CellIteratorType & cell,
871 "User's DataType class should be derived from QPData");
875 const unsigned int n = values_at_qp.
n();
878 std::vector<double> single_qp_data(n);
881 for (
unsigned int q = 0; q < qpd->size(); ++q)
883 for (
unsigned int i = 0; i < n; ++i)
884 single_qp_data[i] = values_at_qp(q, i);
885 (*qpd)[q]->unpack_values(single_qp_data);
891 # ifdef DEAL_II_WITH_P4EST
895 namespace distributed
897 template <
int dim,
typename DataType>
903 std::unique_ptr<const
FiniteElement<dim>>(projection_fe_.clone()))
904 , data_size_in_bytes(0)
905 , n_q_points(rhs_quadrature.size())
906 , project_to_fe_matrix(projection_fe->dofs_per_cell, n_q_points)
907 , project_to_qp_matrix(n_q_points, projection_fe->dofs_per_cell)
909 , data_storage(nullptr)
915 "ContinuousQuadratureDataTransfer requires scalar FiniteElement"));
929 template <
int dim,
typename DataType>
936 Assert(data_storage ==
nullptr,
937 ExcMessage(
"This function can be called only once"));
939 data_storage = &data_storage_;
944 dim>::cell_iterator &cell,
946 status) {
return this->pack_function(cell, status); },
952 template <
int dim,
typename DataType>
960 dim>::cell_iterator &cell,
963 const boost::iterator_range<std::vector<char>::const_iterator>
964 &data_range) { this->unpack_function(cell, status, data_range); });
967 data_storage =
nullptr;
973 template <
int dim,
typename DataType>
974 inline std::vector<char>
981 pack_cell_data(cell, data_storage, matrix_quadrature);
984 const unsigned int number_of_values = matrix_quadrature.n();
985 matrix_dofs.reinit(project_to_fe_matrix.m(), number_of_values);
986 if (number_of_values > 0)
987 project_to_fe_matrix.mmult(matrix_dofs, matrix_quadrature);
994 template <
int dim,
typename DataType>
1001 const boost::iterator_range<std::vector<char>::const_iterator>
1010 Utilities::unpack<FullMatrix<double>>(data_range.begin(),
1013 const unsigned int number_of_values = matrix_dofs.n();
1014 if (number_of_values == 0)
1017 matrix_quadrature.reinit(n_q_points, number_of_values);
1019 if (cell->has_children())
1023 matrix_dofs_child.reinit(projection_fe->dofs_per_cell,
1025 for (
unsigned int child = 0; child < cell->n_children(); ++child)
1026 if (cell->child(child)->is_locally_owned())
1029 ->get_prolongation_matrix(child, cell->refinement_case())
1030 .mmult(matrix_dofs_child, matrix_dofs);
1034 project_to_qp_matrix.mmult(matrix_quadrature,
1038 unpack_to_cell_data(cell->child(child),
1047 project_to_qp_matrix.mmult(matrix_quadrature, matrix_dofs);
1050 unpack_to_cell_data(cell, matrix_quadrature, data_storage);
1058 # endif // DEAL_II_WITH_P4EST