OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE > Class Template Reference

This class is created for the objects handed to the mesh loops. More...

#include <mesh_loop_info_objects.h>

Inheritance diagram for FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >:
Inheritance graph
[legend]
Collaboration diagram for FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >:
Collaboration graph
[legend]

Public Member Functions

Constructors, destructor, and initialization
 IntegrationInfo (const BlockInfo &block_info)
 Constructor. More...
 
 IntegrationInfo (const FEVectors &data, const BlockInfo &block_info)
 Constructor. More...
 
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 fevalv and allocate memory for values, gradients, hessians. More...
 
void initialize_data (const FEVectors &data)
 Initialize global_data. More...
 
Reinitialization
template<typename DHCellIterator >
void reinit (const DHCellIterator &c)
 Reinitialize internal data for use on a cell. More...
 
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 of a cell. More...
 
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 of a face of a cell. More...
 
template<typename TYPE >
void fill_local_data (std::vector< std::vector< std::vector< TYPE > > > &data, bool split_fevalues) const
 Use the finite element functions in global_data and fill the vectors values, gradients, hessians. More...
 
Accessors and info
const FEVALUESBASE & fe () const
 Access to a single actual FEVALUES object. More...
 
const FEVALUESBASE & fe (unsigned int i) const
 Access to a group of actual FEVALUES objects. More...
 
const FEVALUESBASE & get_fe_val_unsplit () const
 Access to fe_val_unsplit. More...
 
void clear ()
 Resize fevalv to 0. More...
 
- Public Member Functions inherited from FuelCell::ApplicationCore::DoFInfo< dim >
 DoFInfo (const BlockInfo &block_info)
 Constructor. More...
 
 DoFInfo (const FEVectors &, const BlockInfo &block_info)
 
void reinit (const DHCellIterator &c)
 Set the current cell and fill indices. More...
 
void reinit (const DHCellIterator &c, const DHFaceIterator &f, const unsigned int fn)
 Set the current cell and face and fill indices. More...
 
void reinit (const DHCellIterator &c, const DHFaceIterator &f, const unsigned int fn, const unsigned int sn)
 Set the current cell, face, and subface and fill indices. More...
 

DATA

bool multigrid
 true if we assemble for multigrid. More...
 
SmartPointer< const FEVectorsglobal_data
 The smart pointer to the FEVectors object called global_data. More...
 
std::vector< std::vector
< std::vector< double > > > 
values
 The vector containing the values of finite element functions in the quadrature points. More...
 
std::vector< std::vector
< std::vector< Tensor< 1, dim > > > > 
gradients
 The vector containing the gradients of finite element functions in the quadrature points. More...
 
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. More...
 
std::vector< boost::shared_ptr
< FEVALUESBASE > > 
fevalv
 The vector of smart pointers to FEVALUESBASE objects. More...
 
boost::shared_ptr< FEVALUESBASE > fe_val_unsplit
 In the MODE1 (or in the main mode), the class splits an FEVALUESBASE object into several sub-objects and fills fevalv. More...
 

Additional Inherited Members

- Public Attributes inherited from FuelCell::ApplicationCore::DoFInfo< dim >
DoFHandler< dim >::cell_iterator dof_cell
 The current DoFHandler<dim> cell. More...
 
DoFHandler< dim >
::active_cell_iterator 
dof_active_cell
 The current DoFHandler<dim> active cell. More...
 
DoFHandler< dim >::face_iterator dof_face
 The current DoFHandler<dim> face. More...
 
Triangulation< dim >::cell_iterator cell
 The current Triangulation<dim> cell. More...
 
Triangulation< dim >::face_iterator face
 The current Triangulation<dim> face. More...
 
unsigned int face_number
 The number of the current face on the current cell. More...
 
unsigned int sub_number
 The number of the current subface on the current face of the current cell. More...
 
std::vector< unsigned int > indices
 The local dof indices associated with the current cell. More...
 
SmartPointer< const BlockInfoblock_info
 The block structure of the problem at hand. More...
 

Detailed Description

template<int dim, typename FEVALUESBASE>
class FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >

This class is created for the objects handed to the mesh loops.

The objects of this class contain an object of type FEVALUESBASE which is either FEValuesBase (base for cell) or FEFaceValuesBase (base for face and subface).

At the same time, the actual type FEVALUES for FEValuesBase is FEValues (cell) and the actual type FEVALUES for FEFaceValuesBase is either FEFaceValues (face) or FESubfaceValues (subface).

Note
FEVALUESBASE is a template parameter of the class itself.
FEVALUES is a template parameter of the initialize() function.

Additionally, this class includes containers to store the values, gradients and hessians of finite element functions stored in global_data in the quadrature points of a mesh entity.

Author
Guido Kanschat
Valentin N. Zingan

Constructor & Destructor Documentation

template<int dim, typename FEVALUESBASE >
FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::IntegrationInfo ( const BlockInfo block_info)
inline

Constructor.

template<int dim, typename FEVALUESBASE >
FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::IntegrationInfo ( const FEVectors data,
const BlockInfo block_info 
)
inline

Constructor.

Member Function Documentation

template<int dim, typename FEVALUESBASE >
void FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::clear ( )
inline

Resize fevalv to 0.

template<int dim, typename FEVALUESBASE >
const FEVALUESBASE & FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::fe ( ) const
inline

Access to a single actual FEVALUES object.

Referenced by FuelCellShop::PostProcessing::WaterSorptionResponse< dim >::compute_responses().

Here is the caller graph for this function:

template<int dim, typename FEVALUESBASE >
const FEVALUESBASE & FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::fe ( unsigned int  i) const
inline

Access to a group of actual FEVALUES objects.

template<int dim, typename FEVALUESBASE >
template<typename TYPE >
void FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::fill_local_data ( std::vector< std::vector< std::vector< TYPE > > > &  data,
bool  split_fevalues 
) const
inline

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

References FuelCell::ApplicationCore::fill_data().

Here is the call graph for this function:

template<int dim, typename FEVALUESBASE >
const FEVALUESBASE & FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::get_fe_val_unsplit ( ) const
inline

Access to fe_val_unsplit.

template<int dim, typename FEVALUESBASE >
template<typename FEVALUES >
void FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::initialize ( const FEVALUES *  FE_VALUES,
const FiniteElement< dim > &  fe,
const Mapping< dim > &  mapping,
const Quadrature< FEVALUES::integral_dimension > &  quadrature,
const UpdateFlags  flags 
)
inline

Build internal structures fevalv and allocate memory for values, gradients, hessians.

Parameters
FE_VALUESis an object of the actual type FEVALUES which is either FEValues (cell) or FEFaceValues (face) or FESubfaceValues (subface).
Note
FE_VALUES is really not needed here and only serves to define the actual type FEVALUES.
Parameters
feis the finite element (for instance, FESystem object) stored in the DoFHandler. It is used in the constructor of FE_VALUES.
mappingis the Mapping object used to map the mesh cells. It is used in the constructor of FE_VALUES.
quadratureis a quadrature formula. It is used in the constructor of FE_VALUES.
flagsare the UpdateFlags. They are used in the constructor of FE_VALUES.
template<int dim, typename FEVALUESBASE >
void FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::initialize_data ( const FEVectors data)
inline

Initialize global_data.

template<int dim, typename FEVALUESBASE >
template<typename DHCellIterator >
void FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::reinit ( const DHCellIterator &  c)
inline

Reinitialize internal data for use on a cell.

References FuelCell::ApplicationCore::DoFInfo< dim, spacedim >::reinit().

Here is the call graph for this function:

template<int dim, typename FEVALUESBASE >
template<typename DHCellIterator , typename DHFaceIterator >
void FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::reinit ( const DHCellIterator &  c,
const DHFaceIterator &  f,
const unsigned int  fn 
)
inline

Reinitialize internal data for use on a face of a cell.

References FuelCell::ApplicationCore::DoFInfo< dim, spacedim >::reinit().

Here is the call graph for this function:

template<int dim, typename FEVALUESBASE >
template<typename DHCellIterator , typename DHFaceIterator >
void FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::reinit ( const DHCellIterator &  c,
const DHFaceIterator &  f,
const unsigned int  fn,
const unsigned int  sn 
)
inline

Reinitialize internal data for use on a subface of a face of a cell.

References FuelCell::ApplicationCore::DoFInfo< dim, spacedim >::reinit().

Here is the call graph for this function:

Member Data Documentation

template<int dim, typename FEVALUESBASE >
std::vector< std::vector< std::vector< Tensor<1,dim> > > >& FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::derivatives
Deprecated:
Reference to gradients for compatibility.
template<int dim, typename FEVALUESBASE >
boost::shared_ptr<FEVALUESBASE> FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::fe_val_unsplit
private

In the MODE1 (or in the main mode), the class splits an FEVALUESBASE object into several sub-objects and fills fevalv.

However, for the variety of problems, there is still a need to keep the unsplit FEVALUESBASE object as well.

template<int dim, typename FEVALUESBASE >
std::vector< boost::shared_ptr<FEVALUESBASE> > FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::fevalv
private

The vector of smart pointers to FEVALUESBASE objects.

template<int dim, typename FEVALUESBASE >
SmartPointer<const FEVectors> FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::global_data

The smart pointer to the FEVectors object called global_data.

Referenced by FuelCellShop::PostProcessing::WaterSorptionResponse< dim >::compute_responses().

template<int dim, typename FEVALUESBASE >
std::vector< std::vector< std::vector< Tensor<1,dim> > > > FuelCell::ApplicationCore::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]

template<int dim, typename FEVALUESBASE >
std::vector< std::vector< std::vector< Tensor<2,dim> > > > FuelCell::ApplicationCore::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]

template<int dim, typename FEVALUESBASE >
bool FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::multigrid

true if we assemble for multigrid.

template<int dim, typename FEVALUESBASE >
std::vector< std::vector< std::vector<double> > > FuelCell::ApplicationCore::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]

Referenced by FuelCellShop::PostProcessing::WaterSorptionResponse< dim >::compute_responses().


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