Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | 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 MappingCollection< dim, spacedim > &mapping_collection, const FECollection< dim, spacedim > &fe_collection, const QCollection< dim > &q_collection, const UpdateFlags update_flags)
 
 FEValues (const FECollection< dim, spacedim > &fe_collection, const QCollection< dim > &q_collection, const UpdateFlags update_flags)
 
template<bool lda>
void reinit (const TriaIterator< DoFCellAccessor< dim, spacedim, 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)
 
void precalculate_fe_values (const std::vector< unsigned int > &fe_indices, const std::vector< unsigned int > &mapping_indices, const std::vector< unsigned int > &q_indices)
 
void precalculate_fe_values ()
 
const FECollection< dim, FEValuesType::space_dimension > & get_fe_collection () const
 
const MappingCollection< dim, FEValuesType::space_dimension > & get_mapping_collection () const
 
const QCollection< q_dim > & get_quadrature_collection () const
 
UpdateFlags get_update_flags () const
 
const ::FEValues< dim, dim > & get_present_fe_values () const
 

Static Public Attributes

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

Protected Member Functions

::FEValues< dim, dim > & select_fe_values (const unsigned int fe_index, const unsigned int mapping_index, const unsigned int q_index)
 

Protected Attributes

const SmartPointer< const FECollection< dim, FEValuesType::space_dimension >, FEValuesBase< dim, q_dim, ::FEValues< dim, dim > > > fe_collection
 
const SmartPointer< const MappingCollection< dim, FEValuesType::space_dimension >, FEValuesBase< dim, q_dim, ::FEValues< dim, dim > > > mapping_collection
 
const QCollection< q_dim > q_collection
 
const std::vector< QCollection< q_dim > > q_collections
 

Private Attributes

Table< 3, std::unique_ptr< ::FEValues< dim, dim > > > fe_values_table
 
TableIndices< 3 > present_fe_values_index
 
const UpdateFlags update_flags
 

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-DoFHandler 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. Alternatively if desired, this can be bypassed by computing all objects in advance with the corresponding hp::FEValuesBase::precalculate_fe_values() function.

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

Definition at line 315 of file fe_values.h.

Constructor & Destructor Documentation

◆ FEValues() [1/2]

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

Constructor. Initialize this object with the given parameters.

Definition at line 267 of file fe_values.cc.

◆ FEValues() [2/2]

template<int dim, int spacedim>
FEValues< dim, spacedim >::FEValues ( const FECollection< dim, spacedim > &  fe_collection,
const 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.

Definition at line 280 of file fe_values.cc.

Member Function Documentation

◆ reinit() [1/2]

template<int dim, int spacedim>
template<bool lda>
void FEValues< dim, spacedim >::reinit ( const TriaIterator< DoFCellAccessor< dim, spacedim, 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 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 293 of file fe_values.cc.

◆ reinit() [2/2]

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-DoFHandler iterators. The reason this function exists is so that one can use hp::FEValues for Triangulation objects too.

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 336 of file fe_values.cc.

◆ precalculate_fe_values() [1/2]

void FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::precalculate_fe_values ( const std::vector< unsigned int > &  fe_indices,
const std::vector< unsigned int > &  mapping_indices,
const std::vector< unsigned int > &  q_indices 
)
inherited

For timing purposes it may be useful to create all required FE*Values objects in advance, rather than computing them on request via lazy allocation as usual in this class.

This function precalculates the FE*Values objects corresponding to the provided parameters: The total of all vector entries corresponding to the same index describes an FE*Values object similarly to select_fe_values().

Definition at line 137 of file fe_values.cc.

◆ precalculate_fe_values() [2/2]

void FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::precalculate_fe_values
inherited

Same as above, geared to the most common use of hp::FEValues objects in which FE, quadrature and mapping indices are similar on each individual cell.

FE*Values objects are created for every FE in the FECollection, with quadrature and mapping corresponding to the same index from the QuadratureCollection and MappingCollection, respectively.

If QuadratureCollection or MappingCollection contains only one object, it is used for all FE*Values objects.

Definition at line 154 of file fe_values.cc.

◆ get_fe_collection()

const FECollection< dim, FEValuesType::space_dimension > & FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::get_fe_collection
inlineinherited

Get a reference to the collection of finite element objects used here.

Definition at line 161 of file fe_values.h.

◆ get_mapping_collection()

const MappingCollection< dim, FEValuesType::space_dimension > & FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::get_mapping_collection
inlineinherited

Get a reference to the collection of mapping objects used here.

Definition at line 167 of file fe_values.h.

◆ get_quadrature_collection()

const QCollection< q_dim > & FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::get_quadrature_collection
inlineinherited

Get a reference to the collection of quadrature objects used here.

Definition at line 173 of file fe_values.h.

◆ get_update_flags()

UpdateFlags FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::get_update_flags
inlineinherited

Get the underlying update flags.

Definition at line 179 of file fe_values.h.

◆ get_present_fe_values()

const ::FEValues< dim, dim > & FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::get_present_fe_values
inlineinherited

Return a reference to the FEValues object selected by the last call to select_fe_values(). select_fe_values() in turn is called when you called the reinit function of the hp::FE*Values class the last time.

Definition at line 188 of file fe_values.h.

◆ select_fe_values()

::FEValues< dim, dim > & FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::select_fe_values ( const unsigned int  fe_index,
const unsigned int  mapping_index,
const unsigned int  q_index 
)
protectedinherited

Select a FEValues object suitable for the given FE, quadrature, and mapping indices. If such an object doesn't yet exist, create one.

The function returns a writable reference so that derived classes can also reinit() the selected FEValues object.

Definition at line 199 of file fe_values.cc.

Member Data Documentation

◆ dimension

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

Definition at line 319 of file fe_values.h.

◆ space_dimension

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

Definition at line 321 of file fe_values.h.

◆ fe_collection

const SmartPointer<const FECollection<dim, FEValuesType::space_dimension>, FEValuesBase<dim, q_dim, ::FEValues< dim, dim > > > hp::FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::fe_collection
protectedinherited

A pointer to the collection of finite elements to be used.

Definition at line 209 of file fe_values.h.

◆ mapping_collection

const SmartPointer< const MappingCollection<dim, FEValuesType::space_dimension>, FEValuesBase<dim, q_dim, ::FEValues< dim, dim > > > hp::FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::mapping_collection
protectedinherited

A pointer to the collection of mappings to be used.

Definition at line 217 of file fe_values.h.

◆ q_collection

const QCollection<q_dim> hp::FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::q_collection
protectedinherited

Copy of the quadrature collection object provided to the constructor.

Definition at line 222 of file fe_values.h.

◆ q_collections

const std::vector<QCollection<q_dim> > hp::FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::q_collections
protectedinherited

Vector of quadrature collections. For hp::FEFaceValues, the ith entry of the quadrature collections are interpreted as the face quadrature rules to be applied the ith face.

The variable q_collection collects the first quadrature rule of each quadrature collection of the vector.

Definition at line 232 of file fe_values.h.

◆ fe_values_table

Table<3, std::unique_ptr<::FEValues< dim, dim > > > hp::FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::fe_values_table
privateinherited

A table in which we store pointers to fe_values objects for different finite element, mapping, and quadrature objects from our collection. The first index indicates the index of the finite element within the fe_collection, the second the index of the mapping within the mapping collection, and the last one the index of the quadrature formula within the q_collection.

Initially, all entries have zero pointers, and we will allocate them lazily as needed in select_fe_values() or precalculate_fe_values().

Definition at line 246 of file fe_values.h.

◆ present_fe_values_index

TableIndices<3> hp::FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::present_fe_values_index
privateinherited

Set of indices pointing at the fe_values object selected last time the select_fe_value() function was called.

Definition at line 252 of file fe_values.h.

◆ update_flags

const UpdateFlags hp::FEValuesBase< dim, q_dim, ::FEValues< dim, dim > >::update_flags
privateinherited

Values of the update flags as given to the constructor.

Definition at line 257 of file fe_values.h.


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