OpenFCST: The open-source Fuel Cell Simulation Toolbox
|
This class is created for the objects handed to the local integration functions. More...
#include <mesh_worker_info.h>
Public Member Functions | |
IntegrationInfo (const BlockInfo &block_info) | |
Constructor. | |
IntegrationInfo (const FEVectors &data, const BlockInfo &block_info) | |
template<typename FEVALUES > | |
void | initialize (const FEVALUES *FE_VALUES, const FiniteElement< dim > &fe, const Mapping< dim > &mapping, const Quadrature< FEVALUES::integral_dimension > &quadrature, const UpdateFlags flags) |
Build internal structures like fevalv and allocate memory for values, gradients, hessians. | |
void | initialize_selector (const VectorSelector &s) |
Initialize the selector. | |
void | initialize_data (const FEVectors &data) |
Initialize the global_data. | |
void | clear () |
Resize fevalv to 0. | |
const FEVALUESBASE & | fe () const |
Access to a single actual FEVALUES object. | |
const FEVALUESBASE & | fe (unsigned int i) const |
Access to a group of actual FEVALUES objects. | |
const FEVALUESBASE & | get_fe_val_unsplit () const |
(ZVN) Access to fe_val_unsplit . | |
template<typename DHCellIterator > | |
void | reinit (const DHCellIterator &c) |
Reinitialize internal data for use on a cell. | |
template<typename DHCellIterator , typename DHFaceIterator > | |
void | reinit (const DHCellIterator &c, const DHFaceIterator &f, const unsigned int fn) |
Reinitialize internal data for use on a face. | |
template<typename DHCellIterator , typename DHFaceIterator > | |
void | reinit (const DHCellIterator &c, const DHFaceIterator &f, const unsigned int fn, const unsigned int sn) |
Reinitialize internal data for use on a subface. | |
template<typename TYPE > | |
void | fill_local_data (std::vector< std::vector< std::vector< TYPE > > > &data, bool split_fevalues) const |
Public Member Functions inherited from MeshWorker::InfoObjects::DoFInfo< dim > | |
DoFInfo (const BlockInfo &block_info) | |
Constructor. | |
DoFInfo (const FEVectors &, const BlockInfo &block_info) | |
void | reinit (const DHCellIterator &c) |
Set the current cell and fill #indices. | |
void | reinit (const DHCellIterator &c, const DHFaceIterator &f, const unsigned int n) |
Set the current cell and face and fill #indices if the #cell changed. | |
void | reinit (const DHCellIterator &c, const DHFaceIterator &f, const unsigned int n, const unsigned int s) |
Set the current cell, face, and subface and fill #indices if the #cell changed. | |
Public Attributes | |
bool | multigrid |
This is true if we assemble for multigrid. | |
SmartPointer< const FEVectors > | global_data |
The smart pointer to the FEVectors object called global_data. | |
VectorSelector | selector |
This object selects the named FEVector vectors in global_data and assigns the types of information extracted. | |
std::vector< std::vector < std::vector< double > > > | values |
The vector containing the values of finite element functions in the quadrature points. | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | gradients |
The vector containing the gradients of finite element functions in the quadrature points. | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > & | derivatives |
std::vector< std::vector < std::vector< Tensor< 2, dim > > > > | hessians |
The vector containing the hessians of finite element functions in the quadrature points. | |
Public Attributes inherited from MeshWorker::InfoObjects::DoFInfo< dim > | |
DoFHandler< dim >::cell_iterator | dof_cell |
The current DoFHandler<dim> cell. | |
DoFHandler< dim > ::active_cell_iterator | dof_active_cell |
The current DoFHandler<dim> active cell. | |
DoFHandler< dim >::face_iterator | dof_face |
The current DoFHandler<dim> face. | |
Triangulation< dim >::cell_iterator | cell |
The current Triangulation<dim> cell. | |
Triangulation< dim >::face_iterator | face |
The current Triangulation<dim> face. | |
unsigned int | face_number |
The number of the current face on the current cell. | |
unsigned int | sub_number |
The number of the current subface on the current face of the current cell. | |
std::vector< unsigned int > | indices |
The local DoF indices associated with the current mesh entity. | |
SmartPointer< const BlockInfo > | block_info |
The block structure of the problem at hand. | |
Private Member Functions | |
template<typename TYPE > | |
void | fill_local_data (std::vector< std::vector< std::vector< TYPE > > > &data, VectorSelector &selector, bool split_fevalues) const |
Use the finite element functions in global_data and fill the vectors values, gradients, hessians. | |
Private Attributes | |
std::vector< boost::shared_ptr < FEVALUESBASE > > | fevalv |
The vector of smart pointers to FEVALUESBASE objects. | |
boost::shared_ptr< FEVALUESBASE > | fe_val_unsplit |
(ZVN) In the mode # 1 (or in the main mode), AppFrame splits an FEVALUESBASE object in several objects. | |
This class is created for the objects handed to the local integration functions.
The objects of this class contain one or more objects of either type FEValues (cell) or FEFaceValues (face) or FESubfaceValues (subface) to be used in the local integration functions. They are stored in an array of pointers either to the base class FEValuesBase (cell) or to the base class FEFaceValuesBase (face, subface).
Therefore, the actual type FEVALUES
must be either FEValues (cell) or FEFaceValues (face) or FESubfaceValues (subface). The template parameter FEVALUESBASE
must be either FEValuesBase (cell) or FEFaceValuesBase (face, subface).
Additionally, this class includes containers to store the values (values), gradients (gradients), and hessians (hessians) of finite element functions, stored in global_data, in the quadrature points of mesh entities. You can use initialize_selector() function to select the FEVector vectors in global_data by names.
This class supports two local integration models. One is the standard model suggested by the use of FESystem. Namely, there is one FEValuseBase object in this class, containing all shape functions of the whole system, and having as many components as the system has. 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 local 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.
|
inline |
Constructor.
|
inline |
|
inline |
|
inline |
Access to a single actual FEVALUES object.
|
inline |
Access to a group of actual FEVALUES objects.
|
inline |
Use the finite element functions in global_data and fill the vectors values, gradients, hessians.
References MeshWorker::InfoObjects::fill_data().
|
inlineprivate |
Use the finite element functions in global_data and fill the vectors values, gradients, hessians.
References MeshWorker::InfoObjects::fill_data(), MeshWorker::selector_index(), and MeshWorker::selector_size().
|
inline |
(ZVN) Access to fe_val_unsplit
.
|
inline |
Build internal structures like fevalv and allocate memory for values, gradients, hessians.
FE_VALUES | is an object of the actual type FEVALUES which is either FEValues (cell) or FEFaceValues (face) or FESubfaceValues (subface). |
fe | is the finite element (for instance, FESystem object) stored in the DoFHandler. It is used in the constructor of #FE_VALUES. |
mapping | is the Mapping object used to map the mesh cells. It is used in the constructor of #FE_VALUES. |
quadrature | is a quadrature formula. It is used in the constructor of #FE_VALUES. |
flags | are the UpdateFlags. They are used in the constructor of #FE_VALUES. |
|
inline |
Initialize the global_data.
Cache the selector.
|
inline |
|
inline |
Reinitialize internal data for use on a cell.
References MeshWorker::InfoObjects::DoFInfo< dim, spacedim >::reinit().
|
inline |
Reinitialize internal data for use on a face.
References MeshWorker::InfoObjects::DoFInfo< dim, spacedim >::reinit().
|
inline |
Reinitialize internal data for use on a subface.
References MeshWorker::InfoObjects::DoFInfo< dim, spacedim >::reinit().
std::vector< std::vector< std::vector< Tensor<1,dim> > > >& MeshWorker::InfoObjects::IntegrationInfo< dim, FEVALUESBASE >::derivatives |
|
private |
(ZVN) In the mode # 1 (or in the main mode), AppFrame splits an FEVALUESBASE object in several objects.
Sometimes there is a need to keep the unsplit FEVALUESBASE object as well, even while in mode # 1.
Referenced by MeshWorker::InfoObjects::IntegrationInfo< dim, FEFaceValuesBase< dim > >::get_fe_val_unsplit().
|
private |
The vector of smart pointers to FEVALUESBASE objects.
Referenced by MeshWorker::InfoObjects::IntegrationInfo< dim, FEFaceValuesBase< dim > >::clear(), and MeshWorker::InfoObjects::IntegrationInfo< dim, FEFaceValuesBase< dim > >::fe().
SmartPointer<const FEVectors> MeshWorker::InfoObjects::IntegrationInfo< dim, FEVALUESBASE >::global_data |
The smart pointer to the FEVectors object called global_data.
It contains the named scalars and named vectors of an application.
Referenced by MeshWorker::InfoObjects::IntegrationInfo< dim, FEFaceValuesBase< dim > >::initialize_data(), and MeshWorker::InfoObjects::IntegrationInfo< dim, FEFaceValuesBase< dim > >::initialize_selector().
std::vector< std::vector< std::vector< Tensor<1,dim> > > > MeshWorker::InfoObjects::IntegrationInfo< dim, FEVALUESBASE >::gradients |
The vector containing the gradients of finite element functions in the quadrature points.
Example: nth solution, mth component, gradient at the qth quadrature point = gradients[n][m][q]
std::vector< std::vector< std::vector< Tensor<2,dim> > > > MeshWorker::InfoObjects::IntegrationInfo< dim, FEVALUESBASE >::hessians |
The vector containing the hessians of finite element functions in the quadrature points.
Example: nth solution, mth component, hessian at the qth quadrature point = hessians[n][m][q]
bool MeshWorker::InfoObjects::IntegrationInfo< dim, FEVALUESBASE >::multigrid |
This is true if we assemble for multigrid.
VectorSelector MeshWorker::InfoObjects::IntegrationInfo< dim, FEVALUESBASE >::selector |
This object selects the named FEVector vectors in global_data and assigns the types of information extracted.
These types are values, gradients, and hessians.
Referenced by MeshWorker::InfoObjects::IntegrationInfo< dim, FEFaceValuesBase< dim > >::initialize_data(), and MeshWorker::InfoObjects::IntegrationInfo< dim, FEFaceValuesBase< dim > >::initialize_selector().
std::vector< std::vector< std::vector<double> > > MeshWorker::InfoObjects::IntegrationInfo< dim, FEVALUESBASE >::values |
The vector containing the values of finite element functions in the quadrature points.
Example: nth solution, mth component, value at the qth quadrature point = values[n][m][q]