imaging::TvFlowStep< fem_types > Class Template Reference
[Equations]

Implements an implicit time step of a total variation flow equation. More...

#include <TvFlowStep.hpp>

Inheritance diagram for imaging::TvFlowStep< fem_types >:

imaging::SimpleEquationInterface< fem_types >

List of all members.

Public Types

typedef ublas::fixed_matrix
< float_t,
fem_types::data_dimension,
fem_types::data_dimension > 
matrix_coefficient_t
typedef ublas::fixed_vector
< float_t,
fem_types::data_dimension > 
vector_coefficient_t

Public Member Functions

const ublas::vector< float_t > & input () const
void set_input (boost::shared_ptr< ublas::vector< float_t > > input_ptr)
float_t step_size () const
void set_step_size (float_t step_size)
float_t epsilon () const
void set_epsilon (float_t epsilon)
void stiffness_matrix (std::size_t integrator_node, const FemKernel< fem_types > &kernel, matrix_coefficient_t &A, ublas::fixed_vector< float_t, fem_types::data_dimension > &a, ublas::fixed_vector< float_t, fem_types::data_dimension > &b, float_t &c) const
void force_vector (std::size_t integrator_node, const FemKernel< fem_types > &kernel, float_t &f, ublas::fixed_vector< float_t, fem_types::data_dimension > &g) const
bool sanity_check_stiffness_matrix (const FemKernel< fem_types > &kernel, std::string &error_message) const
bool sanity_check_force_vector (const FemKernel< fem_types > &kernel, std::string &error_message) const

Static Public Attributes

static const bool a_active = false
static const bool b_active = false
static const bool c_active = true
static const bool f_active = true
static const bool g_active = false
static const size_t boundary_data_type = SimpleEquationInterface<fem_types>::NO_BOUNDARY_DATA


Detailed Description

template<class fem_types>
class imaging::TvFlowStep< fem_types >

Implements an implicit time step of a total variation flow equation.

This class implements the equation

\[ \frac{u - u^0}\delta = \nabla \cdot (\frac{\nabla u}{| \nabla u |})\,. \]

Here $u$ and $u^0$ are scalar functions on the problem domain and $\delta > 0$. The user must provide the initial value $u^0$ and the time step $\delta$ to assembly the system of linear equations which can be solved for $u$.


Member Typedef Documentation

template<class fem_types>
typedef ublas::fixed_matrix<float_t, fem_types::data_dimension, fem_types::data_dimension> imaging::TvFlowStep< fem_types >::matrix_coefficient_t

Defines the type of matrix as passed to A(). If A() returns a diagonal matrix choosing an appropriate type to store this matrix might yield a speed-up (in particular in higher dimensions; in the plane the difference will most probably neglible).

Reimplemented from imaging::SimpleEquationInterface< fem_types >.


Member Function Documentation

template<class fem_types>
const ublas::vector<float_t>& imaging::TvFlowStep< fem_types >::input (  )  const [inline]

Returns a reference to the vector containing the initial value $u^0$.

template<class fem_types>
void imaging::TvFlowStep< fem_types >::set_input ( boost::shared_ptr< ublas::vector< float_t > >  input_ptr  )  [inline]

Set the vector containing the initial value $u^0$.

template<class fem_types>
float_t imaging::TvFlowStep< fem_types >::step_size (  )  const [inline]

Returns the current time step $\delta$.

template<class fem_types>
void imaging::TvFlowStep< fem_types >::set_step_size ( float_t  step_size  )  [inline]

Set the size of the time step.

template<class fem_types>
float_t imaging::TvFlowStep< fem_types >::epsilon (  )  const [inline]

Returns the regularization parameter $\epsilon$ for the norm of the gradient.

template<class fem_types>
void imaging::TvFlowStep< fem_types >::set_epsilon ( float_t  epsilon  )  [inline]

Set the regularization parameter $\epsilon$ for the norm of the gradient.

template<class fem_types>
void imaging::TvFlowStep< fem_types >::stiffness_matrix ( std::size_t  integrator_node,
const FemKernel< fem_types > &  kernel,
matrix_coefficient_t A,
ublas::fixed_vector< float_t, fem_types::data_dimension > &  a,
ublas::fixed_vector< float_t, fem_types::data_dimension > &  b,
float_t c 
) const [inline]

Evaluates $A$, $a$, $b$ and $c$ in integrator_node on the current element of kernel.

As kernel is set to the current element, in the implementation of this function all relevant values of the shape functions and the element transform can be retrieved from kernel via integrator_node. One can also obtain the FE grid from the kernel (FemKernel::grid()) and use its interpolation methods (passing the already correctly initialized kernel to them).

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

References imaging::FemKernel< fem_types >::grid(), and imaging::square().

template<class fem_types>
void imaging::TvFlowStep< fem_types >::force_vector ( std::size_t  integrator_node,
const FemKernel< fem_types > &  kernel,
float_t f,
ublas::fixed_vector< float_t, fem_types::data_dimension > &  g 
) const [inline]

Evaluates $f$ and $g$ in integrator_node on the current element of kernel.

As kernel is set to the current element, in the implementation of this function all relevant values of the shape functions and the element transform can be retrieved from kernel via integrator_node. One can also obtain the FE grid from the kernel (FemKernel::grid()) and use its interpolation methods (passing the already correctly initialized kernel to them).

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

References imaging::FemKernel< fem_types >::grid().

template<class fem_types>
bool imaging::TvFlowStep< fem_types >::sanity_check_stiffness_matrix ( const FemKernel< fem_types > &  kernel,
std::string &  error_message 
) const [inline]

Checks if the dimension of the data corresponds to the dimension of the grid stored in kernel for the assembly of the stiffness matrix.

This function is called by the assemble routine right before the assembly of the stiffness matrix. It should return false if data which is necessary for the assembly of the stiffness matrix is still missing. Otherwise, return true;

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

References imaging::FemKernel< fem_types >::grid().

template<class fem_types>
bool imaging::TvFlowStep< fem_types >::sanity_check_force_vector ( const FemKernel< fem_types > &  kernel,
std::string &  error_message 
) const [inline]

Checks if the dimension of the data corresponds to the dimension of the grid stored in kernel for the assembly of the force vector.

This function is called by the assemble routine right before the assembly of the force vector. It should return false if data which is necessary for the assembly of the force vector is still missing. Otherwise, return true;

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

References imaging::FemKernel< fem_types >::grid().


Member Data Documentation

template<class fem_types>
const bool imaging::TvFlowStep< fem_types >::a_active = false [static]

Must be set to true if $a$ can be non-zero. If it is false, $a$ will not be evaluated in the assembly functions to save computation time.

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

template<class fem_types>
const bool imaging::TvFlowStep< fem_types >::b_active = false [static]

Must be set to true if $b$ can be non-zero. If it is false, $b$ will not be evaluated in the assembly functions to save computation time.

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

template<class fem_types>
const bool imaging::TvFlowStep< fem_types >::c_active = true [static]

Must be set to true if $c$ can be non-zero. If it is false, $c$ will not be evaluated in the assembly functions to save computation time.

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

template<class fem_types>
const bool imaging::TvFlowStep< fem_types >::f_active = true [static]

Must be set to true if $f$ can be non-zero. If it is false, $f$ will not be evaluated in the assembly functions to save computation time.

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

template<class fem_types>
const bool imaging::TvFlowStep< fem_types >::g_active = false [static]

Must be set to true if $g$ can be non-zero. If it is false, $g$ will not be evaluated in the assembly functions to save computation time.

Reimplemented from imaging::SimpleEquationInterface< fem_types >.

template<class fem_types>
const size_t imaging::TvFlowStep< fem_types >::boundary_data_type = SimpleEquationInterface<fem_types>::NO_BOUNDARY_DATA [static]

The type of boundary conditions for this equation as defined in SimpleEquationInterface::boundary_data_types.

Reimplemented from imaging::SimpleEquationInterface< fem_types >.


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

Generated on Tue Feb 10 10:01:31 2009 for imaging2 by  doxygen 1.5.5