Reference documentation for deal.II version 9.5.0
\(\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 | Public Attributes | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
MeshWorker::IntegrationInfo< dim, spacedim > Class Template Reference

#include <deal.II/meshworker/integration_info.h>

Inheritance diagram for MeshWorker::IntegrationInfo< dim, spacedim >:
[legend]

Public Member Functions

 IntegrationInfo ()
 
 IntegrationInfo (const IntegrationInfo< dim, spacedim > &other)
 
template<class FEVALUES >
void initialize (const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const Quadrature< FEVALUES::integral_dimension > &quadrature, const UpdateFlags flags, const BlockInfo *local_block_info=nullptr)
 
void initialize_data (const std::shared_ptr< VectorDataBase< dim, spacedim > > &data)
 
void clear ()
 
const FiniteElement< dim, spacedim > & finite_element () const
 
const FEValuesBase< dim, spacedim > & fe_values () const
 Access to finite element.
 
const FEValuesBase< dim, spacedim > & fe_values (const unsigned int i) const
 Access to finite elements.
 
template<typename number >
void reinit (const DoFInfo< dim, spacedim, number > &i)
 
template<typename number >
void fill_local_data (const DoFInfo< dim, spacedim, number > &info, bool split_fevalues)
 
std::size_t memory_consumption () const
 
template<class FEVALUES >
void initialize (const FiniteElement< dim, sdim > &el, const Mapping< dim, sdim > &mapping, const Quadrature< FEVALUES::integral_dimension > &quadrature, const UpdateFlags flags, const BlockInfo *block_info)
 

Public Attributes

bool multigrid
 This is true if we are assembling for multigrid.
 
std::vector< std::vector< std::vector< double > > > values
 
std::vector< std::vector< std::vector< Tensor< 1, spacedim > > > > gradients
 
std::vector< std::vector< std::vector< Tensor< 2, spacedim > > > > hessians
 
std::shared_ptr< VectorDataBase< dim, spacedim > > global_data
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 
static constexpr unsigned int space_dimension = spacedim
 

Private Member Functions

template<typename TYPE >
void fill_local_data (std::vector< std::vector< std::vector< TYPE > > > &data, VectorSelector &selector, bool split_fevalues) const
 

Private Attributes

std::vector< std::shared_ptr< FEValuesBase< dim, spacedim > > > fevalv
 vector of FEValues objects
 
SmartPointer< const FiniteElement< dim, spacedim >, IntegrationInfo< dim, spacedim > > fe_pointer
 
unsigned int n_components
 

Detailed Description

template<int dim, int spacedim = dim>
class MeshWorker::IntegrationInfo< dim, spacedim >

Class for objects handed to local integration functions.

Objects of this class contain one or more objects of type FEValues, FEFaceValues or FESubfaceValues to be used in local integration. They are stored in an array of pointers to the base classes FEValuesBase. The template parameter VectorType allows the use of different data types for the global system.

Additionally, this function contains space to store the values of finite element functions stored in global_data in the quadrature points. These vectors are initialized automatically on each cell or face. In order to avoid initializing unused vectors, you can use initialize_selector() to select the vectors by name that you actually want to use.

Integration models

This class supports two local integration models, corresponding to the data models in the documentation of the Assembler namespace. One is the standard model suggested by the use of FESystem. Namely, there is one FEValuesBase object in this class, containing all shape functions of the whole system, and having as many components as the system. Using this model involves loops over all system shape functions. It requires to identify the system components for each shape function and to select the correct bilinear form, usually in an if or switch statement.

The second integration model builds one FEValuesBase object per base element of the system. The degrees of freedom on each cell are renumbered by block, such that they represent the same block structure as the global system. Objects performing the integration can then process each block separately, which improves reusability of code considerably.

Note
As described in DoFInfo, the use of the local block model is triggered by calling BlockInfo::initialize_local() before using initialize() in this class.

Definition at line 77 of file integration_info.h.

Constructor & Destructor Documentation

◆ IntegrationInfo() [1/2]

template<int dim, int sdim>
MeshWorker::IntegrationInfo< dim, sdim >::IntegrationInfo
inline

Constructor.

Definition at line 582 of file integration_info.h.

◆ IntegrationInfo() [2/2]

template<int dim, int spacedim = dim>
MeshWorker::IntegrationInfo< dim, spacedim >::IntegrationInfo ( const IntegrationInfo< dim, spacedim > &  other)

Copy constructor, creating a clone to be used by WorkStream::run().

Member Function Documentation

◆ initialize() [1/2]

template<int dim, int spacedim = dim>
template<class FEVALUES >
void MeshWorker::IntegrationInfo< dim, spacedim >::initialize ( const FiniteElement< dim, spacedim > &  el,
const Mapping< dim, spacedim > &  mapping,
const Quadrature< FEVALUES::integral_dimension > &  quadrature,
const UpdateFlags  flags,
const BlockInfo local_block_info = nullptr 
)

Build all internal structures, in particular the FEValuesBase objects and allocate space for data vectors.

Parameters
elis the finite element of the DoFHandler.
mappingis the Mapping object used to map the mesh cells.
quadratureis a Quadrature formula used in the constructor of the FEVALUES objects.
flagsare the UpdateFlags used in the constructor of the FEVALUES objects.
local_block_infois an optional parameter for systems of PDE. If it is provided with reasonable data, then the degrees of freedom on the cells will be re-ordered to reflect the block structure of the system.

◆ initialize_data()

template<int dim, int spacedim = dim>
void MeshWorker::IntegrationInfo< dim, spacedim >::initialize_data ( const std::shared_ptr< VectorDataBase< dim, spacedim > > &  data)

Initialize the data vector and cache the selector.

◆ clear()

template<int dim, int spacedim = dim>
void MeshWorker::IntegrationInfo< dim, spacedim >::clear ( )

Delete the data created by initialize().

◆ finite_element()

template<int dim, int spacedim>
const FiniteElement< dim, spacedim > & MeshWorker::IntegrationInfo< dim, spacedim >::finite_element
inline

Return a reference to the FiniteElement that was used to initialize this object.

Definition at line 669 of file integration_info.h.

◆ fe_values() [1/2]

template<int dim, int spacedim>
const FEValuesBase< dim, spacedim > & MeshWorker::IntegrationInfo< dim, spacedim >::fe_values
inline

Access to finite element.

This is the access function being used, if initialize() for a single element was used (without the BlockInfo argument). It throws an exception, if applied to a vector of elements.

Definition at line 677 of file integration_info.h.

◆ fe_values() [2/2]

template<int dim, int spacedim>
const FEValuesBase< dim, spacedim > & MeshWorker::IntegrationInfo< dim, spacedim >::fe_values ( const unsigned int  i) const
inline

Access to finite elements.

This access function must be used if the initialize() for a group of elements was used (with a valid BlockInfo object).

Definition at line 686 of file integration_info.h.

◆ reinit()

template<int dim, int spacedim>
template<typename number >
void MeshWorker::IntegrationInfo< dim, spacedim >::reinit ( const DoFInfo< dim, spacedim, number > &  i)
inline

Reinitialize internal data structures for use on a cell.

Definition at line 696 of file integration_info.h.

◆ fill_local_data() [1/2]

template<int dim, int spacedim = dim>
template<typename number >
void MeshWorker::IntegrationInfo< dim, spacedim >::fill_local_data ( const DoFInfo< dim, spacedim, number > &  info,
bool  split_fevalues 
)

Use the finite element functions in global_data and fill the vectors values, gradients and hessians.

◆ memory_consumption()

template<int dim, int spacedim = dim>
std::size_t MeshWorker::IntegrationInfo< dim, spacedim >::memory_consumption ( ) const

The memory used by this object.

◆ fill_local_data() [2/2]

template<int dim, int spacedim = dim>
template<typename TYPE >
void MeshWorker::IntegrationInfo< dim, spacedim >::fill_local_data ( std::vector< std::vector< std::vector< TYPE > > > &  data,
VectorSelector selector,
bool  split_fevalues 
) const
private

Use the finite element functions in global_data and fill the vectors values, gradients and hessians with values according to the selector.

◆ initialize() [2/2]

template<int dim, int spacedim = dim>
template<class FEVALUES >
void MeshWorker::IntegrationInfo< dim, spacedim >::initialize ( const FiniteElement< dim, sdim > &  el,
const Mapping< dim, sdim > &  mapping,
const Quadrature< FEVALUES::integral_dimension > &  quadrature,
const UpdateFlags  flags,
const BlockInfo block_info 
)
inline

Definition at line 641 of file integration_info.h.

Member Data Documentation

◆ fevalv

template<int dim, int spacedim = dim>
std::vector<std::shared_ptr<FEValuesBase<dim, spacedim> > > MeshWorker::IntegrationInfo< dim, spacedim >::fevalv
private

vector of FEValues objects

Definition at line 81 of file integration_info.h.

◆ dimension

template<int dim, int spacedim = dim>
constexpr unsigned int MeshWorker::IntegrationInfo< dim, spacedim >::dimension = dim
staticconstexpr

Definition at line 84 of file integration_info.h.

◆ space_dimension

template<int dim, int spacedim = dim>
constexpr unsigned int MeshWorker::IntegrationInfo< dim, spacedim >::space_dimension = spacedim
staticconstexpr

Definition at line 85 of file integration_info.h.

◆ multigrid

template<int dim, int spacedim = dim>
bool MeshWorker::IntegrationInfo< dim, spacedim >::multigrid

This is true if we are assembling for multigrid.

Definition at line 143 of file integration_info.h.

◆ values

template<int dim, int spacedim = dim>
std::vector<std::vector<std::vector<double> > > MeshWorker::IntegrationInfo< dim, spacedim >::values

The vector containing the values of finite element functions in the quadrature points.

There is one vector per selected finite element function, containing one vector for each component, containing vectors with values for each quadrature point.

Definition at line 169 of file integration_info.h.

◆ gradients

template<int dim, int spacedim = dim>
std::vector<std::vector<std::vector<Tensor<1, spacedim> > > > MeshWorker::IntegrationInfo< dim, spacedim >::gradients

The vector containing the derivatives of finite element functions in the quadrature points.

There is one vector per selected finite element function, containing one vector for each component, containing vectors with values for each quadrature point.

Definition at line 179 of file integration_info.h.

◆ hessians

template<int dim, int spacedim = dim>
std::vector<std::vector<std::vector<Tensor<2, spacedim> > > > MeshWorker::IntegrationInfo< dim, spacedim >::hessians

The vector containing the second derivatives of finite element functions in the quadrature points.

There is one vector per selected finite element function, containing one vector for each component, containing vectors with values for each quadrature point.

Definition at line 189 of file integration_info.h.

◆ global_data

template<int dim, int spacedim = dim>
std::shared_ptr<VectorDataBase<dim, spacedim> > MeshWorker::IntegrationInfo< dim, spacedim >::global_data

The global data vector used to compute function values in quadrature points.

Definition at line 211 of file integration_info.h.

◆ fe_pointer

template<int dim, int spacedim = dim>
SmartPointer<const FiniteElement<dim, spacedim>, IntegrationInfo<dim, spacedim> > MeshWorker::IntegrationInfo< dim, spacedim >::fe_pointer
private

The pointer to the (system) element used for initialization.

Definition at line 225 of file integration_info.h.

◆ n_components

template<int dim, int spacedim = dim>
unsigned int MeshWorker::IntegrationInfo< dim, spacedim >::n_components
private

Cache the number of components of the system element.

Definition at line 240 of file integration_info.h.


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