OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | Public Attributes | List of all members
AppFrame::LocalIntegrator< dim > Class Template Referenceabstract

The class using CellInfo and FaceInfo to compute local contributions to a 1form. More...

#include <dof_application.h>

Public Types

typedef
MeshWorker::InfoObjects::IntegrationInfo
< dim, FEValuesBase< dim, dim > > 
CellInfo
 
typedef
MeshWorker::InfoObjects::IntegrationInfo
< dim, FEFaceValuesBase< dim,
dim > > 
FaceInfo
 

Public Member Functions

 LocalIntegrator ()
 Virtual constructor.
 
virtual ~LocalIntegrator ()
 Virtual destructor.
 
virtual void cell (FEVector &cell_vector, const CellInfo &cell)=0
 Integration of local cell form.
 
virtual void bdry (FEVector &face_vector, const FaceInfo &face)=0
 Integration of local boundary form.
 
virtual void face (FEVector &face_vector1, FEVector &face_vector2, const FaceInfo &face1, const FaceInfo &face2)=0
 Integration of local interior face form.
 
virtual void post_loop (FEVector &dst, const FEVectors &srcs)
 An optional function called after the integration loop.
 
void set_gauss_quadrature (unsigned int n_cell_points, unsigned int n_bdry_points, unsigned int n_face_points)
 Assign n-point Gauss quadratures to each of the quadrature rules.
 

Public Attributes

Quadrature< dimcell_quadrature
 The quadrature rule used on cells.
 
Quadrature< dim-1 > bdry_quadrature
 The quadrature rule used on boundary faces.
 
Quadrature< dim-1 > face_quadrature
 The quadrature rule used on interior faces.
 

Detailed Description

template<int dim>
class AppFrame::LocalIntegrator< dim >

The class using CellInfo and FaceInfo to compute local contributions to a 1form.

Objects of this class are handed to the function DoFApplication::integrate_1form(). Derived classes for instance compute local contributions to the residual.

In addition to the local integration functions cell(), bdry() and face(), objects of this class provide control information to the DoFApplication loop:

Author
Guido Kanschat, 2007

Member Typedef Documentation

template<int dim>
typedef MeshWorker::InfoObjects::IntegrationInfo<dim, FEFaceValuesBase<dim, dim> > AppFrame::LocalIntegrator< dim >::FaceInfo

Constructor & Destructor Documentation

template<int dim>
AppFrame::LocalIntegrator< dim >::LocalIntegrator ( )

Virtual constructor.

template<int dim>
virtual AppFrame::LocalIntegrator< dim >::~LocalIntegrator ( )
virtual

Virtual destructor.

Member Function Documentation

template<int dim>
virtual void AppFrame::LocalIntegrator< dim >::bdry ( FEVector face_vector,
const FaceInfo face 
)
pure virtual

Integration of local boundary form.

template<int dim>
virtual void AppFrame::LocalIntegrator< dim >::cell ( FEVector cell_vector,
const CellInfo cell 
)
pure virtual

Integration of local cell form.

template<int dim>
virtual void AppFrame::LocalIntegrator< dim >::face ( FEVector face_vector1,
FEVector face_vector2,
const FaceInfo face1,
const FaceInfo face2 
)
pure virtual

Integration of local interior face form.

template<int dim>
virtual void AppFrame::LocalIntegrator< dim >::post_loop ( FEVector dst,
const FEVectors srcs 
)
virtual

An optional function called after the integration loop.

template<int dim>
void AppFrame::LocalIntegrator< dim >::set_gauss_quadrature ( unsigned int  n_cell_points,
unsigned int  n_bdry_points,
unsigned int  n_face_points 
)

Assign n-point Gauss quadratures to each of the quadrature rules.

Here, a size of zero points means that no loop over these grid entities should be performed.

Member Data Documentation

template<int dim>
Quadrature<dim-1> AppFrame::LocalIntegrator< dim >::bdry_quadrature

The quadrature rule used on boundary faces.

template<int dim>
Quadrature<dim> AppFrame::LocalIntegrator< dim >::cell_quadrature

The quadrature rule used on cells.

template<int dim>
Quadrature<dim-1> AppFrame::LocalIntegrator< dim >::face_quadrature

The quadrature rule used on interior faces.


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