ug4
ug::StdLinConsistentGravity< refDim > Class Template Reference

Class for the computation of the standard version ('Voss-Souza-type') of the consistent gravity. More...

#include <consistent_gravity.h>

+ Inheritance diagram for ug::StdLinConsistentGravity< refDim >:

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 More...
 
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) More...
 
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) More...
 
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) More...
 
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 More...
 
 StdLinConsistentGravity ()
 constructor (sets the 'not init.' flag) More...
 

Protected Member Functions

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 More...
 
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 More...
 
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 More...
 
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 More...
 
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 More...
 
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 More...
 
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 More...
 

Protected Attributes

int m_nCo
 number of corners of the element for which the object is init. (0 if not init) More...
 

Static Private Attributes

static const size_t _X_ = 0
 
static const size_t _Y_ = 1
 
static const size_t _Z_ = 2
 

Detailed Description

template<int refDim>
class ug::StdLinConsistentGravity< refDim >

Class for the computation of the standard version ('Voss-Souza-type') of the consistent gravity.

Density driven flow models involve the convection velocity depending on the gradient of the pressure and the density-dependent gravity force, typically of the form \(\mathbf{v} = \mathbf{K} (- \nabla p + \rho \mathbf{g}) / \mu\), where \(\mathbf{K}\) is a matrix, \(\mu\) a scalar, \(\mathbf{g}\) gravity vector (all independent of \(p\) or \(\rho\)), \(p\) is pressure, \(rho\) is density (depending on the concentration, ...). (Both the pressure and the concentration are unknown grid functions in the PDEs.) Treating \(p\) and \(\rho\) as grid functions of the same class (like the piecewise linear/bilinear interpolations of nodal values) leads to oscillatory solutions because \(\nabla p\) and \(\rho \mathbf{g}\) belong to different classes of grid functions. In particular, if \(p\) is a piecewise linear function, \(\nabla p\) can never cancel the of the gravitaty force if \(\rho\) increases linearly from the top down, although this happens in the analytical solution. The idea of the consistent gravity is to consider a vector function \(\mathbf{h} = (h_x, h_y, \dots)\) that is in some sence a primitive function for \(\rho \mathbf{g}\) so that the velocity can be written as \(\mathbf{v} = - \mathbf{K} (p_x - h_x, p_y - h_y, \dots) / \mu\). In the discretization, the function \(\mathbf{h}\) should belong to the same class of the grid functions as \(p\). In the present implementation, \(\mathbf{h}\) is a piecewise linear/bilinear grid function.

Class method 'prepare' computes the values of \(h\) at the corners of an element. Using these values, method 'compute' computes the consistent gravity force \(\rho \mathbf{g}\) (not \(\mathbf{h}\)!) at any given point.

Remarks
(Consistent gravity in the computation of the Jacobian) The consistent gravity computed by the implemented methods depends linearly on the density. Thus, to compute the derivative of the consistent gravity w.r.t. the density at one of the corners of the element, set the density to 1 at that corner and to 0 at all other corners, then call the functions. Then 'prepare' prepares the nodal values of the derivatives and 'compute' computes the derivative itself.

Alternatively, you can set the density at that corner to the derivative of the density w.r.t. the concentration (instead of 1). Then you get directly the derivative of the consistent gravity w.r.t. the concentration at that node.

Remarks
There is an enhanced version of the consistent gravity.
See also
StdLinConsistentGravityX

References:

  • P. Frolkovic, P. Knabner, Consistent Velocity Approximations in Finite Element or Volume Discretizations of Density Driven Flow, In: Computational Methods in WaterResources XI, Vol. 1 (A.A. Aldama et al., eds.), Computational Mechanics Publication, Southhampten, 1996, p. 93-100
Template Parameters
refDimdimensionality of the reference element (e.g. 2 for triangles, 3 for tetrahedra)

Constructor & Destructor Documentation

◆ StdLinConsistentGravity()

template<int refDim>
ug::StdLinConsistentGravity< refDim >::StdLinConsistentGravity ( )
inline

constructor (sets the 'not init.' flag)

Member Function Documentation

◆ compute()

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::compute ( MathVector< dim > &  ConsistentGravity,
const MathVector< refDim > &  LocalCoord,
const MathMatrix< dim, refDim > &  JTInv,
const MathVector< refDim > *  vLocalGrad,
const MathVector< refDim > *  vConsGravity 
)
inline

computation of the consistent gravity at a given point

Parameters
ConsistentGravitywhere to save the vector
LocalCoordlocal coordinates of the point
JTInvinverse transposed Jacobian
vLocalGradgradients of the shape functions at the given point
vConsGravityprimary function of the consistent gravity at corners

References ug::StdLinConsistentGravity< refDim >::m_nCo, ug::MatVecMult(), UG_ASSERT, and ug::VecSet().

◆ prepare() [1/4]

void ug::StdLinConsistentGravity< 1 >::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)

Parameters
vConsGravitywhere to save the values (n_co vectors)
n_conumber of corners of the element (should be 2)
vCorners(global) coordinates of the corners (n_co vectors)
vDensitycorner density (n_co scalar values)
PhysicalGravitythe gravity vector

References UG_ASSERT.

◆ prepare() [2/4]

void ug::StdLinConsistentGravity< 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)

Parameters
vConsGravitywhere to save the values (n_co vectors)
n_conumber of corners of the element
vCorners(global) coordinates of the corners (n_co vectors)
vDensitycorner density (n_co scalar values)
PhysicalGravitythe gravity vector

References UG_THROW.

◆ prepare() [3/4]

void ug::StdLinConsistentGravity< 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)

Parameters
vConsGravitywhere to save the values (n_co vectors)
n_conumber of corners of the element
vCorners(global) coordinates of the corners (n_co vectors)
vDensitycorner density (n_co scalar values)
PhysicalGravitythe gravity vector

References UG_THROW.

◆ prepare() [4/4]

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::prepare ( MathVector< refDim > *  vConsGravity,
const int  n_co,
const MathVector< dim > *  vCorners,
const number vDensity,
const MathVector< dim > &  PhysicalGravity 
)
inline

computation of the primary function for the consistent gravity at corners, cf. the specializations

Parameters
vConsGravitywhere to save the values (n_co vectors)
n_conumber of corners of the element
vCorners(global) coordinates of the corners (n_co vectors)
vDensitycorner density (n_co scalar values)
PhysicalGravitythe gravity vector

References dim, and UG_THROW.

◆ prepare_edge()

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::prepare_edge ( MathVector< 1 > *  vConsGravity,
const MathVector< dim > *  vCorners,
const number vDensity,
const MathVector< dim > &  PhysicalGravity 
)
inlineprotected

computation of the primary function for the consistent gravity at corners of an edge

Parameters
vConsGravitywhere to save the values (2 vectors)
vCorners(global) coordinates of the corners (2 vectors)
vDensitycorner density (2 scalar values)
PhysicalGravitythe gravity vector

References ug::StdLinConsistentGravity< refDim >::_X_, ug::ReferenceMapping< TRefElem, TWorldDim >::jacobian_transposed(), ug::MatVecMult(), ug::ReferenceMapping< TRefElem, TWorldDim >::update(), and ug::VecSet().

◆ prepare_hexahedron()

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::prepare_hexahedron ( MathVector< 3 > *  vConsGravity,
const MathVector< dim > *  vCorners,
const number vDensity,
const MathVector< dim > &  PhysicalGravity 
)
inlineprotected

computation of the primary function for the consistent gravity at corners of a hexahedron

Parameters
vConsGravitywhere to save the values (8 vectors)
vCorners(global) coordinates of the corners (8 vectors)
vDensitycorner density (8 scalar values)
PhysicalGravitythe gravity vector

References ug::StdLinConsistentGravity< refDim >::_X_, ug::StdLinConsistentGravity< refDim >::_Y_, ug::StdLinConsistentGravity< refDim >::_Z_, ug::ReferenceMapping< TRefElem, TWorldDim >::jacobian_transposed(), ug::MatVecMult(), and ug::ReferenceMapping< TRefElem, TWorldDim >::update().

◆ prepare_prism()

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::prepare_prism ( MathVector< 3 > *  vConsGravity,
const MathVector< dim > *  vCorners,
const number vDensity,
const MathVector< dim > &  PhysicalGravity 
)
inlineprotected

computation of the primary function for the consistent gravity at corners of a prism

Parameters
vConsGravitywhere to save the values (6 vectors)
vCorners(global) coordinates of the corners (6 vectors)
vDensitycorner density (6 scalar values)
PhysicalGravitythe gravity vector

References ug::StdLinConsistentGravity< refDim >::_X_, ug::StdLinConsistentGravity< refDim >::_Y_, ug::StdLinConsistentGravity< refDim >::_Z_, ug::ReferenceMapping< TRefElem, TWorldDim >::jacobian_transposed(), ug::MatVecMult(), and ug::ReferenceMapping< TRefElem, TWorldDim >::update().

◆ prepare_pyramid()

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::prepare_pyramid ( MathVector< 3 > *  vConsGravity,
const MathVector< dim > *  vCorners,
const number vDensity,
const MathVector< dim > &  PhysicalGravity 
)
inlineprotected

computation of the primary function for the consistent gravity at corners of a pyramid

TODO: Verify this implementation! Cf. UG3.

Parameters
vConsGravitywhere to save the values (5 vectors)
vCorners(global) coordinates of the corners (5 vectors)
vDensitycorner density (5 scalar values)
PhysicalGravitythe gravity vector

References ug::StdLinConsistentGravity< refDim >::_X_, ug::StdLinConsistentGravity< refDim >::_Y_, ug::StdLinConsistentGravity< refDim >::_Z_, ug::ReferenceMapping< TRefElem, TWorldDim >::jacobian_transposed(), ug::MatVecMult(), and ug::ReferenceMapping< TRefElem, TWorldDim >::update().

◆ prepare_quadrilateral()

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::prepare_quadrilateral ( MathVector< 2 > *  vConsGravity,
const MathVector< dim > *  vCorners,
const number vDensity,
const MathVector< dim > &  PhysicalGravity 
)
inlineprotected

computation of the primary function for the consistent gravity at corners of a quadrilateral

Parameters
vConsGravitywhere to save the values (4 vectors)
vCorners(global) coordinates of the corners (4 vectors)
vDensitycorner density (4 scalar values)
PhysicalGravitythe gravity vector

References ug::StdLinConsistentGravity< refDim >::_X_, ug::StdLinConsistentGravity< refDim >::_Y_, ug::ReferenceMapping< TRefElem, TWorldDim >::jacobian_transposed(), ug::MatVecMult(), ug::ReferenceMapping< TRefElem, TWorldDim >::update(), and ug::VecSet().

◆ prepare_tetrahedron()

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::prepare_tetrahedron ( MathVector< 3 > *  vConsGravity,
const MathVector< dim > *  vCorners,
const number vDensity,
const MathVector< dim > &  PhysicalGravity 
)
inlineprotected

computation of the primary function for the consistent gravity at corners of a tetrahedron

Parameters
vConsGravitywhere to save the values (4 vectors)
vCorners(global) coordinates of the corners (4 vectors)
vDensitycorner density (4 scalar values)
PhysicalGravitythe gravity vector

References ug::StdLinConsistentGravity< refDim >::_X_, ug::StdLinConsistentGravity< refDim >::_Y_, ug::StdLinConsistentGravity< refDim >::_Z_, ug::ReferenceMapping< TRefElem, TWorldDim >::jacobian_transposed(), ug::MatVecMult(), and ug::ReferenceMapping< TRefElem, TWorldDim >::update().

◆ prepare_triangle()

template<int refDim>
template<int dim>
void ug::StdLinConsistentGravity< refDim >::prepare_triangle ( MathVector< 2 > *  vConsGravity,
const MathVector< dim > *  vCorners,
const number vDensity,
const MathVector< dim > &  PhysicalGravity 
)
inlineprotected

computation of the primary function for the consistent gravity at corners of a triangle

Parameters
vConsGravitywhere to save the values (3 vectors)
vCorners(global) coordinates of the corners (3 vectors)
vDensitycorner density (3 scalar values)
PhysicalGravitythe gravity vector

References ug::StdLinConsistentGravity< refDim >::_X_, ug::StdLinConsistentGravity< refDim >::_Y_, ug::ReferenceMapping< TRefElem, TWorldDim >::jacobian_transposed(), ug::MatVecMult(), ug::ReferenceMapping< TRefElem, TWorldDim >::update(), and ug::VecSet().

Member Data Documentation

◆ _X_

◆ _Y_

◆ _Z_

◆ m_nCo

template<int refDim>
int ug::StdLinConsistentGravity< refDim >::m_nCo
protected

number of corners of the element for which the object is init. (0 if not init)

Referenced by ug::StdLinConsistentGravity< refDim >::compute().


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