Reference documentation for deal.II version 9.0.0
|
#include <deal.II/base/quadrature_point_data.h>
Public Types | |
typedef parallel::distributed::Triangulation< dim >::cell_iterator | CellIteratorType |
Public Member Functions | |
ContinuousQuadratureDataTransfer (const FiniteElement< dim > &projection_fe, const Quadrature< dim > &mass_quadrature, const Quadrature< dim > &data_quadrature) | |
~ContinuousQuadratureDataTransfer () | |
void | prepare_for_coarsening_and_refinement (parallel::distributed::Triangulation< dim > &tria, CellDataStorage< CellIteratorType, DataType > &data_storage) |
void | interpolate () |
Private Member Functions | |
void | pack_function (const typename parallel::distributed::Triangulation< dim, dim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim, dim >::CellStatus status, void *data) |
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) |
Private Attributes | |
const std::unique_ptr< const FiniteElement< dim > > | projection_fe |
std::size_t | data_size_in_bytes |
const unsigned int | n_q_points |
FullMatrix< double > | project_to_fe_matrix |
FullMatrix< double > | project_to_qp_matrix |
FullMatrix< double > | matrix_dofs |
FullMatrix< double > | matrix_dofs_child |
FullMatrix< double > | matrix_quadrature |
unsigned int | offset |
CellDataStorage< CellIteratorType, DataType > * | data_storage |
parallel::distributed::Triangulation< dim > * | triangulation |
A class for the transfer of continuous data stored at quadrature points when performing h-adaptive refinement of parallel::distributed::Triangulation .
This class implements the transfer of the quadrature point data between cells in case of adaptive refinement using L2 projection. That also includes automatic shipping of information between different processors.
To that end, the constructor of the class is provided with three main objects: scalar FiniteElement projection_fe
, mass_quadrature
and data_quadrature
Quadrature rules. First, the data located at data_quadrature
of each cell is L2-projected to the continuous space defined by a single FiniteElement projection_fe
. This is achieved using FETools::compute_projection_from_quadrature_points_matrix(). In doing so the mass matrix of this element is required, which will be calculated with the mass_quadrature
rule . Should the cell now belong to another processor, the data is then sent to this processor. The class makes use of a feature of p4est (and parallel::distributed::Triangulation) that allows one to attach information to cells during mesh refinement and rebalancing. On receiving information on the target cell, the data is projected back to the quadrature points using the matrix calculated by FETools::compute_interpolation_to_quadrature_points_matrix() . In the case that local refinement is performed, this class first project local DoF values of the parent element to each child.
This class is templated by DataType
type, however the user's DataType
class has to be derived from the TransferableQuadraturePointData class. In practice that amounts to implementing the following three functions shown below for a quadrature point data with 2 scalars:
Note that the order of packing and unpacking has to be the same.
This class can then be use with CellDataStorage in the following way:
This approach can be extended to quadrature point data with Tensors of arbitrary order, although with a little bit more work in packing and unpacking of data inside MyQData class.
Definition at line 321 of file quadrature_point_data.h.
typedef parallel::distributed::Triangulation<dim>::cell_iterator parallel::distributed::ContinuousQuadratureDataTransfer< dim, DataType >::CellIteratorType |
A typedef for a cell.
Definition at line 325 of file quadrature_point_data.h.
parallel::distributed::ContinuousQuadratureDataTransfer< dim, DataType >::ContinuousQuadratureDataTransfer | ( | const FiniteElement< dim > & | projection_fe, |
const Quadrature< dim > & | mass_quadrature, | ||
const Quadrature< dim > & | data_quadrature | ||
) |
Constructor which takes the FiniteElement projection_fe
, the quadrature rule mass_quadrature
used to integrate its local mass matrix and finally the quadrature rule data_quadrature
which is used to store DataType
.
projection_fe
has to be scalar-valued.projection_fe
is only required to be continuous within the cell. parallel::distributed::ContinuousQuadratureDataTransfer< dim, DataType >::~ContinuousQuadratureDataTransfer | ( | ) |
Destructor.
void parallel::distributed::ContinuousQuadratureDataTransfer< dim, DataType >::prepare_for_coarsening_and_refinement | ( | parallel::distributed::Triangulation< dim > & | tria, |
CellDataStorage< CellIteratorType, DataType > & | data_storage | ||
) |
Prepare for coarsening and refinement of a triangulation tria
. data_storage
represents the cell data which should be transferred and it should be initialized for each locally owned active cell.
data_storage
contains objects of the same type, more specifically they pack/unpack the same data. void parallel::distributed::ContinuousQuadratureDataTransfer< dim, DataType >::interpolate | ( | ) |
Interpolate the data previously stored in this object before the mesh was refined or coarsened onto the quadrature points of the currently active set of cells.
data_storage
object provided to prepare_for_coarsening_and_refinement() at new cells using CellDataStorage::initialize(). If that is not the case, an exception will be thrown in debug mode.
|
private |
A callback function used to pack the data on the current mesh into objects that can later be retrieved after refinement, coarsening and repartitioning.
|
private |
A callback function used to unpack the data on the current mesh that has been packed up previously on the mesh before refinement, coarsening and repartitioning.
|
private |
FiniteElement used to project data from and to quadrature points.
Definition at line 398 of file quadrature_point_data.h.
|
private |
The size of the data that will be sent, which depends on the DataType class.
Definition at line 403 of file quadrature_point_data.h.
|
private |
Number of quadrature points at which DataType is stored.
Definition at line 408 of file quadrature_point_data.h.
|
private |
Projection matrix from the quadrature points to local DoFs for a single scalar.
Definition at line 414 of file quadrature_point_data.h.
|
private |
Projection matrix from the local DoFs to quadrature points for a single scalar.
Definition at line 420 of file quadrature_point_data.h.
|
private |
Auxiliary matrix which represents projection of each internal value stored at the quadrature point (second index) to the local DoFs of the projection_fe
(first index).
Definition at line 427 of file quadrature_point_data.h.
|
private |
Projection of matrix_dofs
to each child cell in case of adaptive refinement.
Definition at line 432 of file quadrature_point_data.h.
|
private |
Auxiliary matrix which represents data (second index) stored at each quadrature point (first index).
Definition at line 438 of file quadrature_point_data.h.
|
private |
The offset that the parallel::distributed::Triangulation has assigned to this object starting at which we are allowed to read.
Definition at line 444 of file quadrature_point_data.h.
|
private |
A pointer to the CellDataStorage class whose data will be transferred.
Definition at line 449 of file quadrature_point_data.h.
|
private |
A pointer to the distributed triangulation to which cell data is attached.
Definition at line 454 of file quadrature_point_data.h.