FEValues< dim, spacedim > Class Template Reference
[Finite element access/FEValues classes]

Inheritance diagram for FEValues< dim, spacedim >:

Inheritance graph
[legend]

List of all members.

Public Member Functions

 FEValues (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
 FEValues (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
void reinit (const typename DoFHandler< dim, spacedim >::cell_iterator &cell)
void reinit (const typename hp::DoFHandler< dim, spacedim >::cell_iterator &cell)
void reinit (const typename MGDoFHandler< dim, spacedim >::cell_iterator &cell)
void reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell)
const Quadrature< dim > & get_quadrature () const
unsigned int memory_consumption () const
const FEValues< dim, spacedim > & get_present_fe_values () const

Static Public Attributes

static const unsigned int dimension = dim
static const unsigned int space_dimension = spacedim
static const unsigned int integral_dimension = dim

Private Member Functions

void initialize (const UpdateFlags update_flags)
void do_reinit ()

Private Attributes

const Quadrature< dim > quadrature


Detailed Description

template<int dim, int spacedim = dim>
class FEValues< dim, spacedim >

Finite element evaluated in quadrature points of a cell.

This function implements the initialization routines for FEValuesBase, if values in quadrature points of a cell are needed. For further documentation see this class.

Author:
Wolfgang Bangerth, 1998, Guido Kanschat, 2001

Constructor & Destructor Documentation

template<int dim, int spacedim = dim>
FEValues< dim, spacedim >::FEValues ( const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags  update_flags 
)

Constructor. Gets cell independent data from mapping and finite element objects, matching the quadrature rule and update flags.

template<int dim, int spacedim = dim>
FEValues< dim, spacedim >::FEValues ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags  update_flags 
)

Constructor. Uses MappingQ1 implicitly.


Member Function Documentation

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::reinit ( const typename DoFHandler< dim, spacedim >::cell_iterator &  cell  ) 

Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a DoFHandler object", and the finite element associated with this object. It is assumed that the finite element used by the given cell is also the one used by this FEValues object.

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::reinit ( const typename hp::DoFHandler< dim, spacedim >::cell_iterator &  cell  ) 

Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a hp::DoFHandler object", and the finite element associated with this object. It is assumed that the finite element used by the given cell is also the one used by this FEValues object.

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::reinit ( const typename MGDoFHandler< dim, spacedim >::cell_iterator &  cell  ) 

Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a MGDoFHandler object", and the finite element associated with this object. It is assumed that the finite element used by the given cell is also the one used by this FEValues object.

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::reinit ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell  ) 

Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a Triangulation object", and the given finite element. Since iterators into triangulation alone only convey information about the geometry of a cell, but not about degrees of freedom possibly associated with this cell, you will not be able to call some functions of this class if they need information about degrees of freedom. These functions are, above all, the get_function_value/gradients/hessians/laplacians functions. If you want to call these functions, you have to call the reinit variants that take iterators into DoFHandler or other DoF handler type objects.

template<int dim, int spacedim = dim>
const Quadrature<dim>& FEValues< dim, spacedim >::get_quadrature (  )  const

Return a reference to the copy of the quadrature formula stored by this object.

template<int dim, int spacedim = dim>
unsigned int FEValues< dim, spacedim >::memory_consumption (  )  const

Determine an estimate for the memory consumption (in bytes) of this object.

Reimplemented from FEValuesBase< dim, spacedim >.

template<int dim, int spacedim = dim>
const FEValues<dim,spacedim>& FEValues< dim, spacedim >::get_present_fe_values (  )  const

Return a reference to this very object.

Though it seems that it is not very useful, this function is there to provide capability to the hpFEValues class, in which case it provides the FEValues object for the present cell (remember that for hp finite elements, the actual FE object used may change from cell to cell, so we also need different FEValues objects for different cells; once you reinitialize the hpFEValues object for a specific cell, it retrieves the FEValues object for the FE on that cell and returns it through a function of the same name as this one; this function here therefore only provides the same interface so that one can templatize on FEValues/hpFEValues).

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::initialize ( const UpdateFlags  update_flags  )  [private]

Do work common to the two constructors.

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::do_reinit (  )  [private]

The reinit() functions do only that part of the work that requires knowledge of the type of iterator. After setting present_cell(), they pass on to this function, which does the real work, and which is independent of the actual type of the cell iterator.


Member Data Documentation

template<int dim, int spacedim = dim>
const unsigned int FEValues< dim, spacedim >::dimension = dim [static]

Dimension in which this object operates.

template<int dim, int spacedim = dim>
const unsigned int FEValues< dim, spacedim >::space_dimension = spacedim [static]

Dimension of the space in which this object operates.

template<int dim, int spacedim = dim>
const unsigned int FEValues< dim, spacedim >::integral_dimension = dim [static]

Dimension of the object over which we integrate. For the present class, this is equal to dim.

template<int dim, int spacedim = dim>
const Quadrature<dim> FEValues< dim, spacedim >::quadrature [private]

Store a copy of the quadrature formula here.


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

deal.II documentation generated on Sat Aug 15 16:51:56 2009 by doxygen 1.5.9