57 #ifndef TWENTYNODEBRICK_U_P_U_H 58 #define TWENTYNODEBRICK_U_P_U_H 62 #include <domain/mesh/element/ElemWithMaterial.h> 63 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h" 64 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h" 81 static const int Num_IntegrationPts;
82 static const int Num_TotalGaussPts;
83 static const int Num_Nodes;
84 static const int Num_Dim;
85 static const int Num_Dof;
86 static const int Num_ElemDof;
87 static const double pts[3];
88 static const double wts[3];
103 static BJtensor shapeFunction(
double,
double,
double);
104 static BJtensor shapeFunctionDerivative(
double,
double,
double);
112 BJtensor getStiffnessTensorKep(
void)
const;
113 BJtensor getStiffnessTensorG12(
void)
const;
114 BJtensor getStiffnessTensorP(
void)
const;
115 BJtensor getMassTensorMsf(
void)
const;
116 BJtensor getDampTensorC123(
void)
const;
117 const Matrix &getStiff(
int Ki_flag)
const;
118 double getPorePressure(
double,
double,
double)
const;
124 int node_numb_1,
int node_numb_2,
int node_numb_3,
int node_numb_4,
125 int node_numb_5,
int node_numb_6,
int node_numb_7,
int node_numb_8,
126 int node_numb_9,
int node_numb_10,
int node_numb_11,
int node_numb_12,
127 int node_numb_13,
int node_numb_14,
int node_numb_15,
int node_numb_16,
128 int node_numb_17,
int node_numb_18,
int node_numb_19,
int node_numb_20,
130 double nn,
double alf,
double rs,
double rf,
131 double permb_x,
double permb_y,
double permb_z,
132 double kks,
double kkf,
double pp);
146 const Matrix &getInitialStiff(
void)
const;
152 int addInertiaLoadToUnbalance(
const Vector &accel);
159 void Print(std::ostream &s,
int flag =0);
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: TwentyNodeBrick_u_p_U.cpp:173
const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: TwentyNodeBrick_u_p_U.cpp:193
Float vector abstraction.
Definition: Vector.h:93
Element * getCopy(void) const
Virtual constructor.
Definition: TwentyNodeBrick_u_p_U.cpp:150
Definition: Response.h:71
Definition: BJtensor.h:110
void Print(std::ostream &s, int flag=0)
Print stuff.
Definition: TwentyNodeBrick_u_p_U.cpp:491
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: TwentyNodeBrick_u_p_U.cpp:325
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: TwentyNodeBrick_u_p_U.cpp:358
Base class for the finite elements.
Definition: Element.h:109
int sendSelf(CommParameters &)
Send the object.
Definition: TwentyNodeBrick_u_p_U.cpp:398
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TwentyNodeBrick_u_p_U.cpp:238
Body forces over an element.
Definition: BodyForces3D.h:39
int update(void)
Updates the element state.
Definition: TwentyNodeBrick_u_p_U.cpp:535
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: TwentyNodeBrick_u_p_U.cpp:165
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: TwentyNodeBrick_u_p_U.cpp:161
Base class for 2D and 3D materials.
Definition: NDMaterial.h:97
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: TwentyNodeBrick_u_p_U.cpp:440
void zeroLoad(void)
Zeroes the loads over the element.
Definition: TwentyNodeBrick_u_p_U.cpp:265
Response * setResponse(const std::vector< std::string > &argv, Information &eleInfo)
setResponse() is a method invoked to determine if the element will respond to a request for a certain...
Definition: TwentyNodeBrick_u_p_U.cpp:412
twenty node exahedron.
Definition: TwentyNodeBrick_u_p_U.h:73
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
int recvSelf(const CommParameters &)
Receive the object.
Definition: TwentyNodeBrick_u_p_U.cpp:405