9 #if !defined(KRATOS_UPDATED_LAGRANGIAN_U_P_wP_ELEMENT_H_INCLUDED )
10 #define KRATOS_UPDATED_LAGRANGIAN_U_P_wP_ELEMENT_H_INCLUDED
129 void GetValuesVector(
Vector& rValues,
int Step = 0)
const override;
134 void GetFirstDerivativesVector(
Vector& rValues,
int Step = 0)
const override;
139 void GetSecondDerivativesVector(
Vector& rValues,
int Step = 0)
const override;
144 void FinalizeSolutionStep(
ProcessInfo& rCurrentProcessInfo)
override;
155 int Check(
const ProcessInfo& rCurrentProcessInfo)
override;
202 double& rIntegrationWeight)
override;
211 double& rIntegrationWeight)
override;
221 virtual void CalculateAndAddUnconsideredKuuTerms(
MatrixType& rK,
223 double& rIntegrationWeight
229 virtual void CalculateAndAddKuwP(
MatrixType& rK,
231 double& rIntegrationWeight
237 virtual void CalculateAndAddKwPu(
MatrixType& rK,
239 double& rIntegrationWeight
245 virtual void CalculateAndAddKwPP(
MatrixType& rK,
247 double& rIntegrationWeight
253 virtual void CalculateAndAddKwPwP(
MatrixType& rK,
255 double& rIntegrationWeight
261 virtual void CalculateAndAddKwPwPStab(
MatrixType& rK,
263 double& rIntegrationWeight
269 void CalculateAndAddExternalForces(
VectorType& rRightHandSideVector,
272 double& rIntegrationWeight
279 void CalculateAndAddPressureForces(
VectorType& rRightHandSideVector,
281 double& rIntegrationWeight
288 void CalculateAndAddStabilizedPressure(
VectorType& rRightHandSideVector,
290 double& rIntegrationWeight
296 void CalculateAndAddInternalForces(
VectorType& rRightHandSideVector,
298 double& rIntegrationWeight
304 virtual void CalculateAndAddWaterPressureForces(
VectorType& rRightHandSideVector,
306 double& rIntegrationWeight
311 virtual void CalculateAndAddStabilizedWaterPressure(
VectorType& rRightHandSideVector,
313 double& rIntegartionWeight
319 void InitializeSystemMatrices(
MatrixType& rLeftHandSideMatrix,
321 Flags& rCalculationFlags)
override;
327 double& CalculateVolumeChange(
double& rVolumeChange,
ElementDataType& rVariables)
override;
330 void GetConstants(
double& rScalingConstant,
double& rWaterBulk,
double& rDeltaTime,
double& rPermeability);
332 virtual double GetElementSize(
const Matrix& rDN_DX);
376 void save(
Serializer& rSerializer)
const override;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
IntegrationMethod
Definition: geometry_data.h:76
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Large Displacement Lagrangian U-P Element for 3D and 2D geometries. Linear Triangles and Tetrahedra (...
Definition: updated_lagrangian_U_Pressure_element.hpp:46
Updated Lagrangian Large Displacement Lagrangian U-wP Element for 3D and 2D geometries....
Definition: updated_lagrangian_U_P_wP_element.hpp:42
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: updated_lagrangian_U_P_wP_element.hpp:54
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UpdatedLagrangianUPwPElement)
Counted pointer of LargeDisplacementUPElement.
double mTimeStep
Definition: updated_lagrangian_U_P_wP_element.hpp:178
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: updated_lagrangian_U_P_wP_element.hpp:52
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: updated_lagrangian_U_P_wP_element.hpp:50
ConstitutiveLaw ConstitutiveLawType
Definition: updated_lagrangian_U_P_wP_element.hpp:48
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307
Definition: solid_element.hpp:83
Definition: solid_element.hpp:233