16 #ifndef dealii_quadrature_point_data_h 17 #define dealii_quadrature_point_data_h 19 #include <deal.II/base/config.h> 21 #include <deal.II/base/quadrature.h> 22 #include <deal.II/base/subscriptor.h> 23 #include <deal.II/grid/tria.h> 24 #include <deal.II/grid/tria_accessor.h> 25 #include <deal.II/distributed/tria.h> 26 #include <deal.II/fe/fe_tools.h> 27 #include <deal.II/fe/fe.h> 28 #include <deal.II/lac/vector.h> 32 #include <type_traits> 34 DEAL_II_NAMESPACE_OPEN
54 template <
typename CellIteratorType,
typename DataType>
89 template <
typename T=DataType>
91 const unsigned int number_of_data_points_per_cell);
98 template <
typename T=DataType>
99 void initialize(
const CellIteratorType &cell_start,
100 const CellIteratorType &cell_end,
101 const unsigned int number_of_data_points_per_cell);
112 bool erase(
const CellIteratorType &cell);
129 template <
typename T=DataType>
130 std::vector<std::shared_ptr<T> >
get_data(
const CellIteratorType &cell);
142 template <
typename T=DataType>
143 std::vector<std::shared_ptr<const T> >
get_data(
const CellIteratorType &cell)
const;
149 std::map<CellIteratorType, std::vector<std::shared_ptr<DataType>>>
map;
200 virtual void pack_values(std::vector<double> &values)
const = 0;
209 virtual void unpack_values(
const std::vector<double> &values) = 0;
213 #ifdef DEAL_II_WITH_P4EST 216 namespace distributed
320 template <
int dim,
typename DataType>
324 static_assert(std::is_base_of<TransferableQuadraturePointData, DataType>::value,
325 "User's DataType class should be derived from TransferableQuadraturePointData");
383 const typename parallel::distributed::Triangulation<dim,dim>::CellStatus status,
392 const typename parallel::distributed::Triangulation<dim,dim>::CellStatus status,
473 template <
typename CellIteratorType,
typename DataType>
474 template <
typename T>
476 const unsigned int n_q_points)
478 static_assert(std::is_base_of<DataType, T>::value,
479 "User's T class should be derived from user's DataType class");
481 if (map.find(cell) == map.end())
483 map[cell] = std::vector<std::shared_ptr<DataType> >(n_q_points);
486 auto it = map.find(cell);
487 for (
unsigned int q=0; q < n_q_points; q++)
488 it->second[q] = std::make_shared<T>();
494 template <
typename CellIteratorType,
typename DataType>
495 template <
typename T>
497 const CellIteratorType &cell_end,
498 const unsigned int number)
500 for (CellIteratorType it = cell_start; it != cell_end; it++)
501 if (it->is_locally_owned())
502 initialize<T>(it,number);
507 template <
typename CellIteratorType,
typename DataType>
510 const auto it = map.find(cell);
511 for (
unsigned int i = 0; i < it->second.size(); i++)
513 Assert(it->second[i].unique(),
514 ExcMessage(
"Can not erase the cell data multiple objects reference its data."));
517 return (map.erase(cell)==1);
523 template <
typename CellIteratorType,
typename DataType>
529 auto it = map.begin();
530 while (it != map.end())
533 for (
unsigned int i = 0; i < it->second.size(); i++)
535 Assert(it->second[i].unique(),
536 ExcMessage(
"Can not erase the cell data, multiple objects reference it."));
545 template <
typename CellIteratorType,
typename DataType>
546 template <
typename T>
547 std::vector<std::shared_ptr<T> >
550 static_assert(std::is_base_of<DataType, T>::value,
551 "User's T class should be derived from user's DataType class");
553 auto it = map.find(cell);
560 std::vector<std::shared_ptr<T>> res (it->second.size());
561 for (
unsigned int q = 0; q < res.size(); q++)
562 res[q] = std::dynamic_pointer_cast<T>(it->second[q]);
569 template <
typename CellIteratorType,
typename DataType>
570 template <
typename T>
571 std::vector<std::shared_ptr<const T> >
574 static_assert(std::is_base_of<DataType, T>::value,
575 "User's T class should be derived from user's DataType class");
577 auto it = map.find(cell);
583 std::vector<std::shared_ptr<const T>> res (it->second.size());
584 for (
unsigned int q = 0; q < res.size(); q++)
585 res[q] = std::dynamic_pointer_cast<const T>(it->second[q]);
601 template <
typename CellIteratorType,
typename DataType>
603 (
const CellIteratorType &cell,
607 static_assert(std::is_base_of<TransferableQuadraturePointData, DataType>::value,
608 "User's DataType class should be derived from QPData");
610 const std::vector<std::shared_ptr<const DataType> > qpd = data_storage->
get_data(cell);
612 const unsigned int n = matrix_data.
n();
614 std::vector<double> single_qp_data(n);
615 Assert (qpd.size() == matrix_data.
m(),
617 for (
unsigned int q = 0; q < qpd.size(); q++)
619 qpd[q]->pack_values(single_qp_data);
620 Assert (single_qp_data.size() == n,
623 for (
unsigned int i = 0; i < n; i++)
624 matrix_data(q,i) = single_qp_data[i];
634 template <
typename CellIteratorType,
typename DataType>
635 void unpack_to_cell_data
636 (
const CellIteratorType &cell,
640 static_assert(std::is_base_of<TransferableQuadraturePointData, DataType>::value,
641 "User's DataType class should be derived from QPData");
643 std::vector<std::shared_ptr<DataType> > qpd = data_storage->
get_data(cell);
645 const unsigned int n = values_at_qp.
n();
647 std::vector<double> single_qp_data(n);
648 Assert(qpd.size() == values_at_qp.
m(),
651 for (
unsigned int q = 0; q < qpd.size(); q++)
653 for (
unsigned int i = 0; i < n; i++)
654 single_qp_data[i] = values_at_qp(q,i);
655 qpd[q]->unpack_values(single_qp_data);
660 #ifdef DEAL_II_WITH_P4EST 664 namespace distributed
666 template <
int dim,
typename DataType>
672 projection_fe(
std::unique_ptr<const
FiniteElement<dim> >(projection_fe_.clone())),
673 data_size_in_bytes(0),
674 n_q_points(rhs_quadrature.size()),
675 project_to_fe_matrix(projection_fe->dofs_per_cell,n_q_points),
676 project_to_qp_matrix(n_q_points,projection_fe->dofs_per_cell),
678 data_storage(nullptr),
679 triangulation(nullptr)
682 ExcMessage(
"ContinuousQuadratureDataTransfer requires scalar FiniteElement"));
698 template <
int dim,
typename DataType>
704 template <
int dim,
typename DataType>
709 Assert (data_storage ==
nullptr,
710 ExcMessage(
"This function can be called only once"));
711 triangulation = &tr_;
712 data_storage = &data_storage_;
714 unsigned int number_of_values = 0;
718 it != triangulation->end(); it++)
719 if (it->is_locally_owned())
721 std::vector<std::shared_ptr<DataType> > qpd = data_storage->
get_data(it);
722 number_of_values = qpd[0]->number_of_values();
726 number_of_values =
Utilities::MPI::max(number_of_values, triangulation->get_communicator ());
727 Assert (number_of_values > 0,
729 const unsigned int dofs_per_cell = projection_fe->dofs_per_cell;
730 matrix_dofs.reinit(dofs_per_cell,number_of_values);
731 matrix_dofs_child.reinit(dofs_per_cell,number_of_values);
732 matrix_quadrature.reinit(n_q_points,number_of_values);
733 data_size_in_bytes =
sizeof(double) * dofs_per_cell * number_of_values;
735 offset = triangulation->register_data_attach(data_size_in_bytes,
738 std::placeholders::_1,
739 std::placeholders::_2,
740 std::placeholders::_3));
745 template <
int dim,
typename DataType>
748 triangulation->notify_ready_to_unpack(offset,
751 std::placeholders::_1,
752 std::placeholders::_2,
753 std::placeholders::_3));
756 data_storage =
nullptr;
757 triangulation =
nullptr;
762 template <
int dim,
typename DataType>
765 const typename parallel::distributed::Triangulation<dim,dim>::CellStatus ,
768 double *data_store =
reinterpret_cast<double *
>(data);
770 pack_cell_data(cell,data_storage,matrix_quadrature);
773 project_to_fe_matrix.mmult(matrix_dofs, matrix_quadrature);
775 std::memcpy(data_store, &matrix_dofs(0,0), data_size_in_bytes);
780 template <
int dim,
typename DataType>
783 const typename parallel::distributed::Triangulation<dim,dim>::CellStatus status,
790 const double *data_store =
reinterpret_cast<const double *
>(data);
792 std::memcpy(&matrix_dofs(0,0), data_store, data_size_in_bytes);
794 if (cell->has_children())
798 for (
unsigned int child=0; child<cell->n_children(); ++child)
799 if (cell->child(child)->is_locally_owned())
801 projection_fe->get_prolongation_matrix(child, cell->refinement_case())
802 .mmult (matrix_dofs_child, matrix_dofs);
805 project_to_qp_matrix.mmult(matrix_quadrature,matrix_dofs_child);
808 unpack_to_cell_data(cell->child(child),matrix_quadrature,data_storage);
815 project_to_qp_matrix.mmult(matrix_quadrature,matrix_dofs);
818 unpack_to_cell_data(cell,matrix_quadrature,data_storage);
826 #endif // DEAL_II_WITH_P4EST 829 DEAL_II_NAMESPACE_CLOSE
CellDataStorage()=default
std::map< CellIteratorType, std::vector< std::shared_ptr< DataType > > > map
virtual unsigned int number_of_values() const =0
std::vector< std::shared_ptr< T > > get_data(const CellIteratorType &cell)
~CellDataStorage()=default
TransferableQuadraturePointData()=default
bool erase(const CellIteratorType &cell)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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)
void prepare_for_coarsening_and_refinement(parallel::distributed::Triangulation< dim > &tria, CellDataStorage< CellIteratorType, DataType > &data_storage)
virtual void unpack_values(const std::vector< double > &values)=0
std::size_t data_size_in_bytes
void unpack_function(const typename parallel::distributed::Triangulation< dim, dim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim, dim >::CellStatus status, const void *data)
::Triangulation< dim, spacedim >::active_cell_iterator active_cell_iterator
FullMatrix< double > project_to_qp_matrix
parallel::distributed::Triangulation< dim >::cell_iterator CellIteratorType
virtual ~TransferableQuadraturePointData()=default
void pack_function(const typename parallel::distributed::Triangulation< dim, dim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim, dim >::CellStatus status, void *data)
parallel::distributed::Triangulation< dim > * triangulation
::Triangulation< dim, spacedim >::cell_iterator cell_iterator
static ::ExceptionBase & ExcNotImplemented()
CellDataStorage< CellIteratorType, DataType > * data_storage
FullMatrix< double > matrix_quadrature
FullMatrix< double > matrix_dofs_child
FullMatrix< double > project_to_fe_matrix
FullMatrix< double > matrix_dofs
virtual void pack_values(std::vector< double > &values) const =0
T max(const T &t, const MPI_Comm &mpi_communicator)
void initialize(const CellIteratorType &cell, const unsigned int number_of_data_points_per_cell)
~ContinuousQuadratureDataTransfer()
const unsigned int n_q_points
static ::ExceptionBase & ExcInternalError()