KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
List of all members
Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients > Class Template Reference

#include <bingham_fluid.h>

Inheritance diagram for Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >:
Collaboration diagram for Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >:

Public Member Functions

Life Cycle
 BinghamFluid (IndexType NewId=0)
 Default constuctor. More...
 
 BinghamFluid (IndexType NewId, const NodesArrayType &ThisNodes)
 Constructor using an array of nodes. More...
 
 BinghamFluid (IndexType NewId, GeometryType::Pointer pGeometry)
 Constructor using a geometry object. More...
 
 BinghamFluid (IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
 Constuctor using geometry and properties. More...
 
 ~BinghamFluid () override
 Destructor. More...
 
Operations
Element::Pointer Create (IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
 Create a new element of this type. More...
 
Element::Pointer Create (IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
 Create a new element of this type. More...
 
Input and output
std::string Info () const override
 Turn back information as a string. More...
 
void PrintInfo (std::ostream &rOStream) const override
 Print information about this object. More...
 
void PrintData (std::ostream &rOStream) const override
 Print object's data. More...
 

Protected Member Functions

Protected Operations
double EffectiveViscosity (double Density, const TShapeFunctionValues &rN, const TShapeFunctionGradients &rDN_DX, double ElemSize, const ProcessInfo &rProcessInfo) override
 EffectiveViscosity Calculate the effective viscosity at given integration point using the Bingham constitutive model. More...
 

Type Definitions

typedef Node NodeType
 Node type (default is: Node) More...
 
typedef Geometry< NodeTypeGeometryType
 Geometry type (using with given NodeType) More...
 
typedef Geometry< NodeType >::PointsArrayType NodesArrayType
 Definition of nodes container type, redefined from GeometryType. More...
 
typedef Properties PropertiesType
 
typedef Vector VectorType
 Vector type for local contributions to the linear system. More...
 
typedef Matrix MatrixType
 Matrix type for local contributions to the linear system. More...
 
typedef std::size_t IndexType
 
typedef std::size_t SizeType
 
typedef std::vector< std::size_t > EquationIdVectorType
 
typedef std::vector< Dof< double >::Pointer > DofsVectorType
 
typedef PointerVectorSet< Dof< double >, IndexedObjectDofsArrayType
 
typedef Kratos::Vector ShapeFunctionsType
 Type for shape function values container. More...
 
typedef Kratos::Matrix ShapeFunctionDerivativesType
 Type for a matrix containing the shape function gradients. More...
 
typedef GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
 Type for an array of shape function gradient matrices. More...
 
 KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (BinghamFluid)
 

Serialization

class Serializer
 

Member Typedef Documentation

◆ DofsArrayType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef PointerVectorSet<Dof<double>, IndexedObject> Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::DofsArrayType

◆ DofsVectorType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef std::vector< Dof<double>::Pointer > Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::DofsVectorType

◆ EquationIdVectorType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef std::vector<std::size_t> Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::EquationIdVectorType

◆ GeometryType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef Geometry<NodeType> Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::GeometryType

Geometry type (using with given NodeType)

◆ IndexType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef std::size_t Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::IndexType

◆ MatrixType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef Matrix Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::MatrixType

Matrix type for local contributions to the linear system.

◆ NodesArrayType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef Geometry<NodeType>::PointsArrayType Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::NodesArrayType

Definition of nodes container type, redefined from GeometryType.

◆ NodeType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef Node Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::NodeType

Node type (default is: Node)

◆ PropertiesType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef Properties Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::PropertiesType

◆ ShapeFunctionDerivativesArrayType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef GeometryType::ShapeFunctionsGradientsType Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::ShapeFunctionDerivativesArrayType

Type for an array of shape function gradient matrices.

◆ ShapeFunctionDerivativesType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef Kratos::Matrix Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::ShapeFunctionDerivativesType

Type for a matrix containing the shape function gradients.

◆ ShapeFunctionsType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef Kratos::Vector Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::ShapeFunctionsType

Type for shape function values container.

◆ SizeType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef std::size_t Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::SizeType

◆ VectorType

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
typedef Vector Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::VectorType

Vector type for local contributions to the linear system.

Constructor & Destructor Documentation

◆ BinghamFluid() [1/4]

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::BinghamFluid ( IndexType  NewId = 0)
inlineexplicit

Default constuctor.

Parameters
NewIdIndex number of the new element (optional)

◆ BinghamFluid() [2/4]

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::BinghamFluid ( IndexType  NewId,
const NodesArrayType ThisNodes 
)
inline

Constructor using an array of nodes.

Parameters
NewIdIndex of the new element
ThisNodesAn array containing the nodes of the new element

◆ BinghamFluid() [3/4]

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::BinghamFluid ( IndexType  NewId,
GeometryType::Pointer  pGeometry 
)
inline

Constructor using a geometry object.

Parameters
NewIdIndex of the new element
pGeometryPointer to a geometry object

◆ BinghamFluid() [4/4]

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::BinghamFluid ( IndexType  NewId,
GeometryType::Pointer  pGeometry,
PropertiesType::Pointer  pProperties 
)
inline

Constuctor using geometry and properties.

Parameters
NewIdIndex of the new element
pGeometryPointer to a geometry object
pPropertiesPointer to the element's properties

◆ ~BinghamFluid()

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::~BinghamFluid ( )
inlineoverride

Destructor.

Member Function Documentation

◆ Create() [1/2]

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
Element::Pointer Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::Create ( IndexType  NewId,
GeometryType::Pointer  pGeom,
PropertiesType::Pointer  pProperties 
) const
inlineoverride

Create a new element of this type.

Parameters
NewIdIndex of the new element
pGeomA pointer to the geometry of the new element
pPropertiesPointer to the element's properties

◆ Create() [2/2]

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
Element::Pointer Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::Create ( IndexType  NewId,
NodesArrayType const &  ThisNodes,
PropertiesType::Pointer  pProperties 
) const
inlineoverride

Create a new element of this type.

Returns a pointer to a new element, created using given input

Parameters
NewIdthe ID of the new element
ThisNodesthe nodes of the new element
pPropertiesthe properties assigned to the new element
Returns
a Pointer to the new element

◆ EffectiveViscosity()

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
double Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::EffectiveViscosity ( double  Density,
const TShapeFunctionValues &  rN,
const TShapeFunctionGradients &  rDN_DX,
double  ElemSize,
const ProcessInfo rProcessInfo 
)
inlineoverrideprotected

EffectiveViscosity Calculate the effective viscosity at given integration point using the Bingham constitutive model.

The Bingham model is implemented here using a regularized formulation, as the viscosity should be infite while the stress is smaller than the yield stress. The implemented constitutive relation is

stress(i,j) = ( mu + ( 1 - exp(-m*gamma_dot) )* yield_stress / gamma_dot ) * strain_rate(i,j)

where gamma_dot is the second invariant of the strain rate (2*SijSij)^0.5 and m is a regularization coefficient. For details on the formulation see the phd thesis of A. Larese, "A coupled Eulerian-PFEM model for the simulation of overtopping in rockfill dams" (2012) Universitat Politecnica de Catalunya.

mu is interpolated from the nodal VISCOSITY values, yield_stress and m are read from rProcessInfo values YIELD_STRESS and REGULARIZATION_COEFFICIENT respectively.

Note
This assumes that the underlying fluid element provides EquivalentStrainRate(rDN_DX) and EvaluateInPoint(value,VARIABLE,rN) methods.
Parameters
DensityThe fluid's density at the integration point.
rNNodal shape functions evaluated at the integration points (area coordinates for the point).
rDN_DXShape function derivatives at the integration point.
ElemSizeRepresentative length of the element (unused for Bingham).
rProcessInfoProcessInfo instance passed from the ModelPart, containing additional data
Returns
The effective viscosity, in dynamic units (Pa*s or equivalent).

◆ Info()

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
std::string Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::Info ( ) const
inlineoverride

Turn back information as a string.

◆ KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION()

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION ( BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >  )

◆ PrintData()

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
void Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::PrintData ( std::ostream &  rOStream) const
inlineoverride

Print object's data.

◆ PrintInfo()

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
void Kratos::BinghamFluid< TBaseElement, TShapeFunctionValues, TShapeFunctionGradients >::PrintInfo ( std::ostream &  rOStream) const
inlineoverride

Print information about this object.

Friends And Related Function Documentation

◆ Serializer

template<class TBaseElement , class TShapeFunctionValues = typename TBaseElement::ShapeFunctionsType, class TShapeFunctionGradients = typename TBaseElement::ShapeFunctionDerivativesType>
friend class Serializer
friend

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