57 #include <domain/mesh/element/ElemWithMaterial.h> 58 #include <utility/matrix/Vector.h> 59 #include <utility/matrix/Matrix.h> 60 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h" 72 static double xl[][9];
91 void formResidAndTangent(
int tang_flag)
const;
93 void formInertiaTerms(
int tangFlag)
const;
95 const Matrix& computeBbar(
int node,
96 const double natCoor[2],
97 const double shp[3][9],
98 double shpBar[3][9][3] )
const;
100 static void shape2dNine(
double coor[2],
const double x[2][9],
double shp[3][9],
double &xsj );
101 void computeBasis(
void)
const;
102 static double shape1d(
int code,
int node,
double xi );
114 int node1,
int node2,
int node3,
int node4,
int node5,
int node6,
int node7,
int node8,
int node9,
127 void Print(std::ostream &s,
int flag);
131 const Matrix &getInitialStiff(
void)
const;
135 int addInertiaLoadToUnbalance(
const Vector &accel);
Float vector abstraction.
Definition: Vector.h:93
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: NineNodeMixedQuad.cpp:1101
virtual int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: NineNodeMixedQuad.cpp:1123
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: NineNodeMixedQuad.cpp:1112
void setDomain(Domain *)
Sets the domain for the element.
Definition: NineNodeMixedQuad.cpp:110
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: NineNodeMixedQuad.cpp:421
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: NineNodeMixedQuad.cpp:436
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: NineNodeMixedQuad.cpp:117
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: NineNodeMixedQuad.cpp:146
Nine node quad.
Definition: NineNodeMixedQuad.h:69
Base class for the finite elements.
Definition: Element.h:109
virtual ~NineNodeMixedQuad(void)
destructor
Definition: NineNodeMixedQuad.cpp:103
Element * getCopy(void) const
Virtual constructor.
Definition: NineNodeMixedQuad.cpp:99
virtual int recvSelf(const CommParameters &)
Receive the object.
Definition: NineNodeMixedQuad.cpp:1137
Base class for 2D and 3D materials.
Definition: NDMaterial.h:97
void Print(std::ostream &s, int flag)
Print stuff.
Definition: NineNodeMixedQuad.cpp:122
const Matrix & getMass(void) const
return mass matrix
Definition: NineNodeMixedQuad.cpp:366
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
Base class for loads over elements.
Definition: ElementalLoad.h:77
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:116