16#ifndef dealii_quadrature_point_data_h
17#define dealii_quadrature_point_data_h
63template <
typename CellIteratorType,
typename DataType>
102 template <
typename T = DataType>
105 const unsigned int number_of_data_points_per_cell);
112 template <
typename T = DataType>
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>>
168 template <
typename T = DataType>
169 std::vector<std::shared_ptr<const T>>
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.");
324#ifdef DEAL_II_WITH_P4EST
327 namespace distributed
438 template <
int dim,
typename DataType>
443 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
444 "User's DataType class should be derived from TransferableQuadraturePointData");
518 const boost::iterator_range<std::vector<char>::const_iterator>
602template <
typename CellIteratorType,
typename DataType>
606 const CellIteratorType &cell,
607 const unsigned int n_q_points)
609 static_assert(std::is_base_of<DataType, T>::value,
610 "User's T class should be derived from user's DataType class");
615 Assert(&cell->get_triangulation() ==
tria, ExcTriangulationMismatch());
617 const auto key = cell->id();
618 if (map.find(key) == map.end())
620 map[key] = std::vector<std::shared_ptr<DataType>>(n_q_points);
624 const auto it = map.find(key);
625 for (
unsigned int q = 0; q < n_q_points; ++q)
626 it->second[q] = std::make_shared<T>();
632template <
typename CellIteratorType,
typename DataType>
636 const CellIteratorType &cell_start,
637 const CellIteratorType &cell_end,
638 const unsigned int number)
640 for (CellIteratorType it = cell_start; it != cell_end; ++it)
641 if (it->is_locally_owned())
642 initialize<T>(it, number);
647template <
typename CellIteratorType,
typename DataType>
651 const auto key = cell->id();
652 const auto it = map.find(key);
655 Assert(&cell->get_triangulation() ==
tria, ExcTriangulationMismatch());
656 for (
unsigned int i = 0; i < it->second.size(); ++i)
659 it->second[i].unique(),
661 "Can not erase the cell data multiple objects reference its data."));
664 return (map.erase(key) == 1);
669template <
typename CellIteratorType,
typename DataType>
676 auto it = map.begin();
677 while (it != map.end())
680 for (
unsigned int i = 0; i < it->second.size(); ++i)
683 it->second[i].unique(),
685 "Can not erase the cell data, multiple objects reference it."));
693template <
typename CellIteratorType,
typename DataType>
695inline std::vector<std::shared_ptr<T>>
697 const CellIteratorType &cell)
699 static_assert(std::is_base_of<DataType, T>::value,
700 "User's T class should be derived from user's DataType class");
701 Assert(&cell->get_triangulation() ==
tria, ExcTriangulationMismatch());
703 const auto it = map.find(cell->id());
711 std::vector<std::shared_ptr<T>> res(it->second.size());
712 for (
unsigned int q = 0; q < res.size(); ++q)
714 res[q] = std::dynamic_pointer_cast<T>(it->second[q]);
715 Assert(res[q], ExcCellDataTypeMismatch());
722template <
typename CellIteratorType,
typename DataType>
724inline std::vector<std::shared_ptr<const T>>
726 const CellIteratorType &cell)
const
728 static_assert(std::is_base_of<DataType, T>::value,
729 "User's T class should be derived from user's DataType class");
730 Assert(&cell->get_triangulation() ==
tria, ExcTriangulationMismatch());
732 const auto it = map.find(cell->id());
738 std::vector<std::shared_ptr<const T>> res(it->second.size());
739 for (
unsigned int q = 0; q < res.size(); ++q)
741 res[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
742 Assert(res[q], ExcCellDataTypeMismatch());
747template <
typename CellIteratorType,
typename DataType>
749inline std_cxx17::optional<std::vector<std::shared_ptr<T>>>
751 const CellIteratorType &cell)
753 static_assert(std::is_base_of<DataType, T>::value,
754 "User's T class should be derived from user's DataType class");
755 Assert(&cell->get_triangulation() ==
tria, ExcTriangulationMismatch());
757 const auto it = map.find(cell->id());
764 std::vector<std::shared_ptr<T>> result(it->second.size());
765 for (
unsigned int q = 0; q < result.size(); ++q)
767 result[q] = std::dynamic_pointer_cast<T>(it->second[q]);
768 Assert(result[q], ExcCellDataTypeMismatch());
778template <
typename CellIteratorType,
typename DataType>
780inline std_cxx17::optional<std::vector<std::shared_ptr<const T>>>
782 const CellIteratorType &cell)
const
784 static_assert(std::is_base_of<DataType, T>::value,
785 "User's T class should be derived from user's DataType class");
786 Assert(&cell->get_triangulation() ==
tria, ExcTriangulationMismatch());
788 const auto it = map.find(cell->id());
795 std::vector<std::shared_ptr<const T>> result(it->second.size());
796 for (
unsigned int q = 0; q < result.size(); ++q)
798 result[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
799 Assert(result[q], ExcCellDataTypeMismatch());
821template <
typename CellIteratorType,
typename DataType>
823pack_cell_data(
const CellIteratorType & cell,
828 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
829 "User's DataType class should be derived from QPData");
833 const unsigned int m = qpd->size();
835 const unsigned int n = (*qpd)[0]->number_of_values();
836 matrix_data.reinit(m, n);
838 std::vector<double> single_qp_data(n);
839 for (
unsigned int q = 0; q < m; ++q)
841 (*qpd)[q]->pack_values(single_qp_data);
844 for (
unsigned int i = 0; i < n; ++i)
845 matrix_data(q, i) = single_qp_data[i];
850 matrix_data.reinit({0, 0});
859template <
typename CellIteratorType,
typename DataType>
861unpack_to_cell_data(
const CellIteratorType & cell,
866 std::is_base_of<TransferableQuadraturePointData, DataType>::value,
867 "User's DataType class should be derived from QPData");
871 const unsigned int n = values_at_qp.
n();
874 std::vector<double> single_qp_data(n);
877 for (
unsigned int q = 0; q < qpd->size(); ++q)
879 for (
unsigned int i = 0; i < n; ++i)
880 single_qp_data[i] = values_at_qp(q, i);
881 (*qpd)[q]->unpack_values(single_qp_data);
887# ifdef DEAL_II_WITH_P4EST
891 namespace distributed
893 template <
int dim,
typename DataType>
900 , data_size_in_bytes(0)
901 , n_q_points(rhs_quadrature.size())
902 , project_to_fe_matrix(projection_fe->n_dofs_per_cell(), n_q_points)
903 , project_to_qp_matrix(n_q_points, projection_fe->n_dofs_per_cell())
904 , handle(
numbers::invalid_unsigned_int)
905 , data_storage(nullptr)
911 "ContinuousQuadratureDataTransfer requires scalar FiniteElement"));
925 template <
int dim,
typename DataType>
932 Assert(data_storage ==
nullptr,
933 ExcMessage(
"This function can be called only once"));
935 data_storage = &data_storage_;
940 dim>::cell_iterator &cell,
942 status) {
return this->pack_function(cell, status); },
948 template <
int dim,
typename DataType>
956 dim>::cell_iterator &cell,
959 const boost::iterator_range<std::vector<char>::const_iterator>
960 &data_range) { this->unpack_function(cell, status, data_range); });
963 data_storage =
nullptr;
969 template <
int dim,
typename DataType>
970 inline std::vector<char>
977 pack_cell_data(cell, data_storage, matrix_quadrature);
980 const unsigned int number_of_values = matrix_quadrature.n();
981 matrix_dofs.reinit(project_to_fe_matrix.m(), number_of_values);
982 if (number_of_values > 0)
983 project_to_fe_matrix.mmult(matrix_dofs, matrix_quadrature);
990 template <
int dim,
typename DataType>
997 const boost::iterator_range<std::vector<char>::const_iterator>
1006 Utilities::unpack<FullMatrix<double>>(data_range.begin(),
1009 const unsigned int number_of_values = matrix_dofs.n();
1010 if (number_of_values == 0)
1013 matrix_quadrature.reinit(n_q_points, number_of_values);
1015 if (cell->has_children())
1019 matrix_dofs_child.reinit(projection_fe->n_dofs_per_cell(),
1021 for (
unsigned int child = 0; child < cell->n_children(); ++child)
1022 if (cell->child(child)->is_locally_owned())
1025 ->get_prolongation_matrix(child, cell->refinement_case())
1026 .mmult(matrix_dofs_child, matrix_dofs);
1030 project_to_qp_matrix.mmult(matrix_quadrature,
1034 unpack_to_cell_data(cell->child(child),
1043 project_to_qp_matrix.mmult(matrix_quadrature, matrix_dofs);
1046 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