58 #ifndef TWENTYNODEBRICK_U_P_U_H 59 #define TWENTYNODEBRICK_U_P_U_H 63 #include <domain/mesh/element/ElemWithMaterial.h> 64 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h" 65 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h" 88 static const int Num_IntegrationPts;
89 static const int Num_TotalGaussPts;
90 static const int Num_Nodes;
91 static const int Num_Dim;
92 static const int Num_Dof;
93 static const int Num_ElemDof;
94 static const double pts[3];
95 static const double wts[3];
110 static BJtensor shapeFunction(
double,
double,
double);
111 static BJtensor shapeFunctionDerivative(
double,
double,
double);
119 BJtensor getStiffnessTensorKep(
void)
const;
120 BJtensor getStiffnessTensorG12(
void)
const;
121 BJtensor getStiffnessTensorP(
void)
const;
122 BJtensor getMassTensorMsf(
void)
const;
123 BJtensor getDampTensorC123(
void)
const;
124 const Matrix &getStiff(
int Ki_flag)
const;
125 double getPorePressure(
double,
double,
double)
const;
130 int node_numb_1,
int node_numb_2,
int node_numb_3,
int node_numb_4,
131 int node_numb_5,
int node_numb_6,
int node_numb_7,
int node_numb_8,
132 int node_numb_9,
int node_numb_10,
int node_numb_11,
int node_numb_12,
133 int node_numb_13,
int node_numb_14,
int node_numb_15,
int node_numb_16,
134 int node_numb_17,
int node_numb_18,
int node_numb_19,
int node_numb_20,
136 double nn,
double alf,
double rs,
double rf,
137 double permb_x,
double permb_y,
double permb_z,
138 double kks,
double kkf,
double pp);
152 const Matrix &getInitialStiff(
void)
const;
159 int addInertiaLoadToUnbalance(
const Vector &accel);
166 void Print(std::ostream &s,
int flag =0)
const;
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:94
Element * getCopy(void) const
Virtual constructor.
Definition: TwentyNodeBrick_u_p_U.cpp:150
Communication parameters between processes.
Definition: Communicator.h:66
Base class response objects.
Definition: Response.h:81
Boris Jeremic tensor class.
Definition: BJtensor.h:112
int recvSelf(const Communicator &)
Receive the object.
Definition: TwentyNodeBrick_u_p_U.cpp:417
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: TwentyNodeBrick_u_p_U.cpp:337
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: TwentyNodeBrick_u_p_U.cpp:370
Base class for the finite elements.
Definition: Element.h:112
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TwentyNodeBrick_u_p_U.cpp:238
Body forces over an element.
Definition: BodyForces3D.h:40
int update(void)
Updates the element state.
Definition: TwentyNodeBrick_u_p_U.cpp:547
void alive(void)
Reactivates the element.
Definition: TwentyNodeBrick_u_p_U.cpp:265
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 degrees of freedom.
Definition: TwentyNodeBrick_u_p_U.cpp:161
int sendSelf(Communicator &)
Send the object.
Definition: TwentyNodeBrick_u_p_U.cpp:410
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: TwentyNodeBrick_u_p_U.cpp:452
void zeroLoad(void)
Zeroes the loads over the element.
Definition: TwentyNodeBrick_u_p_U.cpp:277
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: TwentyNodeBrick_u_p_U.cpp:503
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:424
Twenty node hexahedral element for three-dimensional coupled problems.
Definition: TwentyNodeBrick_u_p_U.h:80
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
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117