Reference documentation for deal.II version 9.0.0
|
#include <deal.II/hp/fe_values.h>
Additional Inherited Members | |
Protected Member Functions inherited from internal::hp::FEValuesBase< dim, dim-1,::FESubfaceValues< dim, spacedim > > | |
::FESubfaceValues< 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,::FESubfaceValues< dim, spacedim > > | |
const SmartPointer< const ::hp::FECollection< dim, ::FESubfaceValues< dim, spacedim > ::space_dimension >, FEValuesBase< dim, q_dim, ::FESubfaceValues< dim, spacedim > > > | fe_collection |
const SmartPointer< const ::hp::MappingCollection< dim, ::FESubfaceValues< dim, spacedim > ::space_dimension >, FEValuesBase< dim, q_dim, ::FESubfaceValues< dim, spacedim > > > | mapping_collection |
const ::hp::QCollection< q_dim > | q_collection |
This class implements for subfaces what hp::FEFaceValues does for faces. See there for further documentation.
Definition at line 486 of file fe_values.h.
FESubfaceValues< dim, spacedim >::FESubfaceValues | ( | 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 354 of file fe_values.cc.
FESubfaceValues< dim, spacedim >::FESubfaceValues | ( | 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 367 of file fe_values.cc.
void FESubfaceValues< dim, spacedim >::reinit | ( | const TriaIterator< DoFCellAccessor< DoFHandlerType, lda > > | cell, |
const unsigned int | face_no, | ||
const unsigned int | subface_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, face, and subface.
After the call, you can get an FESubfaceValues object using the get_present_fe_values() function that corresponds to the present cell. For this FESubfaceValues 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 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 380 of file fe_values.cc.
void FESubfaceValues< dim, spacedim >::reinit | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
const unsigned int | face_no, | ||
const unsigned int | subface_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 432 of file fe_values.cc.