![]() |
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.
|
This is a triangular prism solid element for the analysis of thin/thick shells undergoing large elastic–plastic strains. More...
#include <solid_shell_element_sprism_3D6N.h>
Classes | |
class | CartesianDerivatives |
Here the cartesian derivatives are defined. More... | |
class | CommonComponents |
Common Components defined in order to compute the Cauchy tensor and the deformation matrix. More... | |
class | EASComponents |
EAS Components. More... | |
struct | GeneralVariables |
struct | LocalSystemComponents |
class | OrthogonalBase |
OrthogonalBase. More... | |
class | StressIntegratedComponents |
Stress integrated Components used during the integration. More... | |
class | TransverseGradient |
TransverseGradient. More... | |
class | TransverseGradientIsoParametric |
TransverseGradientIsoParametric. More... | |
Public Types | |
Enums | |
enum class | Configuration { INITIAL = 0 , CURRENT = 1 } |
This enum is defined in oder to difereniate between initial (TL) and current (UL) configuration. More... | |
enum class | GeometricLevel { LOWER = 0 , CENTER = 5 , UPPER = 9 } |
To differtiate between center, lower part and upper part. More... | |
enum class | OrthogonalBaseApproach { X = 0 , Y = 1 , Z = 2 } |
To differtiate between the different possible orthogonal bases. More... | |
![]() | |
typedef ConstitutiveLaw | ConstitutiveLawType |
Reference type definition for constitutive laws. More... | |
typedef ConstitutiveLawType::Pointer | ConstitutiveLawPointerType |
Pointer type for constitutive laws. More... | |
typedef ConstitutiveLawType::StressMeasure | StressMeasureType |
StressMeasure from constitutive laws. More... | |
typedef GeometryData::IntegrationMethod | IntegrationMethod |
Type definition for integration methods. More... | |
typedef Node | NodeType |
This is the definition of the node. More... | |
typedef Element | BaseType |
The base element type. More... | |
typedef ConstitutiveLaw | ConstitutiveLawType |
Reference type definition for constitutive laws. More... | |
typedef ConstitutiveLawType::Pointer | ConstitutiveLawPointerType |
Pointer type for constitutive laws. More... | |
typedef ConstitutiveLawType::StressMeasure | StressMeasureType |
StressMeasure from constitutive laws. More... | |
typedef GeometryData::IntegrationMethod | IntegrationMethod |
Type definition for integration methods. More... | |
typedef Node | NodeType |
This is the definition of the node. More... | |
typedef Element | BaseType |
The base element type. More... | |
![]() | |
typedef Element | ElementType |
definition of element type More... | |
typedef GeometricalObject | BaseType |
base type: an GeometricalObject that automatically has a unique number More... | |
typedef Node | NodeType |
definition of node type (default is: Node) More... | |
typedef Properties | PropertiesType |
typedef Geometry< NodeType > | GeometryType |
definition of the geometry type with given NodeType More... | |
typedef Geometry< NodeType >::PointsArrayType | NodesArrayType |
definition of nodes container type, redefined from GeometryType More... | |
typedef Vector | VectorType |
typedef Matrix | MatrixType |
typedef std::size_t | IndexType |
typedef std::size_t | SizeType |
typedef Dof< double > | DofType |
typedef std::vector< std::size_t > | EquationIdVectorType |
typedef std::vector< DofType::Pointer > | DofsVectorType |
typedef PointerVectorSet< DofType > | DofsArrayType |
typedef GeometryData::IntegrationMethod | IntegrationMethod |
Type definition for integration methods. More... | |
typedef GeometryData | GeometryDataType |
![]() | |
typedef Node | NodeType |
Definition of the node type. More... | |
typedef Geometry< NodeType > | GeometryType |
The geometry type definition. More... | |
typedef std::size_t | IndexType |
Defines the index type. More... | |
typedef std::size_t | result_type |
Defines the result type. More... | |
![]() | |
typedef std::size_t | IndexType |
The definition of the index type. More... | |
typedef std::size_t | result_type |
The definition of the result_type. More... | |
![]() | |
enum | FlagsList { Flag0 = BlockType(1) , Flag1 = BlockType(1) << 1 , Flag2 = BlockType(1) << 2 , Flag3 = BlockType(1) << 3 , Flag4 = BlockType(1) << 4 , Flag5 = BlockType(1) << 5 , Flag6 = BlockType(1) << 6 , Flag7 = BlockType(1) << 7 , Flag8 = BlockType(1) << 8 , Flag9 = BlockType(1) << 9 , Flag10 = BlockType(1) << 10 , Flag11 = BlockType(1) << 11 , Flag12 = BlockType(1) << 12 , Flag13 = BlockType(1) << 13 , Flag14 = BlockType(1) << 14 , Flag15 = BlockType(1) << 15 , Flag16 = BlockType(1) << 16 , Flag17 = BlockType(1) << 17 , Flag18 = BlockType(1) << 18 , Flag19 = BlockType(1) << 19 , Flag20 = BlockType(1) << 20 , Flag21 = BlockType(1) << 21 , Flag22 = BlockType(1) << 22 , Flag23 = BlockType(1) << 23 , Flag24 = BlockType(1) << 24 , Flag25 = BlockType(1) << 25 , Flag26 = BlockType(1) << 26 , Flag27 = BlockType(1) << 27 , Flag28 = BlockType(1) << 28 , Flag29 = BlockType(1) << 29 , Flag30 = BlockType(1) << 30 } |
typedef int64_t | BlockType |
typedef int64_t | FlagType |
typedef std::size_t | IndexType |
Public Member Functions | |
Life Cycle | |
SolidShellElementSprism3D6N () | |
SolidShellElementSprism3D6N (IndexType NewId, GeometryType::Pointer pGeometry) | |
SolidShellElementSprism3D6N (IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) | |
SolidShellElementSprism3D6N (SolidShellElementSprism3D6N const &rOther) | |
~SolidShellElementSprism3D6N () override | |
Operators | |
SolidShellElementSprism3D6N & | operator= (SolidShellElementSprism3D6N const &rOther) |
Operations | |
Element::Pointer | Create (IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override |
Creates a new element. More... | |
Element::Pointer | Create (IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override |
Creates a new element. More... | |
Element::Pointer | Clone (IndexType NewId, NodesArrayType const &ThisNodes) const override |
void | EquationIdVector (EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override |
Sets on rResult the ID's of the element degrees of freedom. More... | |
void | GetDofList (DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override |
Sets on rElementalDofList the degrees of freedom of the considered element geometry. More... | |
void | GetValuesVector (Vector &rValues, int Step=0) const override |
Sets on rValues the nodal displacements. More... | |
void | GetFirstDerivativesVector (Vector &rValues, int Step=0) const override |
Sets on rValues the nodal velocities. More... | |
void | GetSecondDerivativesVector (Vector &rValues, int Step=0) const override |
Sets on rValues the nodal accelerations. More... | |
void | CalculateRightHandSide (VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override |
This is called during the assembling process in order to calculate the elemental right hand side vector only. More... | |
void | CalculateLeftHandSide (MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override |
This is called during the assembling process in order to calculate the elemental leTransverseGradientFt hand side vector only. More... | |
void | CalculateLocalSystem (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override |
This is called during the assembling process in order to calculate all elemental contributions to the global system matrix and the right hand side. More... | |
void | CalculateMassMatrix (MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override |
This is called during the assembling process in order to calculate the elemental mass matrix. More... | |
void | CalculateDampingMatrix (MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo) override |
This is called during the assembling process in order to calculate the elemental damping matrix. More... | |
void | CalculateDampingMatrix (MatrixType &rDampingMatrix, const MatrixType &rStiffnessMatrix, const MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) |
This is called during the assembling process in order to calculate the elemental damping matrix (Reusing the stiffness matrix and mass matrix) More... | |
void | CalculateOnIntegrationPoints (const Variable< bool > &rVariable, std::vector< bool > &rOutput, const ProcessInfo &rCurrentProcessInfo) override |
Calculate a boolean Variable on the Element Constitutive Law. More... | |
void | CalculateOnIntegrationPoints (const Variable< int > &rVariable, std::vector< int > &rOutput, const ProcessInfo &rCurrentProcessInfo) override |
Calculate a integer Variable on the Element Constitutive Law. More... | |
void | CalculateOnIntegrationPoints (const Variable< double > &rVariable, std::vector< double > &rOutput, const ProcessInfo &rCurrentProcessInfo) override |
Calculate a double Variable on the Element Constitutive Law. More... | |
void | CalculateOnIntegrationPoints (const Variable< array_1d< double, 3 >> &rVariable, std::vector< array_1d< double, 3 >> &rOutput, const ProcessInfo &rCurrentProcessInfo) override |
Calculate a 3 components array_1d on the Element Constitutive Law. More... | |
void | CalculateOnIntegrationPoints (const Variable< array_1d< double, 6 >> &rVariable, std::vector< array_1d< double, 6 >> &rOutput, const ProcessInfo &rCurrentProcessInfo) override |
Calculate a 6 components array_1d on the Element Constitutive Law. More... | |
void | CalculateOnIntegrationPoints (const Variable< Vector > &rVariable, std::vector< Vector > &rOutput, const ProcessInfo &rCurrentProcessInfo) override |
Calculate a Vector Variable on the Element Constitutive Law. More... | |
void | CalculateOnIntegrationPoints (const Variable< Matrix > &rVariable, std::vector< Matrix > &rOutput, const ProcessInfo &rCurrentProcessInfo) override |
Calculate a Matrix Variable on the Element Constitutive Law. More... | |
void | SetValuesOnIntegrationPoints (const Variable< double > &rVariable, const std::vector< double > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a double Value on the Element Constitutive Law. More... | |
void | SetValuesOnIntegrationPoints (const Variable< Vector > &rVariable, const std::vector< Vector > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a Vector Value on the Element Constitutive Law. More... | |
void | SetValuesOnIntegrationPoints (const Variable< Matrix > &rVariable, const std::vector< Matrix > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a Matrix Value on the Element Constitutive Law. More... | |
void | SetValuesOnIntegrationPoints (const Variable< ConstitutiveLaw::Pointer > &rVariable, const std::vector< ConstitutiveLaw::Pointer > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a Constitutive Law Value. More... | |
int | Check (const ProcessInfo &rCurrentProcessInfo) const override |
This function provides the place to perform checks on the completeness of the input. 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... | |
void | InitializeSolutionStep (const ProcessInfo &rCurrentProcessInfo) override |
Called at the beginning of each solution step. More... | |
void | FinalizeSolutionStep (const ProcessInfo &rCurrentProcessInfo) override |
Called at the end of each solution step. More... | |
void | InitializeNonLinearIteration (const ProcessInfo &rCurrentProcessInfo) override |
This is called for non-linear analysis at the beginning of the iteration process. More... | |
void | FinalizeNonLinearIteration (const ProcessInfo &rCurrentProcessInfo) override |
This is called for non-linear analysis at the beginning of the iteration process. More... | |
void | Initialize (const ProcessInfo &rCurrentProcessInfo) override |
Called to initialize the element. More... | |
![]() | |
BaseSolidElement () | |
BaseSolidElement (IndexType NewId, GeometryType::Pointer pGeometry) | |
BaseSolidElement (IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) | |
BaseSolidElement (BaseSolidElement const &rOther) | |
~BaseSolidElement () override | |
BaseSolidElement () | |
BaseSolidElement (IndexType NewId, GeometryType::Pointer pGeometry) | |
BaseSolidElement (IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) | |
BaseSolidElement (BaseSolidElement const &rOther) | |
~BaseSolidElement () override | |
void | ResetConstitutiveLaw () override |
This resets the constitutive law. More... | |
IntegrationMethod | GetIntegrationMethod () const override |
Returns the used integration method. More... | |
void | AddExplicitContribution (const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< double > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo) override |
This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable rDestinationVariable (double version) More... | |
void | AddExplicitContribution (const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< array_1d< double, 3 > > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo) override |
This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable (array_1d<double, 3>) version rDestinationVariable. More... | |
void | SetValuesOnIntegrationPoints (const Variable< bool > &rVariable, const std::vector< bool > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a bool Value on the Element Constitutive Law. More... | |
void | SetValuesOnIntegrationPoints (const Variable< int > &rVariable, const std::vector< int > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a int Value on the Element Constitutive Law. More... | |
void | CalculateOnIntegrationPoints (const Variable< ConstitutiveLaw::Pointer > &rVariable, std::vector< ConstitutiveLaw::Pointer > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Get on rVariable Constitutive Law from the element. More... | |
void | SetValuesOnIntegrationPoints (const Variable< array_1d< double, 3 > > &rVariable, const std::vector< array_1d< double, 3 > > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a array of 3 compoenents Value on the Element. More... | |
void | SetValuesOnIntegrationPoints (const Variable< array_1d< double, 6 > > &rVariable, const std::vector< array_1d< double, 6 > > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a array of 6 compoenents Value on the Element. More... | |
void | CalculateRayleighDampingMatrix (Element &rElement, Element::MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo, const std::size_t MatrixSize) |
void | ResetConstitutiveLaw () override |
This resets the constitutive law. More... | |
IntegrationMethod | GetIntegrationMethod () const override |
Returns the used integration method. More... | |
virtual bool | UseGeometryIntegrationMethod () const |
virtual const GeometryType::IntegrationPointsArrayType | IntegrationPoints () const |
virtual const GeometryType::IntegrationPointsArrayType | IntegrationPoints (IntegrationMethod ThisMethod) const |
virtual const Matrix & | ShapeFunctionsValues (IntegrationMethod ThisMethod) const |
void | AddExplicitContribution (const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< double > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo) override |
This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable rDestinationVariable (double version) More... | |
void | AddExplicitContribution (const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< array_1d< double, 3 > > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo) override |
This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable (array_1d<double, 3>) version rDestinationVariable. More... | |
void | CalculateOnIntegrationPoints (const Variable< ConstitutiveLaw::Pointer > &rVariable, std::vector< ConstitutiveLaw::Pointer > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Get on rVariable Constitutive Law from the element. More... | |
void | SetValuesOnIntegrationPoints (const Variable< bool > &rVariable, const std::vector< bool > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a bool Value on the Element Constitutive Law. More... | |
void | SetValuesOnIntegrationPoints (const Variable< int > &rVariable, const std::vector< int > &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a int Value on the Element Constitutive Law. More... | |
void | SetValuesOnIntegrationPoints (const Variable< array_1d< double, 3 >> &rVariable, const std::vector< array_1d< double, 3 >> &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a array of 3 compoenents Value on the Element. More... | |
void | SetValuesOnIntegrationPoints (const Variable< array_1d< double, 6 >> &rVariable, const std::vector< array_1d< double, 6 >> &rValues, const ProcessInfo &rCurrentProcessInfo) override |
Set a array of 6 compoenents Value on the Element. More... | |
const Parameters | GetSpecifications () const override |
This method provides the specifications/requirements of the element. More... | |
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (BaseSolidElement) | |
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (BaseSolidElement) | |
KRATOS_DEFINE_LOCAL_FLAG (ROTATED) | |
![]() | |
Element (IndexType NewId=0) | |
Element (IndexType NewId, const NodesArrayType &ThisNodes) | |
Element (IndexType NewId, GeometryType::Pointer pGeometry) | |
Element (IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) | |
Element (Element const &rOther) | |
Copy constructor. More... | |
~Element () override | |
Destructor. More... | |
Element & | operator= (Element const &rOther) |
Assignment operator. More... | |
virtual void | CalculateFirstDerivativesContributions (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateFirstDerivativesLHS (MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateFirstDerivativesRHS (VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateSecondDerivativesContributions (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateSecondDerivativesLHS (MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateSecondDerivativesRHS (VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) |
virtual void | AddExplicitContribution (const ProcessInfo &rCurrentProcessInfo) |
virtual void | AddExplicitContribution (const MatrixType &rLHSMatrix, const Variable< MatrixType > &rLHSVariable, const Variable< Matrix > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo) |
This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable rDestinationVariable. (This is the matrix version) More... | |
virtual void | Calculate (const Variable< double > &rVariable, double &Output, const ProcessInfo &rCurrentProcessInfo) |
virtual void | Calculate (const Variable< array_1d< double, 3 > > &rVariable, array_1d< double, 3 > &Output, const ProcessInfo &rCurrentProcessInfo) |
virtual void | Calculate (const Variable< Vector > &rVariable, Vector &Output, const ProcessInfo &rCurrentProcessInfo) |
virtual void | Calculate (const Variable< Matrix > &rVariable, Matrix &Output, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateOnIntegrationPoints (const Variable< array_1d< double, 4 >> &rVariable, std::vector< array_1d< double, 4 >> &rOutput, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateOnIntegrationPoints (const Variable< array_1d< double, 9 >> &rVariable, std::vector< array_1d< double, 9 >> &rOutput, const ProcessInfo &rCurrentProcessInfo) |
virtual void | SetValuesOnIntegrationPoints (const Variable< array_1d< double, 4 >> &rVariable, const std::vector< array_1d< double, 4 >> &rValues, const ProcessInfo &rCurrentProcessInfo) |
virtual void | SetValuesOnIntegrationPoints (const Variable< array_1d< double, 9 >> &rVariable, const std::vector< array_1d< double, 9 >> &rValues, const ProcessInfo &rCurrentProcessInfo) |
virtual void | MassMatrix (MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) |
virtual void | AddMassMatrix (MatrixType &rLeftHandSideMatrix, double coeff, const ProcessInfo &rCurrentProcessInfo) |
virtual void | DampMatrix (MatrixType &rDampMatrix, const ProcessInfo &rCurrentProcessInfo) |
virtual void | AddInertiaForces (VectorType &rRightHandSideVector, double coeff, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateLocalVelocityContribution (MatrixType &rDampingMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateSensitivityMatrix (const Variable< double > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo) |
virtual void | CalculateSensitivityMatrix (const Variable< array_1d< double, 3 > > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo) |
PropertiesType::Pointer | pGetProperties () |
returns the pointer to the property of the element. Does not throw an error, to allow copying of elements which don't have any property assigned. More... | |
const PropertiesType::Pointer | pGetProperties () const |
PropertiesType & | GetProperties () |
PropertiesType const & | GetProperties () const |
void | SetProperties (PropertiesType::Pointer pProperties) |
bool | HasProperties () const |
Check that the Element has a correctly initialized pointer to a Properties instance. More... | |
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (Element) | |
![]() | |
GeometricalObject (IndexType NewId=0) | |
Default constructor. More... | |
GeometricalObject (IndexType NewId, GeometryType::Pointer pGeometry) | |
Default constructor. More... | |
~GeometricalObject () override | |
Destructor. More... | |
GeometricalObject (GeometricalObject const &rOther) | |
Copy constructor. More... | |
GeometricalObject & | operator= (GeometricalObject const &rOther) |
Assignment operator. More... | |
virtual void | SetGeometry (GeometryType::Pointer pGeometry) |
Sets the pointer to the geometry. More... | |
GeometryType::Pointer | pGetGeometry () |
Returns the pointer to the geometry. More... | |
const GeometryType::Pointer | pGetGeometry () const |
Returns the pointer to the geometry (const version) More... | |
GeometryType & | GetGeometry () |
Returns the reference of the geometry. More... | |
GeometryType const & | GetGeometry () const |
Returns the reference of the geometry (const version) More... | |
Flags & | GetFlags () |
Returns the flags of the object. More... | |
Flags const & | GetFlags () const |
Returns the flags of the object (const version) More... | |
void | SetFlags (Flags const &rThisFlags) |
Sets the flags of the object. More... | |
DataValueContainer & | Data () |
DataValueContainer & | GetData () |
DataValueContainer const & | GetData () const |
void | SetData (DataValueContainer const &rThisData) |
template<class TDataType > | |
bool | Has (const Variable< TDataType > &rThisVariable) const |
template<class TVariableType > | |
void | SetValue (const TVariableType &rThisVariable, typename TVariableType::Type const &rValue) |
template<class TVariableType > | |
TVariableType::Type & | GetValue (const TVariableType &rThisVariable) |
template<class TVariableType > | |
TVariableType::Type const & | GetValue (const TVariableType &rThisVariable) const |
unsigned int | use_count () const noexcept |
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (GeometricalObject) | |
Pointer definition of GeometricalObject. More... | |
bool | IsActive () const |
Checks if the GeometricalObject is active. More... | |
![]() | |
IndexedObject (IndexType NewId=0) | |
Default constructor. More... | |
virtual | ~IndexedObject () |
Destructor. More... | |
IndexedObject (IndexedObject const &rOther) | |
Copy constructor. More... | |
IndexedObject & | operator= (IndexedObject const &rOther) |
Assignment operator. More... | |
template<class TObjectType > | |
IndexType | operator() (TObjectType const &rThisObject) const |
IndexType | Id () const |
IndexType | GetId () const |
virtual void | SetId (IndexType NewId) |
IndexType & | DepricatedIdAccess () |
TODO: remove this function when removing data_file_io object. More... | |
KRATOS_CLASS_POINTER_DEFINITION (IndexedObject) | |
Pointer definition of IndexedObject. More... | |
![]() | |
Flags & | operator= (Flags const &rOther) |
Assignment operator. More... | |
operator bool () const | |
Flags | operator~ () const |
bool | operator! () const |
void | AssignFlags (Flags const &rOther) |
void | Set (const Flags ThisFlag) |
void | Set (const Flags ThisFlag, bool Value) |
void | Reset (const Flags ThisFlag) |
void | Flip (const Flags ThisFlag) |
void | SetPosition (IndexType Position, bool Value=true) |
bool | GetPosition (IndexType Position) const |
void | FlipPosition (IndexType Position) |
void | ClearPosition (IndexType Position) |
void | Clear () |
Flags | AsFalse () const |
bool | Is (Flags const &rOther) const |
bool | IsDefined (Flags const &rOther) const |
bool | IsNot (Flags const &rOther) const |
bool | IsNotDefined (Flags const &rOther) const |
KRATOS_CLASS_POINTER_DEFINITION (Flags) | |
Pointer definition of Flags. More... | |
const Flags & | operator|= (const Flags &Other) |
const Flags & | operator&= (const Flags &Other) |
Flags () | |
Default constructor. More... | |
Flags (Flags const &rOther) | |
Copy constructor. More... | |
virtual | ~Flags () |
Destructor. More... | |
Protected Member Functions | |
Protected Operators | |
void | CalculateElementalSystem (LocalSystemComponents &rLocalSystem, const ProcessInfo &rCurrentProcessInfo) |
Calculates the elemental contributions. More... | |
void | PrintElementCalculation (LocalSystemComponents &rLocalSystem, GeneralVariables &rVariables) |
Prints element information for each gauss point (debugging purposes) More... | |
Protected Operations | |
bool | HasNeighbour (const IndexType Index, const NodeType &NeighbourNode) const |
Check if the node has a neighbour: More... | |
std::size_t | NumberOfActiveNeighbours (const GlobalPointersVector< NodeType > &pNeighbourNodes) const |
Calculates the number of active neighbours: More... | |
void | GetNodalCoordinates (BoundedMatrix< double, 12, 3 > &NodesCoord, const GlobalPointersVector< NodeType > &pNeighbourNodes, const Configuration ThisConfiguration) const |
It gets the nodal coordinates, according to the configutaion. More... | |
void | CalculateCartesianDerivatives (CartesianDerivatives &rCartesianDerivatives) |
Calculate the cartesian derivatives. More... | |
void | CalculateCommonComponents (CommonComponents &rCommonComponents, const CartesianDerivatives &rCartesianDerivatives) |
Calculate the necessary components for the Kinematic calculus. More... | |
void | CalculateLocalCoordinateSystem (OrthogonalBase &ThisOrthogonalBase, const OrthogonalBaseApproach ThisOrthogonalBaseApproach, const double ThisAngle) |
Calculates the Local Coordinates System. More... | |
void | CalculateIdVector (array_1d< IndexType, 18 > &rIdVector) |
Calculate the vector of the element Ids. More... | |
void | ComputeLocalDerivatives (BoundedMatrix< double, 6, 3 > &LocalDerivativePatch, const array_1d< double, 3 > &rLocalCoordinates) |
Calculate the local derivatives of the element for a given coordinates. More... | |
void | ComputeLocalDerivativesQuadratic (BoundedMatrix< double, 4, 2 > &rLocalDerivativePatch, const IndexType NodeGauss) |
Calculate the local quadratic derivatives of the element for a given gauss node. More... | |
void | CalculateJacobianCenterGauss (GeometryType::JacobiansType &J, std::vector< Matrix > &Jinv, Vector &detJ, const IndexType rPointNumber, const double ZetaGauss) |
Calculate the Jacobian and his inverse. More... | |
void | CalculateJacobian (double &detJ, BoundedMatrix< double, 3, 3 > &J, BoundedMatrix< double, 6, 3 > &LocalDerivativePatch, const BoundedMatrix< double, 12, 3 > &NodesCoord, const array_1d< double, 3 > &rLocalCoordinates) |
Calculate the Jacobian. More... | |
void | CalculateJacobianAndInv (BoundedMatrix< double, 3, 3 > &J, BoundedMatrix< double, 3, 3 > &Jinv, BoundedMatrix< double, 6, 3 > &LocalDerivativePatch, const BoundedMatrix< double, 3, 6 > &NodesCoord, const array_1d< double, 3 > &rLocalCoordinates) |
Calculate the Jacobian and its inverse. More... | |
void | CalculateJacobianAndInv (BoundedMatrix< double, 3, 3 > &J, BoundedMatrix< double, 3, 3 > &Jinv, const BoundedMatrix< double, 3, 6 > &NodesCoord, const array_1d< double, 3 > &rLocalCoordinates) |
Calculate the Jacobian and his inverse. More... | |
void | CalculateCartesianDerivativesOnCenterPlane (BoundedMatrix< double, 2, 4 > &CartesianDerivativesCenter, const OrthogonalBase &ThisOrthogonalBase, const GeometricLevel Part) |
Calculate the Cartesian derivatives in the Gauss points, for the plane. More... | |
void | CalculateCartesianDerOnGaussPlane (BoundedMatrix< double, 2, 4 > &InPlaneCartesianDerivativesGauss, const BoundedMatrix< double, 12, 3 > &NodesCoord, const OrthogonalBase &ThisOrthogonalBase, const IndexType NodeGauss, const GeometricLevel Part) |
Calculate the Cartesian derivatives in the Gauss points, for the plane. More... | |
void | CalculateCartesianDerOnGaussTrans (BoundedMatrix< double, 6, 1 > &TransversalCartesianDerivativesGauss, const BoundedMatrix< double, 12, 3 > &NodesCoord, const OrthogonalBase &ThisOrthogonalBase, const array_1d< double, 3 > &rLocalCoordinates) |
Calculate the Cartesian derivatives in the Gauss points, for the transversal direction. More... | |
void | CalculateCartesianDerOnCenterTrans (CartesianDerivatives &rCartesianDerivatives, const BoundedMatrix< double, 12, 3 > &NodesCoord, const OrthogonalBase &ThisOrthogonalBase, const GeometricLevel Part) |
Calculate the Cartesian derivatives in the center, for the transversal direction. More... | |
void | CalculateInPlaneGradientFGauss (BoundedMatrix< double, 3, 2 > &InPlaneGradientFGauss, const BoundedMatrix< double, 2, 4 > &InPlaneCartesianDerivativesGauss, const BoundedMatrix< double, 12, 3 > &NodesCoord, const IndexType NodeGauss, const GeometricLevel Part) |
void | CalculateTransverseGradientF (array_1d< double, 3 > &TransverseGradientF, const BoundedMatrix< double, 6, 1 > &TransversalCartesianDerivativesGauss, const BoundedMatrix< double, 12, 3 > &NodesCoord) |
void | CalculateTransverseGradientFinP (TransverseGradientIsoParametric &rTransverseGradientIsoParametric, const BoundedMatrix< double, 12, 3 > &NodesCoord, const GeometricLevel Part) |
Calculate the transversal components of the deformation gradient, in each one of the faces: More... | |
void | CalculateAndAddBMembrane (BoundedMatrix< double, 3, 18 > &BMembrane, BoundedMatrix< double, 3, 1 > &CMembrane, const BoundedMatrix< double, 2, 4 > &InPlaneCartesianDerivativesGauss, const BoundedMatrix< double, 3, 2 > &InPlaneGradientFGauss, const IndexType NodeGauss) |
Construction of the membrane deformation tangent matrix: More... | |
void | CalculateAndAddMembraneKgeometric (BoundedMatrix< double, 36, 36 > &Kgeometricmembrane, const CartesianDerivatives &rCartesianDerivatives, const array_1d< double, 3 > &SMembrane, const GeometricLevel Part) |
Construction of the in-plane geometric stiffness matrix: More... | |
void | CalculateAndAddBShear (BoundedMatrix< double, 2, 18 > &BShear, BoundedMatrix< double, 2, 1 > &CShear, const CartesianDerivatives &rCartesianDerivatives, const TransverseGradient &rTransverseGradient, const TransverseGradientIsoParametric &rTransverseGradientIsoParametric, const GeometricLevel Part) |
Construction of the shear deformation tangent matrix: More... | |
void | CalculateAndAddShearKgeometric (BoundedMatrix< double, 36, 36 > &Kgeometricshear, const CartesianDerivatives &rCartesianDerivatives, const array_1d< double, 2 > &SShear, const GeometricLevel Part) |
Construction of the shear geometric contribution to the stiffness matrix: More... | |
void | CalculateAndAddBNormal (BoundedMatrix< double, 1, 18 > &BNormal, double &CNormal, const BoundedMatrix< double, 6, 1 > &TransversalCartesianDerivativesGaussCenter, const array_1d< double, 3 > &TransversalDeformationGradientF) |
Construction of the transversal deformation tangent matrix: More... | |
void | CalculateAndAddNormalKgeometric (BoundedMatrix< double, 36, 36 > &Kgeometricnormal, const BoundedMatrix< double, 6, 1 > &TransversalCartesianDerivativesGaussCenter, const double SNormal) |
Construction of the transversal geometric contribution to the stiffness matrix: More... | |
BoundedMatrix< double, 36, 1 > | GetVectorCurrentPosition () |
Calculates the vector of current position. More... | |
BoundedMatrix< double, 36, 1 > | GetVectorPreviousPosition () |
void | IntegrateStressesInZeta (GeneralVariables &rVariables, StressIntegratedComponents &rIntegratedStress, const double AlphaEAS, const double ZetaGauss, const double IntegrationWeight) |
Integrates stresses in zeta using the Gauss Quadrature. More... | |
void | IntegrateEASInZeta (GeneralVariables &rVariables, EASComponents &rEAS, const double ZetaGauss, const double IntegrationWeight) |
Integrates the EAS components in zeta using the Gauss Quadrature. More... | |
void | CalculateAndAddLHS (LocalSystemComponents &rLocalSystem, GeneralVariables &rVariables, ConstitutiveLaw::Parameters &rValues, const StressIntegratedComponents &rIntegratedStress, const CommonComponents &rCommonComponents, const CartesianDerivatives &rCartesianDerivatives, const EASComponents &rEAS, double &AlphaEAS) |
Calculation and addition of the matrix of the LHS. More... | |
void | CalculateAndAddRHS (LocalSystemComponents &rLocalSystem, GeneralVariables &rVariables, Vector &rVolumeForce, const StressIntegratedComponents &rIntegratedStress, const CommonComponents &rCommonComponents, const EASComponents &rEAS, double &AlphaEAS) |
Calculation and addition of the vectors of the RHS. More... | |
void | CalculateAndAddKuum (MatrixType &rLeftHandSideMatrix, GeneralVariables &rVariables, const double IntegrationWeight) |
Calculation of the Material Stiffness Matrix. Kuum = BT * C * B. More... | |
void | CalculateAndAddKuug (MatrixType &rLeftHandSideMatrix, const StressIntegratedComponents &rIntegratedStress, const CartesianDerivatives &rCartesianDerivatives) |
Calculation of the Geometric Stiffness Matrix. Kuug = BT * S. More... | |
void | ApplyEASLHS (MatrixType &rLeftHandSideMatrix, const EASComponents &rEAS) |
Update the LHS of the system with the EAS. More... | |
void | ApplyEASRHS (BoundedMatrix< double, 36, 1 > &rRHSFull, const EASComponents &rEAS, double &AlphaEAS) |
Update the RHS of the system with the EAS and the internal variable alpha. More... | |
void | CalculateAndAddExternalForces (VectorType &rRightHandSideVector, GeneralVariables &rVariables, Vector &rVolumeForce) |
Calculation of the External Forces Vector. Fe = N * t + N * b. More... | |
void | CalculateAndAddInternalForces (VectorType &rRightHandSideVector, const StressIntegratedComponents &rIntegratedStress, const CommonComponents &rCommonComponents, const EASComponents &rEAS, double &AlphaEAS) |
Calculation of the Internal Forces Vector. Fi = B * sigma. More... | |
void | SetGeneralVariables (GeneralVariables &rVariables, ConstitutiveLaw::Parameters &rValues, const IndexType rPointNumber) |
Set Variables of the Element to the Parameters of the Constitutive Law. More... | |
void | InitializeSystemMatrices (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, Flags &rCalculationFlags) |
Initialize System Matrices. More... | |
void | CalculateDeltaPosition (Matrix &rDeltaPosition) |
This method computes the delta position matrix necessary for UL formulation. More... | |
void | CalculateKinematics (GeneralVariables &rVariables, const CommonComponents &rCommonComponents, const GeometryType::IntegrationPointsArrayType &rIntegrationPoints, const IndexType rPointNumber, const double AlphaEAS, const double ZetaGauss) |
Calculate Element Kinematics. More... | |
void | CbartoFbar (GeneralVariables &rVariables, const int rPointNumber) |
Calculate Fbar from Cbar. More... | |
void | CalculateDeformationMatrix (Matrix &rB, const CommonComponents &rCommonComponents, const double ZetaGauss, const double AlphaEAS) |
Calculation of the Deformation Matrix BL. More... | |
void | InitializeGeneralVariables (GeneralVariables &rVariables) |
Initialize Element General Variables. More... | |
void | FinalizeStepVariables (GeneralVariables &rVariables, const IndexType rPointNumber) |
Finalize Element Internal Variables. More... | |
void | GetHistoricalVariables (GeneralVariables &rVariables, const IndexType rPointNumber) |
Get the Historical Deformation Gradient to calculate aTransverseGradientFter finalize the step. More... | |
void | CalculateVolumeChange (double &rVolumeChange, GeneralVariables &rVariables) |
This function calculates the variation of the element volume. More... | |
void | CalculateVolumeForce (Vector &rVolumeForce, GeneralVariables &rVariables, const double IntegrationWeight) |
Calculation of the Volume Force of the Element. More... | |
![]() | |
void | SetIntegrationMethod (const IntegrationMethod &ThisIntegrationMethod) |
Sets the used integration method. More... | |
void | SetConstitutiveLawVector (const std::vector< ConstitutiveLaw::Pointer > &ThisConstitutiveLawVector) |
Sets the used constitutive laws. More... | |
virtual void | InitializeMaterial () |
It initializes the material. More... | |
virtual ConstitutiveLaw::StressMeasure | GetStressMeasure () const |
Gives the StressMeasure used. More... | |
virtual bool | UseElementProvidedStrain () const |
This method returns if the element provides the strain. More... | |
virtual void | CalculateAll (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo, const bool CalculateStiffnessMatrixFlag, const bool CalculateResidualVectorFlag) |
This functions calculates both the RHS and the LHS. More... | |
virtual void | CalculateKinematicVariables (KinematicVariables &rThisKinematicVariables, const IndexType PointNumber, const GeometryType::IntegrationMethod &rIntegrationMethod) |
This functions updates the kinematics variables. More... | |
virtual void | SetConstitutiveVariables (KinematicVariables &rThisKinematicVariables, ConstitutiveVariables &rThisConstitutiveVariables, ConstitutiveLaw::Parameters &rValues, const IndexType PointNumber, const GeometryType::IntegrationPointsArrayType &IntegrationPoints) |
This functions updates the data structure passed to the CL. More... | |
virtual void | CalculateConstitutiveVariables (KinematicVariables &rThisKinematicVariables, ConstitutiveVariables &rThisConstitutiveVariables, ConstitutiveLaw::Parameters &rValues, const IndexType PointNumber, const GeometryType::IntegrationPointsArrayType &IntegrationPoints, const ConstitutiveLaw::StressMeasure ThisStressMeasure=ConstitutiveLaw::StressMeasure_PK2) |
This functions updates the constitutive variables. More... | |
Matrix & | CalculateDeltaDisplacement (Matrix &DeltaDisplacement) const |
This methods gives us a matrix with the increment of displacement. More... | |
virtual double | CalculateDerivativesOnReferenceConfiguration (Matrix &rJ0, Matrix &rInvJ0, Matrix &rDN_DX, const IndexType PointNumber, IntegrationMethod ThisIntegrationMethod) const |
This functions calculate the derivatives in the reference frame. More... | |
double | CalculateDerivativesOnCurrentConfiguration (Matrix &rJ, Matrix &rInvJ, Matrix &rDN_DX, const IndexType PointNumber, IntegrationMethod ThisIntegrationMethod) const |
This functions calculate the derivatives in the current frame. More... | |
virtual array_1d< double, 3 > | GetBodyForce (const GeometryType::IntegrationPointsArrayType &IntegrationPoints, const IndexType PointNumber) const |
This function computes the body force. More... | |
virtual void | CalculateAndAddKm (MatrixType &rLeftHandSideMatrix, const Matrix &B, const Matrix &D, const double IntegrationWeight) const |
Calculation of the Material Stiffness Matrix. Km = B^T * D *B. More... | |
void | CalculateAndAddKg (MatrixType &rLeftHandSideMatrix, const Matrix &DN_DX, const Vector &rStressVector, const double IntegrationWeight) const |
Calculation of the Geometric Stiffness Matrix. Kg = dB * S. More... | |
virtual void | CalculateAndAddResidualVector (VectorType &rRightHandSideVector, const KinematicVariables &rThisKinematicVariables, const ProcessInfo &rCurrentProcessInfo, const array_1d< double, 3 > &rBodyForce, const Vector &rStressVector, const double IntegrationWeight) const |
Calculation of the RHS. More... | |
void | CalculateAndAddExtForceContribution (const Vector &rN, const ProcessInfo &rCurrentProcessInfo, const array_1d< double, 3 > &rBodyForce, VectorType &rRightHandSideVector, const double IntegrationWeight) const |
This function add the external force contribution. More... | |
virtual double | GetIntegrationWeight (const GeometryType::IntegrationPointsArrayType &rThisIntegrationPoints, const IndexType PointNumber, const double detJ) const |
This functions computes the integration weight to consider. More... | |
void | CalculateShapeGradientOfMassMatrix (MatrixType &rMassMatrix, ShapeParameter Deriv) const |
This function computes the shape gradient of mass matrix. More... | |
double | GetRayleighAlpha (const Properties &rProperties, const ProcessInfo &rCurrentProcessInfo) |
double | GetRayleighBeta (const Properties &rProperties, const ProcessInfo &rCurrentProcessInfo) |
void | SetIntegrationMethod (const IntegrationMethod &ThisIntegrationMethod) |
Sets the used integration method. More... | |
void | SetConstitutiveLawVector (const std::vector< ConstitutiveLaw::Pointer > &ThisConstitutiveLawVector) |
Sets the used constitutive laws. More... | |
virtual void | InitializeMaterial () |
It initializes the material. More... | |
virtual ConstitutiveLaw::StressMeasure | GetStressMeasure () const |
Gives the StressMeasure used. More... | |
virtual bool | UseElementProvidedStrain () const |
This method returns if the element provides the strain. More... | |
virtual void | CalculateAll (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo, const bool CalculateStiffnessMatrixFlag, const bool CalculateResidualVectorFlag) |
This functions calculates both the RHS and the LHS. More... | |
virtual void | CalculateKinematicVariables (KinematicVariables &rThisKinematicVariables, const IndexType PointNumber, const GeometryType::IntegrationMethod &rIntegrationMethod) |
This functions updates the kinematics variables. More... | |
virtual void | SetConstitutiveVariables (KinematicVariables &rThisKinematicVariables, ConstitutiveVariables &rThisConstitutiveVariables, ConstitutiveLaw::Parameters &rValues, const IndexType PointNumber, const GeometryType::IntegrationPointsArrayType &IntegrationPoints) |
This functions updates the data structure passed to the CL. More... | |
virtual void | CalculateConstitutiveVariables (KinematicVariables &rThisKinematicVariables, ConstitutiveVariables &rThisConstitutiveVariables, ConstitutiveLaw::Parameters &rValues, const IndexType PointNumber, const GeometryType::IntegrationPointsArrayType &IntegrationPoints, const ConstitutiveLaw::StressMeasure ThisStressMeasure=ConstitutiveLaw::StressMeasure_PK2, const bool IsElementRotated=true) |
This functions updates the constitutive variables. More... | |
Matrix & | CalculateDeltaDisplacement (Matrix &DeltaDisplacement) const |
This methods gives us a matrix with the increment of displacement. More... | |
virtual double | CalculateDerivativesOnReferenceConfiguration (Matrix &rJ0, Matrix &rInvJ0, Matrix &rDN_DX, const IndexType PointNumber, IntegrationMethod ThisIntegrationMethod) const |
This functions calculate the derivatives in the reference frame. More... | |
double | CalculateDerivativesOnCurrentConfiguration (Matrix &rJ, Matrix &rInvJ, Matrix &rDN_DX, const IndexType PointNumber, IntegrationMethod ThisIntegrationMethod) const |
This functions calculate the derivatives in the current frame. More... | |
virtual array_1d< double, 3 > | GetBodyForce (const GeometryType::IntegrationPointsArrayType &IntegrationPoints, const IndexType PointNumber) const |
This function computes the body force. More... | |
virtual void | CalculateAndAddKm (MatrixType &rLeftHandSideMatrix, const Matrix &B, const Matrix &D, const double IntegrationWeight) const |
Calculation of the Material Stiffness Matrix. Km = B^T * D *B. More... | |
void | CalculateAndAddKg (MatrixType &rLeftHandSideMatrix, const Matrix &DN_DX, const Vector &rStressVector, const double IntegrationWeight) const |
Calculation of the Geometric Stiffness Matrix. Kg = dB * S. More... | |
virtual void | CalculateAndAddResidualVector (VectorType &rRightHandSideVector, const KinematicVariables &rThisKinematicVariables, const ProcessInfo &rCurrentProcessInfo, const array_1d< double, 3 > &rBodyForce, const Vector &rStressVector, const double IntegrationWeight) const |
Calculation of the RHS. More... | |
void | CalculateAndAddExtForceContribution (const Vector &rN, const ProcessInfo &rCurrentProcessInfo, const array_1d< double, 3 > &rBodyForce, VectorType &rRightHandSideVector, const double IntegrationWeight) const |
This function add the external force contribution. More... | |
virtual double | GetIntegrationWeight (const GeometryType::IntegrationPointsArrayType &rThisIntegrationPoints, const IndexType PointNumber, const double detJ) const |
This functions computes the integration weight to consider. More... | |
void | CalculateShapeGradientOfMassMatrix (MatrixType &rMassMatrix, ShapeParameter Deriv) const |
This function computes the shape gradient of mass matrix. More... | |
virtual bool | IsElementRotated () const |
This method checks is an element has to be rotated according to a set of local axes. More... | |
void | RotateToLocalAxes (ConstitutiveLaw::Parameters &rValues, KinematicVariables &rThisKinematicVariables) |
This method rotates the F or strain according to local axis from global to local coordinates. More... | |
void | RotateToGlobalAxes (ConstitutiveLaw::Parameters &rValues, KinematicVariables &rThisKinematicVariables) |
This method rotates the F or strain according to local axis from local de global. More... | |
Protected Attributes | |
Protected member Variables | |
bool | mFinalizedStep |
std::vector< Matrix > | mAuxContainer |
Flags | mELementalFlags |
![]() | |
IntegrationMethod | mThisIntegrationMethod |
std::vector< ConstitutiveLaw::Pointer > | mConstitutiveLawVector |
Currently selected integration methods. More... | |
Type Definitions | |
typedef ConstitutiveLaw | ConstitutiveLawType |
Reference type definition for constitutive laws. More... | |
typedef ConstitutiveLawType::Pointer | ConstitutiveLawPointerType |
Pointer type for constitutive laws. More... | |
typedef ConstitutiveLawType::StressMeasure | StressMeasureType |
StressMeasure from constitutive laws. More... | |
typedef GeometryData::IntegrationMethod | IntegrationMethod |
Type definition for integration methods. More... | |
typedef Node | NodeType |
This is the definition of the node. More... | |
typedef BaseSolidElement | BaseType |
The base element type. More... | |
typedef std::size_t | IndexType |
The definition of the index type. More... | |
typedef std::size_t | SizeType |
The definition of the sizetype. More... | |
typedef GlobalPointersVector< NodeType > | WeakPointerVectorNodesType |
KRATOS_DEFINE_LOCAL_FLAG (COMPUTE_RHS_VECTOR) | |
Flags related to the element computation. More... | |
KRATOS_DEFINE_LOCAL_FLAG (COMPUTE_LHS_MATRIX) | |
KRATOS_DEFINE_LOCAL_FLAG (COMPUTE_RHS_VECTOR_WITH_COMPONENTS) | |
KRATOS_DEFINE_LOCAL_FLAG (COMPUTE_LHS_MATRIX_WITH_COMPONENTS) | |
KRATOS_DEFINE_LOCAL_FLAG (EAS_IMPLICIT_EXPLICIT) | |
KRATOS_DEFINE_LOCAL_FLAG (TOTAL_UPDATED_LAGRANGIAN) | |
KRATOS_DEFINE_LOCAL_FLAG (QUADRATIC_ELEMENT) | |
KRATOS_DEFINE_LOCAL_FLAG (EXPLICIT_RHS_COMPUTATION) | |
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (SolidShellElementSprism3D6N) | |
Counted pointer of SolidShellElementSprism3D6N. More... | |
Serialization | |
class | Serializer |
Additional Inherited Members | |
![]() | |
static bool | HasSameType (const GeometricalObject &rLHS, const GeometricalObject &rRHS) |
Checks if two GeometricalObject have the same type. More... | |
static bool | HasSameType (const GeometricalObject *rLHS, const GeometricalObject *rRHS) |
Checks if two GeometricalObject have the same type (pointer version) More... | |
static bool | HasSameGeometryType (const GeometricalObject &rLHS, const GeometricalObject &rRHS) |
Checks if two GeometricalObject have the same geometry type. More... | |
static bool | HasSameGeometryType (const GeometricalObject *rLHS, const GeometricalObject *rRHS) |
Checks if two GeometricalObject have the same geometry type (pointer version) More... | |
static bool | IsSame (const GeometricalObject &rLHS, const GeometricalObject &rRHS) |
Checks if two GeometricalObject are the same. More... | |
static bool | IsSame (const GeometricalObject *rLHS, const GeometricalObject *rRHS) |
Checks if two GeometricalObject are the same (pointer version) More... | |
![]() | |
static const Flags | AllDefined () |
static const Flags | AllTrue () |
static Flags | Create (IndexType ThisPosition, bool Value=true) |
This is a triangular prism solid element for the analysis of thin/thick shells undergoing large elastic–plastic strains.
The element is based on a total Lagrangian formulation and uses as strain measure the logarithm of the right stretch tensor (U) obtained from a modified right Cauchy–Green deformation tensor (C). Three are the introduced modifications: (a) a classical assumed strain approach for transverse shear strains (b) an assumed strain approach for the in-plane components using information from neighbor elements and (c) an averaging of the volumetric strain over the element. The objective is to use this type of elements for the simulation of shells avoiding transverse shear locking, improving the membrane behavior of the in-plane triangle and to handle quasi-incompressible materials or materials with isochoric plastic flow.
The base element type.
typedef ConstitutiveLawType::Pointer Kratos::SolidShellElementSprism3D6N::ConstitutiveLawPointerType |
Pointer type for constitutive laws.
Reference type definition for constitutive laws.
typedef std::size_t Kratos::SolidShellElementSprism3D6N::IndexType |
The definition of the index type.
Type definition for integration methods.
This is the definition of the node.
typedef std::size_t Kratos::SolidShellElementSprism3D6N::SizeType |
The definition of the sizetype.
StressMeasure from constitutive laws.
typedef GlobalPointersVector<NodeType> Kratos::SolidShellElementSprism3D6N::WeakPointerVectorNodesType |
Kratos::SolidShellElementSprism3D6N::SolidShellElementSprism3D6N | ( | ) |
Kratos::SolidShellElementSprism3D6N::SolidShellElementSprism3D6N | ( | IndexType | NewId, |
GeometryType::Pointer | pGeometry | ||
) |
Kratos::SolidShellElementSprism3D6N::SolidShellElementSprism3D6N | ( | IndexType | NewId, |
GeometryType::Pointer | pGeometry, | ||
PropertiesType::Pointer | pProperties | ||
) |
Kratos::SolidShellElementSprism3D6N::SolidShellElementSprism3D6N | ( | SolidShellElementSprism3D6N const & | rOther | ) |
|
override |
|
protected |
Update the LHS of the system with the EAS.
rLeftHandSideMatrix | LHS of the system |
rEAS | The components of the EAS stabilization |
|
protected |
Update the RHS of the system with the EAS and the internal variable alpha.
rRHSFull | The full internal forces vector |
rEAS | The components of the EAS stabilization |
AlphaEAS | The internal variable for the EAS |
|
protected |
Construction of the membrane deformation tangent matrix:
BMembrane | Membrane component of the deformation tangent matrix |
CMembrane | Membrane component of the Cauchy tensor |
InPlaneCartesianDerivativesGauss | The in-plane cartesian derivatives of the Gauss points |
InPlaneGradientFGauss | The in-plane deformation gradient components |
NodeGauss | Number of Gauss node calculated |
|
protected |
Construction of the transversal deformation tangent matrix:
BNormal | Transversal deformation tangent matrix |
TransversalCartesianDerivativesGaussCenter | Transversal cartesian derivatives in the central point of the element |
TransversalDeformationGradientF | Transversal components of the deformation gradient in the central point of the element |
|
protected |
Construction of the shear deformation tangent matrix:
BShear | Shear component of the deformation tangent matrix |
CShear | Shear components of the Cauchy tensor |
rCartesianDerivatives | Cartesian derivatives auxiliary struct |
rTransverseGradient | Local deformation gradient components for each Gauss point |
rTransverseGradientIsoParametric | Local deformation gradient components in the isogeometric space |
Part | The enum that indicates upper or lower face |
|
protected |
Calculation of the External Forces Vector. Fe = N * t + N * b.
rRightHandSideVector | RHS of the system |
rVariables | The internal variables in the element |
rVolumeForce | The force due to the acceleration of the body |
|
protected |
Calculation of the Internal Forces Vector. Fi = B * sigma.
rRightHandSideVector | RHS of the system |
rEAS | The components of the EAS stabilization |
|
protected |
Calculation of the Geometric Stiffness Matrix. Kuug = BT * S.
rLeftHandSideMatrix | LHS of the system |
|
protected |
Calculation of the Material Stiffness Matrix. Kuum = BT * C * B.
rLeftHandSideMatrix | LHS of the system |
rVariables | The internal variables in the element |
IntegrationWeight | Contribution in the numerical integration |
|
protected |
Calculation and addition of the matrix of the LHS.
rLocalSystem | The local system of equations |
rVariables | The internal variables in the element |
rValues | Values of the ContstitutiveLaw |
rEAS | The components of the EAS stabilization |
AlphaEAS | The internal variable for the EAS |
|
protected |
Construction of the in-plane geometric stiffness matrix:
Kgeometricmembrane | Membrane component of the stiffness matrix |
rCartesianDerivatives | Cartesian derivatives auxiliary struct |
Part | The enum that indicates upper or lower face |
|
protected |
Construction of the transversal geometric contribution to the stiffness matrix:
Kgeometricnormal | The transversal geometric contribution to the stiffness matrix |
TransversalCartesianDerivativesGaussCenter | Transversal cartesian derivatives in the central point of the element |
SNormal | Enhanced transversal component of the PK2 tensor |
|
protected |
Calculation and addition of the vectors of the RHS.
rLocalSystem | The local system of equations |
rVariables | The internal variables in the element |
rVolumeForce | The force due to the acceleration of the body |
rEAS | The components of the EAS stabilization |
AlphaEAS | The internal variable for the EAS |
|
protected |
Construction of the shear geometric contribution to the stiffness matrix:
Kgeometricshear | The shear geometric contribution to the stiffness matrix |
rCartesianDerivatives | Cartesian derivatives auxiliary struct |
SShear | The shear components of the PK2 tensor |
Part | The enum that indicates upper or lower face |
|
protected |
Calculate the cartesian derivatives.
|
protected |
Calculate the Cartesian derivatives in the Gauss points, for the plane.
CartesianDerivativesCenter | The cartesian derivatives in the plane |
Part | The enum that indicates upper or lower face |
|
protected |
Calculate the Cartesian derivatives in the center, for the transversal direction.
rCartesianDerivatives | The class containing the cartesian derivatives |
NodesCoord | The matrix with the coordinates of the nodes of the element |
Part | 0 for center node of the element, 1 for upper part and 2 for lower part |
|
protected |
Calculate the Cartesian derivatives in the Gauss points, for the plane.
NodeGauss | Number of Gauss node calculated |
Part | The index that indicates upper or lower face |
NodesCoord | The matrix with the coordinates of the nodes of the element |
InPlaneCartesianDerivativesGauss | The cartesian derivatives in the plane |
|
protected |
Calculate the Cartesian derivatives in the Gauss points, for the transversal direction.
NodesCoord | The matrix with the coordinates of the nodes of the element |
TransversalCartesianDerivativesGauss | The cartesian derivatives in the transversal direction |
rLocalCoordinates | The local coordinates |
|
protected |
Calculate the necessary components for the Kinematic calculus.
void Kratos::SolidShellElementSprism3D6N::CalculateDampingMatrix | ( | MatrixType & | rDampingMatrix, |
const MatrixType & | rStiffnessMatrix, | ||
const MatrixType & | rMassMatrix, | ||
const ProcessInfo & | rCurrentProcessInfo | ||
) |
This is called during the assembling process in order to calculate the elemental damping matrix (Reusing the stiffness matrix and mass matrix)
rDampingMatrix | the elemental damping matrix |
rStiffnessMatrix | the elemental stiffness matrix |
rMassMatrix | the elemental mass matrix |
rCurrentProcessInfo | the current process info instance |
|
overridevirtual |
This is called during the assembling process in order to calculate the elemental damping matrix.
rDampingMatrix | the elemental damping matrix |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Calculation of the Deformation Matrix BL.
rB | Deformation matrix |
ZetaGauss | The zeta coordinate for the Gauss Quadrature |
AlphaEAS | The internal variable for the EAS |
|
protected |
This method computes the delta position matrix necessary for UL formulation.
rDeltaPosition | The matrix that contents the increase of displacements |
|
protected |
Calculates the elemental contributions.
rLocalSystem | The local system of equations |
rCurrentProcessInfo | The current process info instance |
|
protected |
Calculate the vector of the element Ids.
|
protected |
Calculate the components of the deformation gradient in the plane, for the Gauss nodes:
InPlaneGradientFGauss | The components of the deformation gradient in the plane, for the gauss node |
InPlaneCartesianDerivativesGauss | The cartesian derivatives of a Gauss node in the plane |
NodesCoord | The coordinates of the nodes of the element |
NodeGauss | Number of Gauss node calculated |
Part | The enum that indicates upper or lower face |
|
protected |
Calculate the Jacobian.
detJ | Determinant of the Jacobian |
J | The Jacobian of the element |
LocalDerivativePatch | The local derivatives of the element |
NodesCoord | The matrix with the coordinates of the nodes of the element |
rLocalCoordinates | The local coordinates |
|
protected |
Calculate the Jacobian and its inverse.
J | The Jacobian of the element |
Jinv | The inverse of the Jacobian |
LocalDerivativePatch | The local derivatives of the element |
NodesCoord | The matrix with the coordinates of the nodes of the element |
rLocalCoordinates | The local coordinates |
|
protected |
Calculate the Jacobian and his inverse.
J | The Jacobian of the element |
Jinv | The inverse of the Jacobian |
NodesCoord | The matrix with the coordinates of the nodes of the element |
rLocalCoordinates | The local coordinates |
|
protected |
Calculate the Jacobian and his inverse.
J | The Jacobian of the element |
Jinv | The inverse of the Jacobian |
detJ | Determinant of the Jacobian |
rPointNumber | The integration point index |
ZetaGauss | The transversal local coordinates |
|
protected |
Calculate Element Kinematics.
rVariables | The internal variables in the element |
rIntegrationPoints | The integration points of the prism |
rPointNumber | The integration point index |
AlphaEAS | The internal variable for the EAS |
ZetaGauss | The zeta coordinate for the Gauss Quadrature |
|
overridevirtual |
This is called during the assembling process in order to calculate the elemental leTransverseGradientFt hand side vector only.
rLeftHandSideMatrix | the elemental leTransverseGradientFt hand side vector |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Calculates the Local Coordinates System.
ThisOrthogonalBaseApproach | The chosen approximation |
ThisAngle | Angle of rotation of the element |
|
overridevirtual |
This is called during the assembling process in order to calculate all elemental contributions to the global system matrix and the right hand side.
rLeftHandSideMatrix | the elemental leTransverseGradientFt hand side matrix |
rRightHandSideVector | the elemental right hand side |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
This is called during the assembling process in order to calculate the elemental mass matrix.
rMassMatrix | the elemental mass matrix |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Calculate a 3 components array_1d on the Element Constitutive Law.
rVariable | The variable we want to get |
rOutput | The values obtained int the integration points |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Calculate a 6 components array_1d on the Element Constitutive Law.
rVariable | The variable we want to get |
rOutput | The values obtained int the integration points |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Calculate a boolean Variable on the Element Constitutive Law.
rVariable | The internal variables in the element |
rOutput | The solution (boolean) |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Calculate a double Variable on the Element Constitutive Law.
rVariable | The internal variables in the element |
rOutput | The solution (double) |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Calculate a integer Variable on the Element Constitutive Law.
rVariable | The internal variables in the element |
rOutput | The solution (integer) |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Calculate a Matrix Variable on the Element Constitutive Law.
rVariable | The internal variables in the element |
rOutput | The matrix solution |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Calculate a Vector Variable on the Element Constitutive Law.
rVariable | The internal variables in the element |
rOutput | The vector solution |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
This is called during the assembling process in order to calculate the elemental right hand side vector only.
rRightHandSideVector | the elemental right hand side vector |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Calculate the transversal components of the deformation gradient, in the Gauss points:
TransverseGradientF | The transversal components of the deformation gradient |
TransversalCartesianDerivativesGauss | The transversal cartesian derivatives |
NodesCoord | The coordinates of the nodes of the element |
|
protected |
Calculate the transversal components of the deformation gradient, in each one of the faces:
rTransverseGradientIsoParametric | Auxilar components of the deformation gradient |
NodesCoord | The coordinates of the nodes of the element |
Part | The enum that indicates if calculate upper or lower components |
|
protected |
This function calculates the variation of the element volume.
rVolumeChange | Volume variation of the element |
rVariables | The internal variables in the element |
|
protected |
Calculation of the Volume Force of the Element.
rVolumeForce | The volume forces of the element |
rVariables | The internal variables in the element |
IntegrationWeight | Contribution in the numerical integration |
|
protected |
Calculate Fbar from Cbar.
Assuming that the rotation matrix of the polar decomposition of the F_bar is the same of the polar decomposition of F
rVariables | The internal variables in the element |
rPointNumber | The integration point index |
|
overridevirtual |
This function provides the place to perform checks on the completeness of the input.
It is designed to be called only once (or anyway, not oTransverseGradientFten) typically at the beginning of the calculations, so to verify that nothing is missing from the input or that no common error is found.
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Clones the selected element variables, creating a new one
NewId | the ID of the new element |
ThisNodes | the nodes of the new element |
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Calculate the local derivatives of the element for a given coordinates.
LocalDerivativePatch | The local derivatives of the element |
rLocalCoordinates | The local coordinates |
|
protected |
Calculate the local quadratic derivatives of the element for a given gauss node.
rLocalDerivativePatch | The local derivatives of the element |
NodeGauss | The Gauss node index |
|
overridevirtual |
Creates a new element.
NewId | The Id of the new created element |
pGeom | The pointer to the geometry of the element |
pProperties | The pointer to property |
Reimplemented from Kratos::Element.
|
overridevirtual |
Creates a new element.
NewId | The Id of the new created element |
ThisNodes | The array containing nodes |
pProperties | The pointer to property |
Reimplemented from Kratos::Element.
|
overridevirtual |
Sets on rResult the ID's of the element degrees of freedom.
rResult | The result vector with the ID's of the DOF |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
This is called for non-linear analysis at the beginning of the iteration process.
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Called at the end of each solution step.
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Finalize Element Internal Variables.
rVariables | The internal variables in the element |
rPointNumber | The integration point index |
|
overridevirtual |
Sets on rElementalDofList the degrees of freedom of the considered element geometry.
rElementalDofList | The list of the degrees of freedom from the element |
rCurrentProcessInfo | the current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Sets on rValues the nodal velocities.
Step | The calculation step |
rValues | The velocities vector |
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Get the Historical Deformation Gradient to calculate aTransverseGradientFter finalize the step.
rVariables | The internal variables in the element |
rPointNumber | The integration point index |
|
protected |
It gets the nodal coordinates, according to the configutaion.
|
overridevirtual |
Sets on rValues the nodal accelerations.
Step | The calculation step |
rValues | The accelerations vector |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Sets on rValues the nodal displacements.
Step | The calculation step |
rValues | The displacements vector |
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Calculates the vector of current position.
|
protected |
Calculates the vector of previous position
|
protected |
Check if the node has a neighbour:
Index | The index of the node |
NeighbourNode | The neighbours nodes |
|
inlineoverridevirtual |
Turn back information as a string.
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Called to initialize the element.
Must be called before any calculation is done
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Initialize Element General Variables.
rVariables | The internal variables in the element |
|
overridevirtual |
This is called for non-linear analysis at the beginning of the iteration process.
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Called at the beginning of each solution step.
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Initialize System Matrices.
rLeftHandSideMatrix | LHS of the system |
rRightHandSideVector | RHS of the system |
rCalculationFlags | Calculation flags |
|
protected |
Integrates the EAS components in zeta using the Gauss Quadrature.
rVariables | The internal variables in the element |
rEAS | The components of the EAS stabilization |
ZetaGauss | The zeta coordinate for the Gauss Quadrature |
IntegrationWeight | Contribution in the numerical integration |
|
protected |
Integrates stresses in zeta using the Gauss Quadrature.
rVariables | The internal variables in the element |
AlphaEAS | The internal variable for the EAS |
ZetaGauss | The zeta coordinate for the Gauss Quadrature |
IntegrationWeight | Contribution in the numerical integration |
Kratos::SolidShellElementSprism3D6N::KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION | ( | SolidShellElementSprism3D6N | ) |
Counted pointer of SolidShellElementSprism3D6N.
Kratos::SolidShellElementSprism3D6N::KRATOS_DEFINE_LOCAL_FLAG | ( | COMPUTE_LHS_MATRIX | ) |
Kratos::SolidShellElementSprism3D6N::KRATOS_DEFINE_LOCAL_FLAG | ( | COMPUTE_LHS_MATRIX_WITH_COMPONENTS | ) |
Kratos::SolidShellElementSprism3D6N::KRATOS_DEFINE_LOCAL_FLAG | ( | COMPUTE_RHS_VECTOR | ) |
Flags related to the element computation.
Kratos::SolidShellElementSprism3D6N::KRATOS_DEFINE_LOCAL_FLAG | ( | COMPUTE_RHS_VECTOR_WITH_COMPONENTS | ) |
Kratos::SolidShellElementSprism3D6N::KRATOS_DEFINE_LOCAL_FLAG | ( | EAS_IMPLICIT_EXPLICIT | ) |
Kratos::SolidShellElementSprism3D6N::KRATOS_DEFINE_LOCAL_FLAG | ( | EXPLICIT_RHS_COMPUTATION | ) |
Kratos::SolidShellElementSprism3D6N::KRATOS_DEFINE_LOCAL_FLAG | ( | QUADRATIC_ELEMENT | ) |
Kratos::SolidShellElementSprism3D6N::KRATOS_DEFINE_LOCAL_FLAG | ( | TOTAL_UPDATED_LAGRANGIAN | ) |
|
protected |
Calculates the number of active neighbours:
pNeighbourNodes | The neighbours nodes |
SolidShellElementSprism3D6N & Kratos::SolidShellElementSprism3D6N::operator= | ( | SolidShellElementSprism3D6N const & | rOther | ) |
|
inlineoverridevirtual |
Print object's data.
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Prints element information for each gauss point (debugging purposes)
rLocalSystem | The local system of equations |
rVariables | The internal variables in the element |
|
inlineoverridevirtual |
Print information about this object.
Reimplemented from Kratos::BaseSolidElement.
|
protected |
Set Variables of the Element to the Parameters of the Constitutive Law.
rVariables | The internal variables in the element |
rValues | Values of the ContstitutiveLaw |
rPointNumber | The integration point index |
|
overridevirtual |
Set a Constitutive Law Value.
rVariable | The internal variables in the element |
rValues | Values of the ContstitutiveLaw |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Set a double Value on the Element Constitutive Law.
Access for variables on Integration points. This gives access to variables stored in the constitutive law on each integration point. Specialisations of element.h (e.g. the TotalLagrangian) must specify the actual interface to the constitutive law! Note, that these functions expect a std::vector of values for the specified variable type that contains a value for each integration point! SetValuesOnIntegrationPoints: Set the values for given Variable.
rVariable | The internal variables in the element |
rValues | Values of the ContstitutiveLaw |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Set a Matrix Value on the Element Constitutive Law.
rVariable | The internal variables in the element |
rValues | Values of the ContstitutiveLaw |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
overridevirtual |
Set a Vector Value on the Element Constitutive Law.
rVariable | The internal variables in the element |
rValues | Values of the ContstitutiveLaw |
rCurrentProcessInfo | The current process info instance |
Reimplemented from Kratos::BaseSolidElement.
|
friend |
|
protected |
|
protected |
Container for historical total Jacobians for Total Lagrangian Container for historical total elastic deformation measure F0 = dx/dX for Updated Lagrangian
|
protected |