Reference documentation for deal.II version 9.0.0
quadrature_point_data.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2018 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 at
12 // the top level of the deal.II distribution.
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 
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>
29 
30 #include <vector>
31 #include <map>
32 #include <type_traits>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
38 
54 template <typename CellIteratorType, typename DataType>
56 {
57 public:
61  CellDataStorage() = default;
62 
66  ~CellDataStorage() = default;
67 
89  template <typename T=DataType>
90  void initialize(const CellIteratorType &cell,
91  const unsigned int number_of_data_points_per_cell);
92 
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);
102 
112  bool erase(const CellIteratorType &cell);
113 
117  void clear();
118 
129  template <typename T=DataType>
130  std::vector<std::shared_ptr<T> > get_data(const CellIteratorType &cell);
131 
142  template <typename T=DataType>
143  std::vector<std::shared_ptr<const T> > get_data(const CellIteratorType &cell) const;
144 
145 private:
149  std::map<CellIteratorType, std::vector<std::shared_ptr<DataType>>> map;
150 
151 };
152 
153 
172 {
173 public:
178 
182  virtual ~TransferableQuadraturePointData() = default;
183 
189  virtual unsigned int number_of_values() const = 0;
190 
200  virtual void pack_values(std::vector<double> &values) const = 0;
201 
209  virtual void unpack_values(const std::vector<double> &values) = 0;
210 };
211 
212 
213 #ifdef DEAL_II_WITH_P4EST
214 namespace parallel
215 {
216  namespace distributed
217  {
320  template <int dim, typename DataType>
322  {
323  public:
324  static_assert(std::is_base_of<TransferableQuadraturePointData, DataType>::value,
325  "User's DataType class should be derived from TransferableQuadraturePointData");
326 
331 
343  const Quadrature<dim> &mass_quadrature,
344  const Quadrature<dim> &data_quadrature);
345 
350 
363 
374  void interpolate ();
375 
376  private:
383  const typename parallel::distributed::Triangulation<dim,dim>::CellStatus status,
384  void *data);
385 
392  const typename parallel::distributed::Triangulation<dim,dim>::CellStatus status,
393  const void *data);
394 
398  const std::unique_ptr<const FiniteElement<dim> > projection_fe;
399 
403  std::size_t data_size_in_bytes;
404 
408  const unsigned int n_q_points;
409 
415 
421 
428 
433 
439 
444  unsigned int offset;
445 
450 
455  };
456 
457  } // distributed
458 
459 } // parallel
460 
461 #endif
462 
465 #ifndef DOXYGEN
466 
467 // ------------------- inline and template functions ----------------
468 
469 //--------------------------------------------------------------------
470 // CellDataStorage
471 //--------------------------------------------------------------------
472 
473 template <typename CellIteratorType, typename DataType>
474 template <typename T>
475 void CellDataStorage<CellIteratorType,DataType>::initialize(const CellIteratorType &cell,
476  const unsigned int n_q_points)
477 {
478  static_assert(std::is_base_of<DataType, T>::value,
479  "User's T class should be derived from user's DataType class");
480 
481  if (map.find(cell) == map.end())
482  {
483  map[cell] = std::vector<std::shared_ptr<DataType> >(n_q_points);
484  // we need to initialize one-by-one as the std::vector<>(q, T())
485  // will end with a single same T object stored in each element of the vector:
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>();
489  }
490 }
491 
492 
493 
494 template <typename CellIteratorType, typename DataType>
495 template <typename T>
496 void CellDataStorage<CellIteratorType,DataType>::initialize(const CellIteratorType &cell_start,
497  const CellIteratorType &cell_end,
498  const unsigned int number)
499 {
500  for (CellIteratorType it = cell_start; it != cell_end; it++)
501  if (it->is_locally_owned())
502  initialize<T>(it,number);
503 }
504 
505 
506 
507 template <typename CellIteratorType, typename DataType>
508 bool CellDataStorage<CellIteratorType,DataType>::erase(const CellIteratorType &cell)
509 {
510  const auto it = map.find(cell);
511  for (unsigned int i = 0; i < it->second.size(); i++)
512  {
513  Assert(it->second[i].unique(),
514  ExcMessage("Can not erase the cell data multiple objects reference its data."));
515  }
516 
517  return (map.erase(cell)==1);
518 }
519 
520 
521 
522 
523 template <typename CellIteratorType, typename DataType>
525 {
526  // Do not call
527  // map.clear();
528  // as we want to be sure no one uses the stored objects. Loop manually:
529  auto it = map.begin();
530  while (it != map.end())
531  {
532  // loop over all objects and see if no one is using them
533  for (unsigned int i = 0; i < it->second.size(); i++)
534  {
535  Assert(it->second[i].unique(),
536  ExcMessage("Can not erase the cell data, multiple objects reference it."));
537  }
538  it = map.erase(it);
539  }
540 }
541 
542 
543 
544 
545 template <typename CellIteratorType, typename DataType>
546 template <typename T>
547 std::vector<std::shared_ptr<T> >
548 CellDataStorage<CellIteratorType,DataType>::get_data(const CellIteratorType &cell)
549 {
550  static_assert(std::is_base_of<DataType, T>::value,
551  "User's T class should be derived from user's DataType class");
552 
553  auto it = map.find(cell);
554  Assert(it != map.end(), ExcMessage("Could not find data for the cell"));
555 
556  // It would be nice to have a specialized version of this function for T==DataType.
557  // However explicit (i.e full) specialization of a member template is only allowed
558  // when the enclosing class is also explicitly (i.e fully) specialized.
559  // Thus, stick with copying of shared pointers even when the T==DataType:
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]);
563  return res;
564 }
565 
566 
567 
568 
569 template <typename CellIteratorType, typename DataType>
570 template <typename T>
571 std::vector<std::shared_ptr<const T> >
572 CellDataStorage<CellIteratorType,DataType>::get_data(const CellIteratorType &cell) const
573 {
574  static_assert(std::is_base_of<DataType, T>::value,
575  "User's T class should be derived from user's DataType class");
576 
577  auto it = map.find(cell);
578  Assert(it != map.end(), ExcMessage("Could not find QP data for the cell"));
579 
580  // Cast base class to the desired class. This has to be done irrespectively of
581  // T==DataType as we need to return shared_ptr<const T> to make sure the user
582  // does not modify the content of QP objects
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]);
586  return res;
587 }
588 
589 //--------------------------------------------------------------------
590 // ContinuousQuadratureDataTransfer
591 //--------------------------------------------------------------------
592 
593 
594 /*
595  * Pack cell data of type @p DataType stored using @p data_storage in @p cell
596  * at each quadrature point to @p matrix_data. Here @p matrix_data is a matrix
597  * whose first index corresponds to different quadrature points on the cell
598  * whereas the second index represents different values stored at each
599  * quadrature point in the DataType class.
600  */
601 template <typename CellIteratorType, typename DataType>
602 void pack_cell_data
603 (const CellIteratorType &cell,
604  const CellDataStorage<CellIteratorType,DataType> *data_storage,
605  FullMatrix<double> &matrix_data)
606 {
607  static_assert(std::is_base_of<TransferableQuadraturePointData, DataType>::value,
608  "User's DataType class should be derived from QPData");
609 
610  const std::vector<std::shared_ptr<const DataType> > qpd = data_storage->get_data(cell);
611 
612  const unsigned int n = matrix_data.n();
613 
614  std::vector<double> single_qp_data(n);
615  Assert (qpd.size() == matrix_data.m(),
616  ExcDimensionMismatch(qpd.size(), matrix_data.m()));
617  for (unsigned int q = 0; q < qpd.size(); q++)
618  {
619  qpd[q]->pack_values(single_qp_data);
620  Assert (single_qp_data.size() == n,
621  ExcDimensionMismatch(single_qp_data.size(),n));
622 
623  for (unsigned int i = 0; i < n; i++)
624  matrix_data(q,i) = single_qp_data[i];
625  }
626 }
627 
628 
629 
630 
631 /*
632  * the opposite of the pack function above.
633  */
634 template <typename CellIteratorType, typename DataType>
635 void unpack_to_cell_data
636 (const CellIteratorType &cell,
637  const FullMatrix<double> &values_at_qp,
639 {
640  static_assert(std::is_base_of<TransferableQuadraturePointData, DataType>::value,
641  "User's DataType class should be derived from QPData");
642 
643  std::vector<std::shared_ptr<DataType> > qpd = data_storage->get_data(cell);
644 
645  const unsigned int n = values_at_qp.n();
646 
647  std::vector<double> single_qp_data(n);
648  Assert(qpd.size() == values_at_qp.m(),
649  ExcDimensionMismatch(qpd.size(), values_at_qp.m()));
650 
651  for (unsigned int q = 0; q < qpd.size(); q++)
652  {
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);
656  }
657 }
658 
659 
660 #ifdef DEAL_II_WITH_P4EST
661 
662 namespace parallel
663 {
664  namespace distributed
665  {
666  template <int dim, typename DataType>
668  (const FiniteElement<dim> &projection_fe_,
669  const Quadrature<dim> &lhs_quadrature,
670  const Quadrature<dim> &rhs_quadrature)
671  :
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),
677  offset(0),
678  data_storage(nullptr),
679  triangulation(nullptr)
680  {
681  Assert(projection_fe->n_components() == 1,
682  ExcMessage("ContinuousQuadratureDataTransfer requires scalar FiniteElement"));
683 
685  (*projection_fe.get(),
686  lhs_quadrature,
687  rhs_quadrature,
689 
691  (*projection_fe.get(),
692  rhs_quadrature,
694  }
695 
696 
697 
698  template <int dim, typename DataType>
700  {}
701 
702 
703 
704  template <int dim, typename DataType>
708  {
709  Assert (data_storage == nullptr,
710  ExcMessage("This function can be called only once"));
711  triangulation = &tr_;
712  data_storage = &data_storage_;
713  // get the number from the first active cell
714  unsigned int number_of_values = 0;
715  // if triangulation has some active cells locally owned cells on this processor we can expect
716  // data to be initialized. Do that to get the number:
717  for (typename parallel::distributed::Triangulation<dim>::active_cell_iterator it = triangulation->begin_active();
718  it != triangulation->end(); it++)
719  if (it->is_locally_owned())
720  {
721  std::vector<std::shared_ptr<DataType> > qpd = data_storage->get_data(it);
722  number_of_values = qpd[0]->number_of_values();
723  break;
724  }
725  // some processors may have no data stored, thus get the maximum among all processors:
726  number_of_values = Utilities::MPI::max(number_of_values, triangulation->get_communicator ());
727  Assert (number_of_values > 0,
728  ExcInternalError());
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;
734 
735  offset = triangulation->register_data_attach(data_size_in_bytes,
737  this,
738  std::placeholders::_1,
739  std::placeholders::_2,
740  std::placeholders::_3));
741  }
742 
743 
744 
745  template <int dim, typename DataType>
747  {
748  triangulation->notify_ready_to_unpack(offset,
750  this,
751  std::placeholders::_1,
752  std::placeholders::_2,
753  std::placeholders::_3));
754 
755  // invalidate the pointers
756  data_storage = nullptr;
757  triangulation = nullptr;
758  }
759 
760 
761 
762  template <int dim, typename DataType>
765  const typename parallel::distributed::Triangulation<dim,dim>::CellStatus /*status*/,
766  void *data)
767  {
768  double *data_store = reinterpret_cast<double *>(data);
769 
770  pack_cell_data(cell,data_storage,matrix_quadrature);
771 
772  // project to FE
773  project_to_fe_matrix.mmult(matrix_dofs, matrix_quadrature);
774 
775  std::memcpy(data_store, &matrix_dofs(0,0), data_size_in_bytes);
776  }
777 
778 
779 
780  template <int dim, typename DataType>
783  const typename parallel::distributed::Triangulation<dim,dim>::CellStatus status,
784  const void *data)
785  {
788  (void) status;
789 
790  const double *data_store = reinterpret_cast<const double *>(data);
791 
792  std::memcpy(&matrix_dofs(0,0), data_store, data_size_in_bytes);
793 
794  if (cell->has_children())
795  {
796  // we need to first use prolongation matrix to get dofvalues on child
797  // cells based on dofvalues stored in the parent's data_store
798  for (unsigned int child=0; child<cell->n_children(); ++child)
799  if (cell->child(child)->is_locally_owned())
800  {
801  projection_fe->get_prolongation_matrix(child, cell->refinement_case())
802  .mmult (matrix_dofs_child, matrix_dofs);
803 
804  // now we do the usual business of evaluating FE on quadrature points:
805  project_to_qp_matrix.mmult(matrix_quadrature,matrix_dofs_child);
806 
807  // finally, put back into the map:
808  unpack_to_cell_data(cell->child(child),matrix_quadrature,data_storage);
809  }
810  }
811  else
812  {
813  // if there are no children, evaluate FE field at
814  // rhs_quadrature points.
815  project_to_qp_matrix.mmult(matrix_quadrature,matrix_dofs);
816 
817  // finally, put back into the map:
818  unpack_to_cell_data(cell,matrix_quadrature,data_storage);
819  }
820  }
821 
822  } // distributed
823 
824 } // parallel
825 
826 #endif // DEAL_II_WITH_P4EST
827 
828 #endif // DOXYGEN
829 DEAL_II_NAMESPACE_CLOSE
830 
831 #endif
size_type m() const
void compute_interpolation_to_quadrature_points_matrix(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, FullMatrix< double > &I_q)
CellDataStorage()=default
STL namespace.
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
bool erase(const CellIteratorType &cell)
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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)
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
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
Definition: tria.h:285
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
Definition: tria.h:265
static ::ExceptionBase & ExcNotImplemented()
CellDataStorage< CellIteratorType, DataType > * data_storage
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)
static ::ExceptionBase & ExcInternalError()