ug4
|
Class for the computation of the enhanced version ('Frolkovic-type') of the consistent gravity. More...
#include <consistent_gravity.h>
Public Member Functions | |
template<int dim> | |
void | compute (MathVector< dim > &ConsistentGravity, const MathVector< refDim > &LocalCoord, const MathMatrix< dim, refDim > &JTInv, const MathVector< refDim > *vLocalGrad, const MathVector< refDim > *vConsGravity) |
computation of the consistent gravity at a given point | |
void | prepare (MathVector< 2 > *vConsGravity, const int n_co, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
spacialization of the method for faces (reference dimension 2) | |
void | prepare (MathVector< 3 > *vConsGravity, const int n_co, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
spacialization of the method for volumes (reference dimension 3) | |
template<int dim> | |
void | prepare (MathVector< refDim > *vConsGravity, const int n_co, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners, cf. the specializations | |
StdLinConsistentGravityX () | |
constructor | |
![]() | |
template<int dim> | |
void | compute (MathVector< dim > &ConsistentGravity, const MathVector< refDim > &LocalCoord, const MathMatrix< dim, refDim > &JTInv, const MathVector< refDim > *vLocalGrad, const MathVector< refDim > *vConsGravity) |
computation of the consistent gravity at a given point | |
void | prepare (MathVector< 1 > *vConsGravity, const int n_co, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
spacialization of the method for edges (reference dimension 1) | |
void | prepare (MathVector< 2 > *vConsGravity, const int n_co, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
spacialization of the method for faces (reference dimension 2) | |
void | prepare (MathVector< 3 > *vConsGravity, const int n_co, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
spacialization of the method for volumes (reference dimension 3) | |
template<int dim> | |
void | prepare (MathVector< refDim > *vConsGravity, const int n_co, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners, cf. the specializations | |
StdLinConsistentGravity () | |
constructor (sets the 'not init.' flag) | |
Protected Member Functions | |
template<typename refElem , int dim> | |
void | prepare_simplex (MathVector< refDim > *vConsGravity, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the extended version of the corner consistent gravity for simplices (only in full dimension!) | |
![]() | |
template<int dim> | |
void | prepare_edge (MathVector< 1 > *vConsGravity, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners of an edge | |
template<int dim> | |
void | prepare_hexahedron (MathVector< 3 > *vConsGravity, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners of a hexahedron | |
template<int dim> | |
void | prepare_prism (MathVector< 3 > *vConsGravity, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners of a prism | |
template<int dim> | |
void | prepare_pyramid (MathVector< 3 > *vConsGravity, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners of a pyramid | |
template<int dim> | |
void | prepare_quadrilateral (MathVector< 2 > *vConsGravity, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners of a quadrilateral | |
template<int dim> | |
void | prepare_tetrahedron (MathVector< 3 > *vConsGravity, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners of a tetrahedron | |
template<int dim> | |
void | prepare_triangle (MathVector< 2 > *vConsGravity, const MathVector< dim > *vCorners, const number *vDensity, const MathVector< dim > &PhysicalGravity) |
computation of the primary function for the consistent gravity at corners of a triangle | |
Private Types | |
typedef StdLinConsistentGravity< refDim > | base_type |
Additional Inherited Members | |
![]() | |
int | m_nCo |
number of corners of the element for which the object is init. (0 if not init) | |
Class for the computation of the enhanced version ('Frolkovic-type') of the consistent gravity.
Implementation of the enhanced ('Frolkovic-type') version of the consistent gravity for simplices (triangles and tetrahedra) in the case of the gravity force parallel to the z-axis. For all other elements and other gravities, the same method as in StdLinConsistentGravity is used.
References:
refDim | dimensionality of the reference element (e.g. 2 for triangles, 3 for tetrahedra) |
|
private |
|
inline |
constructor
|
inline |
computation of the consistent gravity at a given point
ConsistentGravity | where to save the vector |
LocalCoord | local coordinates of the point |
JTInv | inverse transposed Jacobian |
vLocalGrad | gradients of the shape functions at the given point |
vConsGravity | primary function of the consistent gravity at corners |
References ug::MatVecMult(), UG_ASSERT, and ug::VecSet().
void ug::StdLinConsistentGravityX< 2 >::prepare | ( | MathVector< 2 > * | vConsGravity, |
const int | n_co, | ||
const MathVector< dim > * | vCorners, | ||
const number * | vDensity, | ||
const MathVector< dim > & | PhysicalGravity | ||
) |
spacialization of the method for faces (reference dimension 2)
vConsGravity | where to save the values (n_co vectors) |
n_co | number of corners of the element |
vCorners | (global) coordinates of the corners (n_co vectors) |
vDensity | corner density (n_co scalar values) |
PhysicalGravity | the gravity vector |
void ug::StdLinConsistentGravityX< 3 >::prepare | ( | MathVector< 3 > * | vConsGravity, |
const int | n_co, | ||
const MathVector< dim > * | vCorners, | ||
const number * | vDensity, | ||
const MathVector< dim > & | PhysicalGravity | ||
) |
spacialization of the method for volumes (reference dimension 3)
vConsGravity | where to save the values (n_co vectors) |
n_co | number of corners of the element |
vCorners | (global) coordinates of the corners (n_co vectors) |
vDensity | corner density (n_co scalar values) |
PhysicalGravity | the gravity vector |
|
inline |
computation of the primary function for the consistent gravity at corners, cf. the specializations
vConsGravity | where to save the values (n_co vectors) |
n_co | number of corners of the element |
vCorners | (global) coordinates of the corners (n_co vectors) |
vDensity | corner density (n_co scalar values) |
PhysicalGravity | the gravity vector |
|
inlineprotected |
computation of the extended version of the corner consistent gravity for simplices (only in full dimension!)
vConsGravity | where to save the values (dim+1 vectors) |
vCorners | (global) coordinates of the corners (dim+1 vectors) |
vDensity | corner density (dim+1 scalar values) |
PhysicalGravity | the gravity vector (MUST BE (0, ..., 0, g)) |
References ug::ReferenceMapping< TRefElem, TWorldDim >::jacobian_transposed_inverse(), ug::TransposedMatVecMult(), ug::ReferenceMapping< TRefElem, TWorldDim >::update(), and ug::VecSubtract().