Reference documentation for deal.II version GIT relicensing-480-geae235577e 2024-04-25 01:00:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Private Types | Private Attributes | List of all members
parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim > Class Template Reference

#include <deal.II/distributed/field_transfer.h>

Public Member Functions

 FieldTransfer (const DoFHandler< dim, spacedim > &dof_handler)
 
void prepare_for_coarsening_and_refinement (const VectorType &in, const unsigned int fe_nothing_index)
 
void interpolate (const Number &new_value, const AffineConstraints< Number > &affine_constraints, VectorType &out)
 

Private Types

using Number = typename VectorType::value_type
 

Private Attributes

const DoFHandler< dim, spacedim > & dof_handler
 
std::vector< Vector< Number > > data_to_transfer
 
std::unique_ptr< CellDataTransfer< dim, spacedim, std::vector< Vector< Number > > > > cell_data_transfer
 

Detailed Description

template<int dim, typename VectorType, int spacedim = dim>
class parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim >

This class is similar to SolutionTransfer but it supports the case where elements have been activated during refinement, i.e., FE_Nothing elements have been associated with a finite elements during refinement.

Definition at line 44 of file field_transfer.h.

Member Typedef Documentation

◆ Number

template<int dim, typename VectorType , int spacedim = dim>
using parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim >::Number = typename VectorType::value_type
private

Definition at line 47 of file field_transfer.h.

Constructor & Destructor Documentation

◆ FieldTransfer()

template<int dim, typename VectorType , int spacedim>
parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim >::FieldTransfer ( const DoFHandler< dim, spacedim > &  dof_handler)

Constructor.

Parameters
[in]dof_handlerThe DoFHandler on which all the operations will happen. This constructor must be called before the underlying Triangulation is coarsened/refined.

Definition at line 32 of file field_transfer.cc.

Member Function Documentation

◆ prepare_for_coarsening_and_refinement()

template<int dim, typename VectorType , int spacedim>
void parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim >::prepare_for_coarsening_and_refinement ( const VectorType &  in,
const unsigned int  fe_nothing_index 
)

Prepare the current object for coarsening and refinement.

Parameters
[in]inThe vector that will be interpolated
[in]fe_nothing_indexThe finite element index associated with FE_Nothing

Definition at line 132 of file field_transfer.cc.

◆ interpolate()

template<int dim, typename VectorType , int spacedim>
void parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim >::interpolate ( const Number new_value,
const AffineConstraints< Number > &  affine_constraints,
VectorType &  out 
)

Interpolate the data previously stored in this object before the mesh was refined or coarsened onto the current set of cells. new_value is the value associated to the new degrees of freedom that where created during the element activation. affine_constraints is the AffineConstraints after refinement.

Definition at line 173 of file field_transfer.cc.

Member Data Documentation

◆ dof_handler

template<int dim, typename VectorType , int spacedim = dim>
const DoFHandler<dim, spacedim>& parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim >::dof_handler
private

DoFHandler associated with the object.

Definition at line 86 of file field_transfer.h.

◆ data_to_transfer

template<int dim, typename VectorType , int spacedim = dim>
std::vector<Vector<Number> > parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim >::data_to_transfer
private

Data transferred by cell_data_transfer.

Definition at line 91 of file field_transfer.h.

◆ cell_data_transfer

template<int dim, typename VectorType , int spacedim = dim>
std::unique_ptr< CellDataTransfer<dim, spacedim, std::vector<Vector<Number> > > > parallel::distributed::experimental::FieldTransfer< dim, VectorType, spacedim >::cell_data_transfer
private

CellDataTransfer used to perform the field transfer.

Definition at line 98 of file field_transfer.h.


The documentation for this class was generated from the following files: