43 #ifndef NineFourNodeQuadUP_h 44 #define NineFourNodeQuadUP_h 46 #include <domain/mesh/element/ElemWithMaterial.h> 47 #include <utility/matrix/Matrix.h> 48 #include <utility/matrix/Vector.h> 49 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h" 50 #include "domain/mesh/element/utils/body_forces/BodyForces2D.h" 71 static const int nintu;
72 static const int nintp;
73 static const int nenu;
74 static const int nenp;
77 static double shgu[3][9][9];
78 static double shgp[3][4][4];
79 static double shgq[3][9][4];
80 static double shlu[3][9][9];
81 static double shlp[3][4][4];
82 static double shlq[3][9][4];
85 static double dvolu[9];
86 static double dvolp[4];
87 static double dvolq[4];
90 double mixtureRho(
int ipt)
const;
91 void shapeFunction(
double *w,
int nint,
int nen,
int mode)
const;
92 void globalShapeFunction(
double *dvol,
double *w,
int nint,
int nen,
int mode)
const;
94 double *initNodeDispl;
101 int nd5,
int nd6,
int nd7,
int nd8,
int nd9,
103 double t,
double bulk,
double rhof,
double perm1,
double perm2,
117 const Matrix &getInitialStiff(
void)
const;
123 int addInertiaLoadToUnbalance(
const Vector &);
130 void Print(std::ostream &s,
int flag =0);
Body forces over an element.
Definition: BodyForces2D.h:39
Float vector abstraction.
Definition: Vector.h:93
int update(void)
Updates the element state.
Definition: Nine_Four_Node_QuadUP.cpp:197
Definition: Response.h:71
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Nine_Four_Node_QuadUP.cpp:519
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Nine_Four_Node_QuadUP.cpp:1076
Element * getCopy(void) const
Virtual constructor.
Definition: Nine_Four_Node_QuadUP.cpp:144
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: Nine_Four_Node_QuadUP.cpp:1018
void Print(std::ostream &s, int flag=0)
Print stuff.
Definition: Nine_Four_Node_QuadUP.cpp:1003
Uniaxial p-y material that incorporates liquefaction effects.
Definition: PyLiq1.h:60
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: Nine_Four_Node_QuadUP.cpp:247
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: Nine_Four_Node_QuadUP.cpp:949
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: Nine_Four_Node_QuadUP.cpp:163
Base class for the finite elements.
Definition: Element.h:109
Uniaxial t-z material that incorporates liquefaction effects.
Definition: TzLiq1.h:60
virtual ~NineFourNodeQuadUP(void)
Destructor.
Definition: Nine_Four_Node_QuadUP.cpp:148
int setParameter(const std::vector< std::string > &argv, Parameter ¶m)
Sets the value param to the parameter argv.
Definition: Nine_Four_Node_QuadUP.cpp:1115
Base class for 2D and 3D materials.
Definition: NDMaterial.h:97
const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: Nine_Four_Node_QuadUP.cpp:388
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: Nine_Four_Node_QuadUP.cpp:834
int recvData(const CommParameters &)
Receives object members through the channel.
Definition: Nine_Four_Node_QuadUP.cpp:959
int recvSelf(const CommParameters &)
Receive the object.
Definition: Nine_Four_Node_QuadUP.cpp:982
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
Element with material.
Definition: ElemWithMaterial.h:40
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
Definition: Parameter.h:65
void zeroLoad(void)
Zeroes the loads over the element.
Definition: Nine_Four_Node_QuadUP.cpp:611
Base class for loads over elements.
Definition: ElementalLoad.h:77
Definition: Nine_Four_Node_QuadUP.h:57
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: Nine_Four_Node_QuadUP.cpp:719
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:116
int sendSelf(CommParameters &)
Send the object.
Definition: Nine_Four_Node_QuadUP.cpp:969
int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Nine_Four_Node_QuadUP.cpp:1165
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: Nine_Four_Node_QuadUP.cpp:158