61 #include "domain/mesh/element/plane/TriBase3N.h" 62 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h" 63 #include "domain/mesh/element/utils/body_forces/BodyForces2D.h" 84 static double matrixData[36];
88 static double shp[3][3];
91 double shapeFunction(
const GaussPoint &)
const;
92 void setPressureLoadAtNodes(
void);
99 Tri31(
int tag,
int nd1,
int nd2,
int nd3,
101 double t,
double pressure = 0.0,
107 virtual ~
Tri31(
void);
117 const Matrix &getInitialStiff(
void)
const;
122 inline double getRho(
void)
const 124 void setRho(
const double &r)
126 double getThickness(
void)
const 128 void setThickness(
const double &t)
133 int addInertiaLoadToUnbalance(
const Vector &accel);
141 void Print(std::ostream &s,
int flag =0)
const;
153 const Vector &getAvgStress(
void)
const;
154 const Vector &getAvgStrain(
void)
const;
int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Tri31.cpp:700
int update(void)
Updates the element state.
Definition: Tri31.cpp:140
Body forces over an element.
Definition: BodyForces2D.h:40
Element * getCopy(void) const
Virtual constructor.
Definition: Tri31.cpp:97
Float vector abstraction.
Definition: Vector.h:94
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:673
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Tri31.cpp:302
Communication parameters between processes.
Definition: Communicator.h:66
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: Tri31.cpp:184
Base class response objects.
Definition: Response.h:81
double getThickness(void) const
Return the material thickness.
Definition: SolidMech2D.h:55
int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: Tri31.cpp:522
Uniaxial p-y material that incorporates liquefaction effects.
Definition: PyLiq1.h:61
SolidMech2D physicalProperties
pointers to the material objects and physical properties.
Definition: ElemWithMaterial.h:50
3D position of Gauss points.
Definition: GaussPoint.h:38
void setRho(const double &)
Set the density for all materials.
Definition: NDMaterialPhysicalProperties.cc:197
Base class for the finite elements.
Definition: Element.h:112
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: Tri31.cpp:109
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Tri31.cpp:649
int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: Tri31.cpp:537
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: Tri31.cpp:555
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: Tri31.cpp:112
Base class for 3 node triangles.
Definition: TriBase3N.h:46
Uniaxial t-z material that incorporates liquefaction effects.
Definition: TzLiq1.h:61
3 node triangle.
Definition: Tri31.h:74
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:604
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: Tri31.cpp:405
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: Tri31.cpp:453
int sendData(Communicator &)
Send object members through the communicator argument.
Definition: Tri31.cpp:503
Base class for Gauss integration models.
Definition: GaussModel.h:41
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
void alive(void)
Reactivates the element.
Definition: Tri31.cpp:343
double getRho(void) const
Returns the average of the densities for each material.
Definition: NDMaterialPhysicalProperties.cc:186
void setThickness(const double &t)
Set the material thickness.
Definition: SolidMech2D.h:58
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
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
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: Tri31.cpp:513