Reference documentation for deal.II version 9.1.0-pre
Public Member Functions | List of all members
hp::FEValues< dim, spacedim > Class Template Reference

#include <deal.II/hp/fe_values.h>

Inheritance diagram for hp::FEValues< dim, spacedim >:
[legend]

Public Member Functions

 FEValues (const ::hp::MappingCollection< dim, spacedim > &mapping_collection, const ::hp::FECollection< dim, spacedim > &fe_collection, const ::hp::QCollection< dim > &q_collection, const UpdateFlags update_flags)
 
 FEValues (const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim > &q_collection, const UpdateFlags update_flags)
 
template<typename DoFHandlerType , bool lda>
void reinit (const TriaIterator< DoFCellAccessor< DoFHandlerType, lda >> cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
 
void reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
 
- Public Member Functions inherited from internal::hp::FEValuesBase< dim, dim,::FEValues< dim, spacedim > >
 FEValuesBase (const ::hp::MappingCollection< dim,::FEValues< dim, spacedim >::space_dimension > &mapping_collection, const ::hp::FECollection< dim,::FEValues< dim, spacedim >::space_dimension > &fe_collection, const ::hp::QCollection< q_dim > &q_collection, const ::UpdateFlags update_flags)
 
 FEValuesBase (const ::hp::FECollection< dim,::FEValues< dim, spacedim >::space_dimension > &fe_collection, const ::hp::QCollection< q_dim > &q_collection, const UpdateFlags update_flags)
 
const ::hp::FECollection< dim,::FEValues< dim, spacedim >::space_dimension > & get_fe_collection () const
 
const ::hp::MappingCollection< dim,::FEValues< dim, spacedim >::space_dimension > & get_mapping_collection () const
 
const ::hp::QCollection< q_dim > & get_quadrature_collection () const
 
UpdateFlags get_update_flags () const
 
const ::FEValues< dim, spacedim > & get_present_fe_values () const
 

Additional Inherited Members

- Protected Member Functions inherited from internal::hp::FEValuesBase< dim, dim,::FEValues< dim, spacedim > >
::FEValues< dim, spacedim > & select_fe_values (const unsigned int fe_index, const unsigned int mapping_index, const unsigned int q_index)
 
- Protected Attributes inherited from internal::hp::FEValuesBase< dim, dim,::FEValues< dim, spacedim > >
const SmartPointer< const ::hp::FECollection< dim,::FEValues< dim, spacedim >::space_dimension >, FEValuesBase< dim, q_dim,::FEValues< dim, spacedim > > > fe_collection
 
const SmartPointer< const ::hp::MappingCollection< dim,::FEValues< dim, spacedim >::space_dimension >, FEValuesBase< dim, q_dim,::FEValues< dim, spacedim > > > mapping_collection
 
const ::hp::QCollection< q_dim > q_collection
 

Detailed Description

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

An hp equivalent of the FEValues class. See the step-27 tutorial program for examples of use.

The idea of this class is as follows: when one assembled matrices in the hp finite element method, there may be different finite elements on different cells, and consequently one may also want to use different quadrature formulas for different cells. On the other hand, the FEValues efficiently handles pre-evaluating whatever information is necessary for a single finite element and quadrature object. This class brings these concepts together: it provides a "collection" of FEValues objects.

Upon construction, one passes not one finite element and quadrature object (and possible a mapping), but a whole collection of type hp::FECollection and hp::QCollection. Later on, when one sits on a concrete cell, one would call the reinit() function for this particular cell, just as one does for a regular FEValues object. The difference is that this time, the reinit() function looks up the active_fe_index of that cell, if necessary creates a FEValues object that matches the finite element and quadrature formulas with that particular index in their collections, and then re-initializes it for the current cell. The FEValues object that then fits the finite element and quadrature formula for the current cell can then be accessed using the get_present_fe_values() function, and one would work with it just like with any FEValues object for non-hp DoF handler objects.

The reinit() functions have additional arguments with default values. If not specified, the function takes the index into the hp::FECollection, hp::QCollection, and hp::MappingCollection objects from the active_fe_index of the cell, as explained above. However, one can also select different indices for a current cell. For example, by specifying a different index into the hp::QCollection class, one does not need to sort the quadrature objects in the quadrature collection so that they match one-to-one the order of finite element objects in the FE collection (even though choosing such an order is certainly convenient).

Note that FEValues objects are created on the fly, i.e. only as they are needed. This ensures that we do not create objects for every combination of finite element, quadrature formula and mapping, but only those that will actually be needed.

This class has not yet been implemented for the use in the codimension one case (spacedim != dim ).

Author
Wolfgang Bangerth, 2003

Definition at line 239 of file fe_values.h.

Constructor & Destructor Documentation

template<int dim, int spacedim = dim>
hp::FEValues< dim, spacedim >::FEValues ( const ::hp::MappingCollection< dim, spacedim > &  mapping_collection,
const ::hp::FECollection< dim, spacedim > &  fe_collection,
const ::hp::QCollection< dim > &  q_collection,
const UpdateFlags  update_flags 
)

Constructor. Initialize this object with the given parameters.

The finite element collection parameter is actually ignored, but is in the signature of this function to make it compatible with the signature of the respective constructor of the usual FEValues object, with the respective parameter in that function also being the return value of the DoFHandler::get_fe() function.

template<int dim, int spacedim>
FEValues< dim, spacedim >::FEValues ( const hp::FECollection< dim, spacedim > &  fe_collection,
const hp::QCollection< dim > &  q_collection,
const UpdateFlags  update_flags 
)

Constructor. This constructor is equivalent to the other one except that it makes the object use a \(Q_1\) mapping (i.e., an object of type MappingQGeneric(1)) implicitly.

The finite element collection parameter is actually ignored, but is in the signature of this function to make it compatible with the signature of the respective constructor of the usual FEValues object, with the respective parameter in that function also being the return value of the DoFHandler::get_fe() function.

Definition at line 128 of file fe_values.cc.

Member Function Documentation

template<int dim, int spacedim>
template<typename DoFHandlerType , bool lda>
void FEValues< dim, spacedim >::reinit ( const TriaIterator< DoFCellAccessor< DoFHandlerType, lda >>  cell,
const unsigned int  q_index = numbers::invalid_unsigned_int,
const unsigned int  mapping_index = numbers::invalid_unsigned_int,
const unsigned int  fe_index = numbers::invalid_unsigned_int 
)

Reinitialize the object for the given cell.

After the call, you can get an FEValues object using the get_present_fe_values() function that corresponds to the present cell. For this FEValues object, we use the additional arguments described below to determine which finite element, mapping, and quadrature formula to use. They are order in such a way that the arguments one may want to change most frequently come first. The rules for these arguments are as follows:

If the fe_index argument to this function is left at its default value, then we use that finite element within the hp::FECollection passed to the constructor of this class with index given by cell->active_fe_index(). Consequently, the hp::FECollection argument given to this object should really be the same as that used in the construction of the hp::DoFHandler associated with the present cell. On the other hand, if a value is given for this argument, it overrides the choice of cell->active_fe_index().

If the q_index argument is left at its default value, then we use that quadrature formula within the hp::QCollection passed to the constructor of this class with index given by cell->active_fe_index(), i.e. the same index as that of the finite element. In this case, there should be a corresponding quadrature formula for each finite element in the hp::FECollection. As a special case, if the quadrature collection contains only a single element (a frequent case if one wants to use the same quadrature object for all finite elements in an hp discretization, even if that may not be the most efficient), then this single quadrature is used unless a different value for this argument is specified. On the other hand, if a value is given for this argument, it overrides the choice of cell->active_fe_index() or the choice for the single quadrature.

If the mapping_index argument is left at its default value, then we use that mapping object within the hp::MappingCollection passed to the constructor of this class with index given by cell->active_fe_index(), i.e. the same index as that of the finite element. As above, if the mapping collection contains only a single element (a frequent case if one wants to use a \(Q_1\) mapping for all finite elements in an hp discretization), then this single mapping is used unless a different value for this argument is specified.

Definition at line 142 of file fe_values.cc.

template<int dim, int spacedim>
void FEValues< dim, spacedim >::reinit ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  q_index = numbers::invalid_unsigned_int,
const unsigned int  mapping_index = numbers::invalid_unsigned_int,
const unsigned int  fe_index = numbers::invalid_unsigned_int 
)

Like the previous function, but for non-hp iterators. The reason this (and the other non-hp iterator) function exists is so that one can use hp::FEValues not only for hp::DoFhandler objects, but for all sorts of DoFHandler objects, and triangulations not associated with DoFHandlers in general.

Since cell->active_fe_index() doesn't make sense for triangulation iterators, this function chooses the zero-th finite element, mapping, and quadrature object from the relevant constructions passed to the constructor of this object. The only exception is if you specify a value different from the default value for any of these last three arguments.

Definition at line 193 of file fe_values.cc.


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