52 #include "domain/mesh/element/plane/SolidMech4N.h" 53 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h" 54 #include <utility/matrix/Vector.h> 55 #include <utility/matrix/Matrix.h> 75 static const double sg[4];
76 static const double tg[4];
77 static const double wg[4];
80 static double stressData[][4];
81 static double tangentData[][3][4];
86 static double xl[][4];
89 static void saveData(
int gp,
const Vector &stress,
const Matrix &tangent);
92 void getData(
int gp,
Vector &stress,
Matrix &tangent)
const;
95 const Matrix& computeBenhanced(
int node,
double L1,
double L2,
double j,
const Matrix &Jinv)
const;
99 void computeBasis(
void)
const;
102 void formResidAndTangent(
int tang_flag)
const;
105 void formInertiaTerms(
int tangFlag)
const;
109 static void computeJacobian(
double L1,
double L2,
const double x[2][4],
Matrix &JJ,
Matrix &JJinv );
112 const Matrix& computeB(
int node,
const double shp[3][4])
const;
118 static void shape2d(
double ss,
double tt,
const double x[2][4],
double shp[3][4],
double &xsj );
142 void Print(std::ostream &s,
int flag)
const;
151 int addInertiaLoadToUnbalance(
const Vector &accel);
const Matrix & getMass(void) const
return mass matrix
Definition: EnhancedQuad.cpp:442
void setDomain(Domain *)
Set domain.
Definition: EnhancedQuad.cpp:125
int update(void)
Updates the element state.
Definition: EnhancedQuad.cpp:1018
Float vector abstraction.
Definition: Vector.h:94
int getNumDOF(void) const
return number of dofs
Definition: EnhancedQuad.cpp:132
Communication parameters between processes.
Definition: Communicator.h:66
void Print(std::ostream &s, int flag) const
Print stuff.
Definition: EnhancedQuad.cpp:137
void alive(void)
Reactivates the element.
Definition: EnhancedQuad.cpp:453
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: EnhancedQuad.cpp:154
virtual int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: EnhancedQuad.cpp:1336
int sendData(Communicator &)
Send members through the communicator argument.
Definition: EnhancedQuad.cpp:1280
Element * getCopy(void) const
Virtual constructor.
Definition: EnhancedQuad.cpp:116
const Vector & getResistingForce(void) const
Get resisting force (residual vector).
Definition: EnhancedQuad.cpp:507
Four-node quadrilateral element, which uses a bilinear isoparametric formulation with enhanced strain...
Definition: EnhancedQuad.h:62
EnhancedQuad(int tag=0, const NDMaterial *ptr_mat=nullptr)
Constructor.
Definition: EnhancedQuad.cpp:90
Base class for the finite elements.
Definition: Element.h:112
int updateParameter(int parameterID, Information &info)
Update parameter value.
Definition: EnhancedQuad.cpp:1316
int recvData(const Communicator &)
Receives members through the communicator argument.
Definition: EnhancedQuad.cpp:1288
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: EnhancedQuad.cpp:1351
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
Four node quadrilateral element for two-dimensional problems.
Definition: SolidMech4N.h:76
const Matrix & getInitialStiff(void) const
Return secant matrix.
Definition: EnhancedQuad.cpp:166
Parameter.
Definition: Parameter.h:68
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
const Vector & getResistingForceIncInertia(void) const
Get resisting force with inertia terms.
Definition: EnhancedQuad.cpp:521
int setParameter(const std::vector< std::string > &argv, Parameter ¶m)
Set parameter value.
Definition: EnhancedQuad.cpp:1296