60 #include "domain/mesh/element/plane/TriBase3N.h" 61 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h" 62 #include "domain/mesh/element/utils/body_forces/BodyForces2D.h" 83 static double matrixData[36];
87 static double shp[3][3];
90 double shapeFunction(
const GaussPoint &)
const;
91 void setPressureLoadAtNodes(
void);
98 Tri31(
int tag,
int nd1,
int nd2,
int nd3,
100 double t,
double pressure = 0.0,
106 virtual ~
Tri31(
void);
116 const Matrix &getInitialStiff(
void)
const;
121 inline double getRho(
void)
const 123 void setRho(
const double &r)
125 double getThickness(
void)
const 127 void setThickness(
const double &t)
131 int addInertiaLoadToUnbalance(
const Vector &accel);
139 void Print(std::ostream &s,
int flag =0);
151 const Vector &getAvgStress(
void)
const;
152 const Vector &getAvgStrain(
void)
const;
int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Tri31.cpp:658
int update(void)
Updates the element state.
Definition: Tri31.cpp:140
Body forces over an element.
Definition: BodyForces2D.h:39
Element * getCopy(void) const
Virtual constructor.
Definition: Tri31.cpp:97
Float vector abstraction.
Definition: Vector.h:93
const GaussModel & getGaussModel(void) const
Return the Gauss points of the element.
Definition: Tri31.cpp:339
int setParameter(const std::vector< std::string > &argv, Parameter ¶m)
Sets the value param to the parameter argv.
Definition: Tri31.cpp:631
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Tri31.cpp:302
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: Tri31.cpp:486
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: Tri31.cpp:184
Definition: Response.h:71
Uniaxial p-y material that incorporates liquefaction effects.
Definition: PyLiq1.h:60
SolidMech2D physicalProperties
pointers to the material objects and physical properties.
Definition: ElemWithMaterial.h:43
3D position of Gauss points.
Definition: GaussPoint.h:37
Base class for the finite elements.
Definition: Element.h:109
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: Tri31.cpp:109
int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: Tri31.cpp:505
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Tri31.cpp:607
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: Tri31.cpp:112
Base class for 3 node triangles.
Definition: TriBase3N.h:45
Uniaxial t-z material that incorporates liquefaction effects.
Definition: TzLiq1.h:60
3 node triangle.
Definition: Tri31.h:73
Response * setResponse(const std::vector< std::string > &argv, Information &eleInformation)
setResponse() is a method invoked to determine if the element will respond to a request for a certain...
Definition: Tri31.cpp:587
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: Tri31.cpp:388
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: Tri31.cpp:436
Base class for Gauss integration models.
Definition: GaussModel.h:40
Base class for 2D and 3D materials.
Definition: NDMaterial.h:97
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: Tri31.cpp:496
int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: Tri31.cpp:520
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
void Print(std::ostream &s, int flag=0)
Print stuff.
Definition: Tri31.cpp:538
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
Definition: Parameter.h:65
Base class for loads over elements.
Definition: ElementalLoad.h:77
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:116