Electrowetting on dielectric
Public Types | Public Member Functions | Public Attributes | Protected Attributes | Private Member Functions
AsFunction< dim > Class Template Reference
Collaboration diagram for AsFunction< dim >:
Collaboration graph
[legend]

List of all members.

Public Types

enum  FunctionType { isScalarValued, isVectorValued }

Public Member Functions

 AsFunction (const hp::DoFHandler< dim > &_dh, const TrilinosWrappers::Vector &_sol, const char *_name, const FunctionType _type=isScalarValued, const unsigned fc=0)
 AsFunction (const AsFunction< dim > &_asfunction)
 ~AsFunction ()
void reset (const hp::QCollection< dim > &q, const UpdateFlags f)
void reset ()
void reinit (const typename hp::DoFHandler< dim >::active_cell_iterator &cell)
void reset_face (const hp::QCollection< dim-1 > &q, const UpdateFlags f)
void reset_face ()
void reinit (const typename hp::DoFHandler< dim >::active_cell_iterator &cell, const unsigned face)
void get_function_values (std::vector< double > &values) const
void get_function_gradients (std::vector< Tensor< 1, dim > > &values) const
void get_function_values (std::vector< Tensor< 1, dim > > &values) const
void get_function_gradients (std::vector< Tensor< 2, dim > > &values) const
void get_function_divergences (std::vector< double > &values) const
void get_function_face_values (std::vector< double > &values) const
void get_function_face_gradients (std::vector< Tensor< 1, dim > > &values) const
void get_function_face_values (std::vector< Tensor< 1, dim > > &values) const
void get_function_face_gradients (std::vector< Tensor< 2, dim > > &values) const
const TrilinosWrappers::Vectorget_vec () const
const double operator() (unsigned i) const
const hp::DoFHandler< dim > & get_dh () const
hp::DoFHandler< dim >
::active_cell_iterator 
get_begin () const
hp::DoFHandler< dim >
::active_cell_iterator 
get_end () const

Public Attributes

const std::string name
const
DataComponentInterpretation::DataComponentInterpretation 
interpretation

Protected Attributes

hp::QCollection< dim > * quad
UpdateFlags flags
hp::FEValues< dim > * fe_val
hp::QCollection< dim-1 > * f_quad
UpdateFlags f_flags
hp::FEFaceValues< dim > * fe_face_val
SmartPointer< const
hp::DoFHandler< dim > > 
dh
SmartPointer< const
TrilinosWrappers::Vector
sol
const FunctionType type
const unsigned first_component

Private Member Functions

 AsFunction ()

Detailed Description

template<int dim>
class AsFunction< dim >

When we have more than one coupled PDEs, we usually need to use the computed solution of a subproblem as right hand side or coefficient of the other subproblems. This class implements this technology, by storing a FEValues object, being able to loop over all cells and, more importantly, call the get_function_* methods with a fixed global vector which, usually, is the solution of a given problem or an easily derived quantity.

Definition at line 22 of file AsFunction.h.


Member Enumeration Documentation

template<int dim>
enum AsFunction::FunctionType

In deal.II the FEValuesExtractors do not have a common type therefore we will need to keep track of the type of function (scalar/vector valued) and treat the get_function_* methods accordingly

Enumerator:
isScalarValued 

The function is scalar valued.

isVectorValued 

The function is vector valued.

Definition at line 29 of file AsFunction.h.


Constructor & Destructor Documentation

template<int dim>
AsFunction< dim >::AsFunction ( const hp::DoFHandler< dim > &  _dh,
const TrilinosWrappers::Vector _sol,
const char *  _name,
const FunctionType  _type = isScalarValued,
const unsigned  fc = 0 
)
Parameters:
_dh: The DoFHandler object that is needed to interpret vectors as functions.
_sol: The coefficient vector that represent the discrete function.
_name: The name of this function. Used for plotting.
_type: Indicates if the function must be interpreted as scalar or vector valued
fc: First component. This is used if, for instance, the finite element space is vector valued and we want only one of the components of the space.
template<int dim>
AsFunction< dim >::AsFunction ( const AsFunction< dim > &  _asfunction)

Copy constructor

template<int dim>
AsFunction< dim >::~AsFunction ( )

Clear allocated memory

template<int dim>
AsFunction< dim >::AsFunction ( ) [private]

Disable default constructor


Member Function Documentation

template<int dim>
void AsFunction< dim >::reset ( const hp::QCollection< dim > &  q,
const UpdateFlags  f 
)

Reset the object,

Parameters:
q: is the new quadrature
f: the update flags
template<int dim>
void AsFunction< dim >::reset ( )

Clear the object

template<int dim>
void AsFunction< dim >::reinit ( const typename hp::DoFHandler< dim >::active_cell_iterator &  cell)

Reinit to a new cell

template<int dim>
void AsFunction< dim >::reset_face ( const hp::QCollection< dim-1 > &  q,
const UpdateFlags  f 
)

Reset the face object

Parameters:
q: is the new quadrature
f: the update flags
template<int dim>
void AsFunction< dim >::reset_face ( )

Clears the face object

template<int dim>
void AsFunction< dim >::reinit ( const typename hp::DoFHandler< dim >::active_cell_iterator &  cell,
const unsigned  face 
)

Reinit to a new face.

Parameters:
cell: the cell
face: the face
template<int dim>
void AsFunction< dim >::get_function_values ( std::vector< double > &  values) const

Maps to FEValues.get_function_values for scalars

template<int dim>
void AsFunction< dim >::get_function_gradients ( std::vector< Tensor< 1, dim > > &  values) const

Maps to FEValues.get_function_gradients for scalars

template<int dim>
void AsFunction< dim >::get_function_values ( std::vector< Tensor< 1, dim > > &  values) const

Maps to FEValues.get_function_values for vectors

template<int dim>
void AsFunction< dim >::get_function_gradients ( std::vector< Tensor< 2, dim > > &  values) const

Maps to FEValues.get_function_gradients for vectors

template<int dim>
void AsFunction< dim >::get_function_divergences ( std::vector< double > &  values) const

Maps to FEValues.get_function_divergences (valid only for vectors)

template<int dim>
void AsFunction< dim >::get_function_face_values ( std::vector< double > &  values) const

Maps to FEFaceValues.get_function_values for scalars

template<int dim>
void AsFunction< dim >::get_function_face_gradients ( std::vector< Tensor< 1, dim > > &  values) const
template<int dim>
void AsFunction< dim >::get_function_face_values ( std::vector< Tensor< 1, dim > > &  values) const

Maps to FEFaceValues.get_function_values for vectors

template<int dim>
void AsFunction< dim >::get_function_face_gradients ( std::vector< Tensor< 2, dim > > &  values) const
template<int dim>
const TrilinosWrappers::Vector& AsFunction< dim >::get_vec ( ) const

Raw access to the data vector

template<int dim>
const double AsFunction< dim >::operator() ( unsigned  i) const

Raw access to an entry of the data vector

template<int dim>
const hp::DoFHandler<dim>& AsFunction< dim >::get_dh ( ) const

Get the associated DoFHandler object

template<int dim>
hp::DoFHandler<dim>::active_cell_iterator AsFunction< dim >::get_begin ( ) const

Get the beginning of the cell iterator associated with this DoFHandler

template<int dim>
hp::DoFHandler<dim>::active_cell_iterator AsFunction< dim >::get_end ( ) const

Get the end of the cell iterator associated with this DoFHandler


Member Data Documentation

template<int dim>
hp::QCollection<dim>* AsFunction< dim >::quad [mutable, protected]

The internal cell quadrature rule.

Definition at line 141 of file AsFunction.h.

template<int dim>
UpdateFlags AsFunction< dim >::flags [mutable, protected]

The internal cell update flags

Definition at line 145 of file AsFunction.h.

template<int dim>
hp::FEValues<dim>* AsFunction< dim >::fe_val [mutable, protected]

The internal FEValues object

Definition at line 149 of file AsFunction.h.

template<int dim>
hp::QCollection<dim-1>* AsFunction< dim >::f_quad [mutable, protected]

The internal face quadrature rule

Definition at line 153 of file AsFunction.h.

template<int dim>
UpdateFlags AsFunction< dim >::f_flags [mutable, protected]

The internal face update flags

Definition at line 157 of file AsFunction.h.

template<int dim>
hp::FEFaceValues<dim>* AsFunction< dim >::fe_face_val [mutable, protected]

The internal FEFaceValues object

Definition at line 161 of file AsFunction.h.

template<int dim>
SmartPointer< const hp::DoFHandler<dim> > AsFunction< dim >::dh [protected]

Points to DoFHandler

Definition at line 165 of file AsFunction.h.

template<int dim>
SmartPointer< const TrilinosWrappers::Vector > AsFunction< dim >::sol [protected]

Points to the coefficient vector

Definition at line 169 of file AsFunction.h.

template<int dim>
const FunctionType AsFunction< dim >::type [protected]

The type of function that we have (scalar/vector)

Definition at line 173 of file AsFunction.h.

template<int dim>
const unsigned AsFunction< dim >::first_component [protected]

The first component of the function (meaningful if the DoFHandler is vector valued)

Definition at line 177 of file AsFunction.h.

template<int dim>
const std::string AsFunction< dim >::name

The function name

Definition at line 182 of file AsFunction.h.

The interpetation of the function (scalar/vector)

Definition at line 186 of file AsFunction.h.


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