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

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

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

Public Member Functions

 FEFaceValues (const hp::MappingCollection< dim, spacedim > &mapping_collection, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim-1 > &q_collection, const UpdateFlags update_flags)
 
 FEFaceValues (const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim-1 > &q_collection, const UpdateFlags update_flags)
 
template<typename DoFHandlerType , bool lda>
void reinit (const TriaIterator< DoFCellAccessor< DoFHandlerType, lda > > cell, const unsigned int face_no, 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 face_no, 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-1,::FEFaceValues< dim, spacedim > >
 FEValuesBase (const ::hp::MappingCollection< dim, ::FEFaceValues< dim, spacedim > ::space_dimension > &mapping_collection, const ::hp::FECollection< dim, ::FEFaceValues< dim, spacedim > ::space_dimension > &fe_collection, const ::hp::QCollection< q_dim > &q_collection, const ::UpdateFlags update_flags)
 
 FEValuesBase (const ::hp::FECollection< dim, ::FEFaceValues< dim, spacedim > ::space_dimension > &fe_collection, const ::hp::QCollection< q_dim > &q_collection, const UpdateFlags update_flags)
 
const ::hp::FECollection< dim, ::FEFaceValues< dim, spacedim > ::space_dimension > & get_fe_collection () const
 
const ::hp::MappingCollection< dim, ::FEFaceValues< dim, spacedim > ::space_dimension > & get_mapping_collection () const
 
const ::hp::QCollection< q_dim > & get_quadrature_collection () const
 
UpdateFlags get_update_flags () const
 
const ::FEFaceValues< dim, spacedim > & get_present_fe_values () const
 

Additional Inherited Members

- Protected Member Functions inherited from internal::hp::FEValuesBase< dim, dim-1,::FEFaceValues< dim, spacedim > >
::FEFaceValues< 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-1,::FEFaceValues< dim, spacedim > >
const SmartPointer< const ::hp::FECollection< dim, ::FEFaceValues< dim, spacedim > ::space_dimension >, FEValuesBase< dim, q_dim, ::FEFaceValues< dim, spacedim > > > fe_collection
 
const SmartPointer< const ::hp::MappingCollection< dim, ::FEFaceValues< dim, spacedim > ::space_dimension >, FEValuesBase< dim, q_dim, ::FEFaceValues< dim, spacedim > > > mapping_collection
 
const ::hp::QCollection< q_dim > q_collection
 

Detailed Description

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

This is the equivalent of the hp::FEValues class but for face integrations, i.e. it is to hp::FEValues what FEFaceValues is to FEValues.

The same comments apply as in the documentation of the hp::FEValues class. However, it is important to note that it is here more common that one would want to explicitly specify an index to a particular quadrature formula in the reinit() functions. This is because the default index corresponds to the finite element index on the current function. On the other hand, integration on faces will typically have to happen with a quadrature formula that is adjusted to the finite elements used on both sides of a face. If one sorts the elements of the hp::FECollection with ascending polynomial degree, and matches these finite elements with corresponding quadrature formulas in the hp::QCollection passed to the constructor, then the quadrature index passed to the reinit() function should typically be something like std::max (cell->active_fe_index(), neighbor->active_fe_index() to ensure that a quadrature formula is chosen that is sufficiently accurate for both finite elements.

Author
Wolfgang Bangerth, 2003

Definition at line 368 of file fe_values.h.

Constructor & Destructor Documentation

◆ FEFaceValues() [1/2]

template<int dim, int spacedim>
FEFaceValues< dim, spacedim >::FEFaceValues ( const hp::MappingCollection< dim, spacedim > &  mapping_collection,
const hp::FECollection< dim, spacedim > &  fe_collection,
const hp::QCollection< dim-1 > &  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.

Definition at line 235 of file fe_values.cc.

◆ FEFaceValues() [2/2]

template<int dim, int spacedim>
FEFaceValues< dim, spacedim >::FEFaceValues ( const hp::FECollection< dim, spacedim > &  fe_collection,
const hp::QCollection< dim-1 > &  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 248 of file fe_values.cc.

Member Function Documentation

◆ reinit() [1/2]

template<int dim, int spacedim>
template<typename DoFHandlerType , bool lda>
void FEFaceValues< dim, spacedim >::reinit ( const TriaIterator< DoFCellAccessor< DoFHandlerType, lda > >  cell,
const unsigned int  face_no,
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 and face.

After the call, you can get an FEFaceValues object using the get_present_fe_values() function that corresponds to the present cell. For this FEFaceValues 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 261 of file fe_values.cc.

◆ reinit() [2/2]

template<int dim, int spacedim>
void FEFaceValues< dim, spacedim >::reinit ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
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 312 of file fe_values.cc.


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