Reference documentation for deal.II version 8.1.0
Classes | Public Types | Public Member Functions | Private Attributes | List of all members

#include <fe_values.h>

Classes

struct  ShapeFunctionData
 

Public Types

typedef double value_type
 
typedef ::Tensor< 1, spacedim > gradient_type
 
typedef ::Tensor< 2, spacedim > hessian_type
 

Public Member Functions

 Scalar ()
 
 Scalar (const FEValuesBase< dim, spacedim > &fe_values_base, const unsigned int component)
 
Scalaroperator= (const Scalar< dim, spacedim > &)
 
value_type value (const unsigned int shape_function, const unsigned int q_point) const
 
gradient_type gradient (const unsigned int shape_function, const unsigned int q_point) const
 
hessian_type hessian (const unsigned int shape_function, const unsigned int q_point) const
 
template<class InputVector >
void get_function_values (const InputVector &fe_function, std::vector< value_type > &values) const
 
template<class InputVector >
void get_function_gradients (const InputVector &fe_function, std::vector< gradient_type > &gradients) const
 
template<class InputVector >
void get_function_hessians (const InputVector &fe_function, std::vector< hessian_type > &hessians) const
 
template<class InputVector >
void get_function_laplacians (const InputVector &fe_function, std::vector< value_type > &laplacians) const
 

Private Attributes

const FEValuesBase< dim, spacedim > & fe_values
 
const unsigned int component
 
std::vector< ShapeFunctionDatashape_function_data
 

Detailed Description

template<int dim, int spacedim = dim>
class FEValuesViews::Scalar< dim, spacedim >

A class representing a view to a single scalar component of a possibly vector-valued finite element. Views are discussed in the Handling vector valued problems module.

You get an object of this type if you apply a FEValuesExtractors::Scalar to an FEValues, FEFaceValues or FESubfaceValues object.

Definition at line 142 of file fe_values.h.

Member Typedef Documentation

template<int dim, int spacedim = dim>
typedef double FEValuesViews::Scalar< dim, spacedim >::value_type

A typedef for the data type of values of the view this class represents. Since we deal with a single components, the value type is a scalar double.

Definition at line 150 of file fe_values.h.

template<int dim, int spacedim = dim>
typedef ::Tensor<1,spacedim> FEValuesViews::Scalar< dim, spacedim >::gradient_type

A typedef for the type of gradients of the view this class represents. Here, for a scalar component of the finite element, the gradient is a Tensor<1,dim>.

Definition at line 157 of file fe_values.h.

template<int dim, int spacedim = dim>
typedef ::Tensor<2,spacedim> FEValuesViews::Scalar< dim, spacedim >::hessian_type

A typedef for the type of second derivatives of the view this class represents. Here, for a scalar component of the finite element, the Hessian is a Tensor<2,dim>.

Definition at line 164 of file fe_values.h.

Constructor & Destructor Documentation

template<int dim, int spacedim = dim>
FEValuesViews::Scalar< dim, spacedim >::Scalar ( )

Default constructor. Creates an invalid object.

template<int dim, int spacedim = dim>
FEValuesViews::Scalar< dim, spacedim >::Scalar ( const FEValuesBase< dim, spacedim > &  fe_values_base,
const unsigned int  component 
)

Constructor for an object that represents a single scalar component of a FEValuesBase object (or of one of the classes derived from FEValuesBase).

Member Function Documentation

template<int dim, int spacedim = dim>
Scalar& FEValuesViews::Scalar< dim, spacedim >::operator= ( const Scalar< dim, spacedim > &  )

Copy operator. This is not a lightweight object so we don't allow copying and generate an exception if this function is called.

template<int dim, int spacedim = dim>
value_type FEValuesViews::Scalar< dim, spacedim >::value ( const unsigned int  shape_function,
const unsigned int  q_point 
) const

Return the value of the vector component selected by this view, for the shape function and quadrature point selected by the arguments.

Parameters
shape_functionNumber of the shape function to be evaluated. Note that this number runs from zero to dofs_per_cell, even in the case of an FEFaceValues or FESubfaceValues object.
q_pointNumber of the quadrature point at which function is to be evaluated
template<int dim, int spacedim = dim>
gradient_type FEValuesViews::Scalar< dim, spacedim >::gradient ( const unsigned int  shape_function,
const unsigned int  q_point 
) const

Return the gradient (a tensor of rank 1) of the vector component selected by this view, for the shape function and quadrature point selected by the arguments.

Note
The meaning of the arguments is as documented for the value() function.
template<int dim, int spacedim = dim>
hessian_type FEValuesViews::Scalar< dim, spacedim >::hessian ( const unsigned int  shape_function,
const unsigned int  q_point 
) const

Return the Hessian (the tensor of rank 2 of all second derivatives) of the vector component selected by this view, for the shape function and quadrature point selected by the arguments.

Note
The meaning of the arguments is as documented for the value() function.
template<int dim, int spacedim = dim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
std::vector< value_type > &  values 
) const

Return the values of the selected scalar component of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

This function is the equivalent of the FEValuesBase::get_function_values function but it only works on the selected scalar component.

template<int dim, int spacedim = dim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
std::vector< gradient_type > &  gradients 
) const

Return the gradients of the selected scalar component of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

This function is the equivalent of the FEValuesBase::get_function_gradients function but it only works on the selected scalar component.

template<int dim, int spacedim = dim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
std::vector< hessian_type > &  hessians 
) const

Return the Hessians of the selected scalar component of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

This function is the equivalent of the FEValuesBase::get_function_hessians function but it only works on the selected scalar component.

template<int dim, int spacedim = dim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
std::vector< value_type > &  laplacians 
) const

Return the Laplacians of the selected scalar component of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called. The Laplacians are the trace of the Hessians.

This function is the equivalent of the FEValuesBase::get_function_laplacians function but it only works on the selected scalar component.

Member Data Documentation

template<int dim, int spacedim = dim>
const FEValuesBase<dim,spacedim>& FEValuesViews::Scalar< dim, spacedim >::fe_values
private

A reference to the FEValuesBase object we operate on.

Definition at line 312 of file fe_values.h.

template<int dim, int spacedim = dim>
const unsigned int FEValuesViews::Scalar< dim, spacedim >::component
private

The single scalar component this view represents of the FEValuesBase object.

Definition at line 318 of file fe_values.h.

template<int dim, int spacedim = dim>
std::vector<ShapeFunctionData> FEValuesViews::Scalar< dim, spacedim >::shape_function_data
private

Store the data about shape functions.

Definition at line 323 of file fe_values.h.


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