Reference documentation for deal.II version 9.1.1
|
#include <deal.II/meshworker/integration_info.h>
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. More... | |
const FEValuesBase< dim, spacedim > & | fe_values (const unsigned int i) const |
Access to finite elements. More... | |
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 |
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 |
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 |
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.
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.
Definition at line 78 of file integration_info.h.
|
inline |
Constructor.
Definition at line 584 of file integration_info.h.
|
inline |
Copy constructor, creating a clone to be used by WorkStream::run().
Definition at line 593 of file integration_info.h.
|
inline |
Build all internal structures, in particular the FEValuesBase objects and allocate space for data vectors.
el | is the finite element of the DoFHandler. |
mapping | is the Mapping object used to map the mesh cells. |
quadrature | is a Quadrature formula used in the constructor of the FEVALUES objects. |
flags | are the UpdateFlags used in the constructor of the FEVALUES objects. |
local_block_info | is 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. |
Definition at line 642 of file integration_info.h.
void MeshWorker::IntegrationInfo< dim, spacedim >::initialize_data | ( | const std::shared_ptr< VectorDataBase< dim, spacedim >> & | data | ) |
Initialize the data vector and cache the selector.
void MeshWorker::IntegrationInfo< dim, spacedim >::clear | ( | ) |
Delete the data created by initialize().
|
inline |
Return a reference to the FiniteElement that was used to initialize this object.
Definition at line 670 of file integration_info.h.
|
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 678 of file integration_info.h.
|
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 687 of file integration_info.h.
|
inline |
Reinitialize internal data structures for use on a cell.
Definition at line 697 of file integration_info.h.
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.
std::size_t MeshWorker::IntegrationInfo< dim, spacedim >::memory_consumption | ( | ) | const |
The memory used by this object.
|
private |
Use the finite element functions in global_data and fill the vectors values, gradients and hessians with values according to the selector.
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 170 of file integration_info.h.
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 180 of file integration_info.h.
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 190 of file integration_info.h.
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 212 of file integration_info.h.
|
private |
The pointer to the (system) element used for initialization.
Definition at line 226 of file integration_info.h.
|
private |
Cache the number of components of the system element.
Definition at line 241 of file integration_info.h.