58 #include "domain/mesh/element/ElemWithMaterial.h"    59 #include <utility/matrix/Vector.h>    60 #include <utility/matrix/Matrix.h>    61 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h"    74     static double xl[][9]; 
    93     void formResidAndTangent(
int tang_flag) 
const;
    95     void formInertiaTerms(
int tangFlag) 
const;
    97     const Matrix& computeBbar( 
int node, 
    98                    const double natCoor[2], 
    99                    const double shp[3][9], 
   100                    double shpBar[3][9][3] ) 
const;
   102     static void shape2dNine( 
double coor[2], 
const double x[2][9], 
double shp[3][9], 
double &xsj );
   103     void computeBasis(
void) 
const;
   104     static double shape1d( 
int code, 
int node, 
double xi );
   115                int node1, 
int node2, 
int node3, 
int node4, 
int node5, 
int node6, 
int node7, 
int node8, 
int node9,
   128     void Print(std::ostream &s, 
int flag) 
const;
   132     const Matrix &getInitialStiff(
void) 
const;     
   137     int addInertiaLoadToUnbalance(
const Vector &accel);
 Float vector abstraction. 
Definition: Vector.h:94
Communication parameters between processes. 
Definition: Communicator.h:66
int recvData(const Communicator &)
Receives object members through the communicator argument. 
Definition: NineNodeMixedQuad.cpp:1126
void alive(void)
Reactivates the element. 
Definition: NineNodeMixedQuad.cpp:376
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:435
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces. 
Definition: NineNodeMixedQuad.cpp:450
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
Mixed pressure/volume nine node quadrilateral element for two-dimensional plane strain problems...
Definition: NineNodeMixedQuad.h:71
Base class for the finite elements. 
Definition: Element.h:112
virtual ~NineNodeMixedQuad(void)
destructor 
Definition: NineNodeMixedQuad.cpp:103
void Print(std::ostream &s, int flag) const
Print stuff. 
Definition: NineNodeMixedQuad.cpp:122
Element * getCopy(void) const
Virtual constructor. 
Definition: NineNodeMixedQuad.cpp:99
Base class for 2D and 3D materials. 
Definition: NDMaterial.h:101
const Matrix & getMass(void) const
return mass matrix 
Definition: NineNodeMixedQuad.cpp:366
int sendData(Communicator &)
Send object members through the communicator argument. 
Definition: NineNodeMixedQuad.cpp:1115
Open source finite element program for structural analysis. 
Definition: ContinuaReprComponent.h:35
Element with material. 
Definition: ElemWithMaterial.h:45
Matrix of floats. 
Definition: Matrix.h:111
Base class for loads over elements. 
Definition: ElementalLoad.h:79
virtual int sendSelf(Communicator &)
Sends object through the communicator argument. 
Definition: NineNodeMixedQuad.cpp:1137
Domain (mesh and boundary conditions) of the finite element model. 
Definition: Domain.h:117
virtual int recvSelf(const Communicator &)
Receive the object. 
Definition: NineNodeMixedQuad.cpp:1151