16#ifndef dealii_quadrature_point_data_h
17#define dealii_quadrature_point_data_h
61template <
typename CellIteratorType,
typename DataType>
100 template <
typename T = DataType>
103 const unsigned int number_of_data_points_per_cell);
110 template <
typename T = DataType>
113 const CellIteratorType &cell_end,
114 const unsigned int number_of_data_points_per_cell);
126 erase(
const CellIteratorType &cell);
148 template <
typename T = DataType>
149 std::vector<std::shared_ptr<T>>
166 template <
typename T = DataType>
167 std::vector<std::shared_ptr<const T>>
186 template <
typename T = DataType>
187 std_cxx17::optional<std::vector<std::shared_ptr<T>>>
206 template <
typename T = DataType>
207 std_cxx17::optional<std::vector<std::shared_ptr<const T>>>
215 CellIteratorType::AccessorType::dimension;
221 CellIteratorType::AccessorType::space_dimension;
237 std::map<CellId, std::vector<std::shared_ptr<DataType>>>
map;
244 "Cell data is being retrieved with a type which is different than the type used to initialize it");
251 "The provided cell iterator does not belong to the triangulation that corresponds to the CellDataStorage object.");
322#ifdef DEAL_II_WITH_P4EST
325 namespace distributed
436 template <
int dim,
typename DataType>
441 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
442 "User's DataType class should be derived from TransferableQuadraturePointData");
516 const boost::iterator_range<std::vector<char>::const_iterator>
600template <
typename CellIteratorType,
typename DataType>
604 const CellIteratorType &cell,
605 const unsigned int n_q_points)
607 static_assert(std::is_base_of<DataType, T>::value,
608 "User's T class should be derived from user's DataType class");
613 Assert(&cell->get_triangulation() ==
tria, ExcTriangulationMismatch());
615 const auto key = cell->id();
616 if (map.find(key) == map.end())
618 map[key] = std::vector<std::shared_ptr<DataType>>(n_q_points);
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>();
630template <
typename CellIteratorType,
typename DataType>
634 const CellIteratorType &cell_start,
635 const CellIteratorType &cell_end,
636 const unsigned int number)
638 for (CellIteratorType it = cell_start; it != cell_end; ++it)
639 if (it->is_locally_owned())
640 initialize<T>(it, number);
645template <
typename CellIteratorType,
typename DataType>
649 const auto key = cell->id();
650 const auto it = map.find(key);
653 Assert(&cell->get_triangulation() ==
tria, ExcTriangulationMismatch());
654 for (
unsigned int i = 0; i < it->second.size(); ++i)
657 it->second[i].unique(),
659 "Can not erase the cell data multiple objects reference its data."));
662 return (map.erase(key) == 1);
667template <
typename CellIteratorType,
typename DataType>
674 auto it = map.begin();
675 while (it != map.end())
678 for (
unsigned int i = 0; i < it->second.size(); ++i)
681 it->second[i].unique(),
683 "Can not erase the cell data, multiple objects reference it."));
691template <
typename CellIteratorType,
typename DataType>
693inline std::vector<std::shared_ptr<T>>
695 const CellIteratorType &cell)
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());
701 const auto it = map.find(cell->id());
709 std::vector<std::shared_ptr<T>> res(it->second.size());
710 for (
unsigned int q = 0; q < res.size(); ++q)
712 res[q] = std::dynamic_pointer_cast<T>(it->second[q]);
713 Assert(res[q], ExcCellDataTypeMismatch());
720template <
typename CellIteratorType,
typename DataType>
722inline std::vector<std::shared_ptr<const T>>
724 const CellIteratorType &cell)
const
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());
730 const auto it = map.find(cell->id());
736 std::vector<std::shared_ptr<const T>> res(it->second.size());
737 for (
unsigned int q = 0; q < res.size(); ++q)
739 res[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
740 Assert(res[q], ExcCellDataTypeMismatch());
745template <
typename CellIteratorType,
typename DataType>
747inline std_cxx17::optional<std::vector<std::shared_ptr<T>>>
749 const CellIteratorType &cell)
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());
755 const auto it = map.find(cell->id());
762 std::vector<std::shared_ptr<T>> result(it->second.size());
763 for (
unsigned int q = 0; q < result.size(); ++q)
765 result[q] = std::dynamic_pointer_cast<T>(it->second[q]);
766 Assert(result[q], ExcCellDataTypeMismatch());
776template <
typename CellIteratorType,
typename DataType>
778inline std_cxx17::optional<std::vector<std::shared_ptr<const T>>>
780 const CellIteratorType &cell)
const
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());
786 const auto it = map.find(cell->id());
793 std::vector<std::shared_ptr<const T>> result(it->second.size());
794 for (
unsigned int q = 0; q < result.size(); ++q)
796 result[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
797 Assert(result[q], ExcCellDataTypeMismatch());
819template <
typename CellIteratorType,
typename DataType>
821pack_cell_data(
const CellIteratorType & cell,
826 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
827 "User's DataType class should be derived from QPData");
831 const unsigned int m = qpd->size();
833 const unsigned int n = (*qpd)[0]->number_of_values();
834 matrix_data.reinit(m, n);
836 std::vector<double> single_qp_data(n);
837 for (
unsigned int q = 0; q < m; ++q)
839 (*qpd)[q]->pack_values(single_qp_data);
842 for (
unsigned int i = 0; i < n; ++i)
843 matrix_data(q, i) = single_qp_data[i];
848 matrix_data.reinit({0, 0});
857template <
typename CellIteratorType,
typename DataType>
859unpack_to_cell_data(
const CellIteratorType & cell,
864 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
865 "User's DataType class should be derived from QPData");
869 const unsigned int n = values_at_qp.
n();
872 std::vector<double> single_qp_data(n);
875 for (
unsigned int q = 0; q < qpd->size(); ++q)
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);
885# ifdef DEAL_II_WITH_P4EST
889 namespace distributed
891 template <
int dim,
typename DataType>
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)
909 "ContinuousQuadratureDataTransfer requires scalar FiniteElement"));
923 template <
int dim,
typename DataType>
930 Assert(data_storage ==
nullptr,
931 ExcMessage(
"This function can be called only once"));
933 data_storage = &data_storage_;
938 dim>::cell_iterator &cell,
940 status) {
return this->pack_function(cell, status); },
946 template <
int dim,
typename DataType>
954 dim>::cell_iterator &cell,
957 const boost::iterator_range<std::vector<char>::const_iterator>
958 &data_range) { this->unpack_function(cell, status, data_range); });
961 data_storage =
nullptr;
967 template <
int dim,
typename DataType>
968 inline std::vector<char>
975 pack_cell_data(cell, data_storage, matrix_quadrature);
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);
988 template <
int dim,
typename DataType>
995 const boost::iterator_range<std::vector<char>::const_iterator>
1004 Utilities::unpack<FullMatrix<double>>(data_range.begin(),
1007 const unsigned int number_of_values = matrix_dofs.n();
1008 if (number_of_values == 0)
1011 matrix_quadrature.reinit(n_q_points, number_of_values);
1013 if (cell->has_children())
1017 matrix_dofs_child.reinit(projection_fe->n_dofs_per_cell(),
1019 for (
unsigned int child = 0; child < cell->n_children(); ++child)
1020 if (cell->child(child)->is_locally_owned())
1023 ->get_prolongation_matrix(child, cell->refinement_case())
1024 .mmult(matrix_dofs_child, matrix_dofs);
1028 project_to_qp_matrix.mmult(matrix_quadrature,
1032 unpack_to_cell_data(cell->child(child),
1041 project_to_qp_matrix.mmult(matrix_quadrature, matrix_dofs);
1044 unpack_to_cell_data(cell, matrix_quadrature, data_storage);
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
virtual ~TransferableQuadraturePointData()=default
virtual unsigned int number_of_values() const =0
virtual void unpack_values(const std::vector< double > &values)=0
TransferableQuadraturePointData()=default
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)
std::size_t data_size_in_bytes
parallel::distributed::Triangulation< dim > * triangulation
typename parallel::distributed::Triangulation< dim >::cell_iterator CellIteratorType
FullMatrix< double > matrix_quadrature
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)
FullMatrix< double > matrix_dofs
FullMatrix< double > project_to_qp_matrix
void prepare_for_coarsening_and_refinement(parallel::distributed::Triangulation< dim > &tria, CellDataStorage< CellIteratorType, DataType > &data_storage)
const unsigned int n_q_points
FullMatrix< double > matrix_dofs_child
FullMatrix< double > project_to_fe_matrix
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
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcCellDataTypeMismatch()
#define AssertDimension(dim1, dim2)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcTriangulationMismatch()
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
const ::Triangulation< dim, spacedim > & tria