Evolving surface finite element method  v0.3.0-14-g3598512
Numerical experiments for my papers
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
Esfem::Impl::MCF_op Class Reference

#include <secOrd_op_solutionDriven_impl.h>

Inheritance diagram for Esfem::Impl::MCF_op:
Inheritance graph
[legend]
Collaboration diagram for Esfem::Impl::MCF_op:
Collaboration graph
[legend]

Public Types

using Scalar_fef = Esfem::Grid::Scal_FEfun::Dune_FEfun
 $ f\colon \R^3 \to \R $
 
using Vector_fef = Esfem::Grid::Vec_FEfun::Dune_FEfun
 $ f\colon \R^3 \to \R^3 $
 
template<typename T >
using Local_function = typename T::LocalFunctionType
 Finite element function restricted to a triangle.
 
template<typename T >
using Domain = typename Local_function< T >::DomainType
 
template<typename T >
using Range = typename Local_function< T >::RangeType
 
using RangeField = typename Local_function< Vector_fef >::RangeFieldType
 Simply $ \R $.
 
using Geometry = Vector_fef::DiscreteFunctionSpaceType::IteratorType::Entity::Geometry
 Geometry means access to the finite element.
 
using Grid_part = Vector_fef::GridPartType
 Auxiliary template argument.
 
using Quadrature = Dune::Fem::CachingQuadrature< Grid_part, 0 >
 
template<typename T >
using Jacobian_range = typename Local_function< T >::JacobianRangeType
 Jacobian of the local finite element function. More...
 

Public Member Functions

 MCF_op (const Io::Parameter &, const Grid::Grid_and_time &, const Scalar_fef &u_input)
 
 MCF_op (const MCF_op &)=delete
 
MCF_opoperator= (const MCF_op &)=delete
 
void operator() (const Vector_fef &rhs, Vector_fef &lhs) const override
 Applying the mean curvature operator. The cg solver uses this. More...
 
void brusselator_rhs (const Vector_fef &rhs, Vector_fef &lhs) const
 Generates rhs for the linear system. More...
 

Static Public Attributes

Variable containing dimensions
static constexpr int dim_vec_domain = Domain<Vector_fef>::dimension
 
static constexpr int dim_vec_range = Range<Vector_fef>::dimension
 
static constexpr int dim_scalar_domain = Domain<Scalar_fef>::dimension
 
static constexpr int dim_scalar_range = Range<Scalar_fef>::dimension
 

Private Member Functions

void mcf_lhs_matrixFree_assembly (const Geometry &, const Quadrature &, const Local_function< Vector_fef > &, Local_function< Vector_fef > &) const
 Used by operator() More...
 
void mcf_rhs_matrixFree_assembly (const Geometry &, const Quadrature &, const Local_function< Vector_fef > &, const Local_function< Scalar_fef > &u_loc, Local_function< Vector_fef > &) const
 Used by rhs() More...
 
Range< Vector_fefsurface_normal (const Geometry &) const
 Calculates the outwards pointing normal vector.
 

Private Attributes

const double alpha
 Velocity Laplace regularization parameter by Lubich.
 
const double delta
 Tumor growth parameter.
 
const double epsilon
 Mean curvature regularization parameter by Elliott.
 
const double eps
 Generic tolerance.
 
const Dune::Fem::TimeProviderBase & tp
 Time step provider.
 
const Scalar_fefu
 Concentration of the growth promoting chemical substance.
 

Detailed Description

See also
Esfem::SecOrd_op::Solution_driven

Definition at line 40 of file secOrd_op_solutionDriven_impl.h.

Member Typedef Documentation

template<typename T >
using Esfem::Impl::MCF_op::Jacobian_range = typename Local_function<T>::JacobianRangeType

Jacobian of the local finite element function.

Warning
Jacobian_range is used like this: jac[0][1].

Definition at line 116 of file secOrd_op_solutionDriven_impl.h.

Constructor & Destructor Documentation

MCF_op::MCF_op ( const Io::Parameter p,
const Grid::Grid_and_time &  g,
const Scalar_fef u_input 
)
See also
Constructor of Esfem::SecOrd_op::Solution_driven

Definition at line 58 of file secOrd_op_solutionDriven_impl.cpp.

Member Function Documentation

void MCF_op::brusselator_rhs ( const Vector_fef rhs,
Vector_fef lhs 
) const

Generates rhs for the linear system.

The new value of the finite element function will be

\begin{equation*} (M_3^n + \alpha A_3^n) \nodalValue{X}^n + \tau \delta M_3^n(\nodalValue{u}^n, \nodalValue{\surfaceNormal}). \end{equation*}

Definition at line 83 of file secOrd_op_solutionDriven_impl.cpp.

void MCF_op::mcf_lhs_matrixFree_assembly ( const Geometry g,
const Quadrature &  q,
const Local_function< Vector_fef > &  cf,
Local_function< Vector_fef > &  f 
) const
private

Used by operator()

The matrix vector formulation reads as

\begin{equation*} \parentheses[\big]{M + (\alpha + \varepsilon \tau) A} \nodalValue{X} \end{equation*}

Definition at line 113 of file secOrd_op_solutionDriven_impl.cpp.

void MCF_op::mcf_rhs_matrixFree_assembly ( const Geometry g,
const Quadrature &  q,
const Local_function< Vector_fef > &  cf,
const Local_function< Scalar_fef > &  u_loc,
Local_function< Vector_fef > &  f 
) const
private

Used by rhs()

The matrix vector formulation reads as

\begin{equation*} (M + \alpha A) \nodalValue{X} + \tau \delta M(\nodalValue{u}, \surfaceNormal) \end{equation*}

Definition at line 131 of file secOrd_op_solutionDriven_impl.cpp.

void MCF_op::operator() ( const Vector_fef rhs,
Vector_fef lhs 
) const
override

Applying the mean curvature operator. The cg solver uses this.

The precise formulation reads as

\begin{equation*} \parentheses[\big]{M_3^n + (\alpha + \varepsilon\tau) A_3^n} \nodalValue{X}^{n+1} \end{equation*}

Definition at line 69 of file secOrd_op_solutionDriven_impl.cpp.


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