DUNE-ACFEM (2.5.1)

Dune::ACFem::DeformationTensorOperatorParts< FunctionSpace > Class Template Reference

Define the a model where the flux-part is formed from the symmetric gradient. More...

#include <dune/acfem/models/modules/deformationtensormodel.hh>

Public Types

enum  StructureFlags
 
enum  ConstituentFlags
 
typedef Expression ExpressionType
 The type of the underlying expression.
 

Public Member Functions

template<class Entity , class Point >
void linearizedFlux (const RangeType &uBar, const JacobianRangeType &DuBar, const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, JacobianRangeType &flux) const
 Evaluate the linearized flux in local coordinates. More...
 
template<class Entity , class Point >
void fluxDivergence (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, const HessianRangeType &hessian, RangeType &result) const
 Compute the point-wise value of the flux-part of the operator, meaning the part of the differential operator which is multiplied by the derivative of the test function. More...
 
void setEntity (const Entity &entity) const
 
bool setIntersection (const Intersection &intersection) const
 
void flux (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, JacobianRangeType &flux) const
 
void source (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, RangeType &result) const
 
void linearizedSource (const RangeType &uBar, const JacobianRangeType &DuBar, const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, RangeType &result) const
 
void robinFlux (const Intersection &intersection, const Point &x, const DomainType &unitOuterNormal, const RangeType &value, RangeType &result) const
 
void linearizedRobinFlux (const RangeType &uBar, const Intersection &intersection, const Point &x, const DomainType &unitOuterNormal, const RangeType &value, RangeType &result) const
 
const ExpressionTypeexpression () const
 Return a const reference to the underlying expression.
 
ExpressionTypeexpression ()
 Return a mutable reference to the underlying expression.
 
ExpressionType operator* () const
 Return a copy from of the underlying expression.
 

Detailed Description

template<class FunctionSpace>
class Dune::ACFem::DeformationTensorOperatorParts< FunctionSpace >

Define the a model where the flux-part is formed from the symmetric gradient.

The flux implement the second-order contribution

\[ \int_\Omega (\partial_i U_j + \partial_j U_i)\,\partial_i\phi_j = \frac12 \int_\Omega (\partial_i U_j + \partial_j U_i)\,(\partial_i\phi_j+\partial_j\phi_i) \]

Of course, the above formula uses sum-convention. This is indeed symmetric. Mind the factor 1/2. The DeformationTensorModel::fluxDivergence() takes the form

\[ -\partial_i(\partial_i U_\alpha + \partial_\alpha U_i) = -\Delta U_\alpha - \partial_\alpha\nabla\cdot U. \]

So for divergence free vector fields this would only be the component-wise Laplacian.

Note
Obviously, this model make only sense for when dimDomain = dimRange = GridPart::dimensionworld. For GridPart::dimensionworld != GridPart::dimension one would also have to discuss what the symmetric gradient means in terms of tangential derivatives. In some sense also the vectorfield U would have to be tangential.

Member Enumeration Documentation

◆ ConstituentFlags

◆ StructureFlags

Member Function Documentation

◆ flux()

void Dune::ACFem::DefaultOperatorParts< Expression >::flux ( const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
JacobianRangeType &  flux 
) const
inlineinherited

◆ linearizedRobinFlux()

void Dune::ACFem::DefaultOperatorParts< Expression >::linearizedRobinFlux ( const RangeType &  uBar,
const Intersection &  intersection,
const Point &  x,
const DomainType &  unitOuterNormal,
const RangeType &  value,
RangeType &  result 
) const
inlineinherited

◆ linearizedSource()

void Dune::ACFem::DefaultOperatorParts< Expression >::linearizedSource ( const RangeType &  uBar,
const JacobianRangeType &  DuBar,
const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
RangeType &  result 
) const
inlineinherited

◆ robinFlux()

void Dune::ACFem::DefaultOperatorParts< Expression >::robinFlux ( const Intersection &  intersection,
const Point &  x,
const DomainType &  unitOuterNormal,
const RangeType &  value,
RangeType &  result 
) const
inlineinherited

◆ setEntity()

void Dune::ACFem::DefaultOperatorParts< Expression >::setEntity ( const Entity &  entity) const
inlineinherited

◆ setIntersection()

bool Dune::ACFem::DefaultOperatorParts< Expression >::setIntersection ( const Intersection &  intersection) const
inlineinherited

◆ source()

void Dune::ACFem::DefaultOperatorParts< Expression >::source ( const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
RangeType &  result 
) const
inlineinherited


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 3, 22:32, 2024)