62 #ifndef EIGHTNODEBRICK_U_P_U_H    63 #define EIGHTNODEBRICK_U_P_U_H    65 #include <domain/mesh/element/volumetric/BrickBase.h>    66 #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[2]; 
    95     static const double wts[2]; 
   111                    int node_numb_1, 
int node_numb_2, 
int node_numb_3, 
int node_numb_4,
   112                    int node_numb_5, 
int node_numb_6, 
int node_numb_7, 
int node_numb_8,
   114                    double nn, 
double alf, 
double rs,
double rf,
   115                    double permb_x, 
double permb_y, 
double permb_z,
   116                    double kks, 
double kkf, 
double pp);
   129     const Matrix &getInitialStiff(
void) 
const;
   136     int addInertiaLoadToUnbalance(
const Vector &accel);
   143     void Print(std::ostream &s, 
int flag =0) 
const;
   152     BJtensor shapeFunction(
double, 
double, 
double) 
const;
   153     BJtensor shapeFunctionDerivative(
double, 
double, 
double) 
const;
   161     BJtensor getStiffnessTensorKep() 
const;
   162     BJtensor getStiffnessTensorG12() 
const;
   163     BJtensor getStiffnessTensorP() 
const;
   165     BJtensor getDampTensorC123(
void) 
const;
   166     const Matrix &getStiff(
int Ki_flag) 
const;
   167     double getPorePressure(
double, 
double, 
double);
 const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix. 
Definition: EightNodeBrick_u_p_U.cpp:159
Float vector abstraction. 
Definition: Vector.h:94
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis. 
Definition: EightNodeBrick_u_p_U.cpp:440
int sendSelf(Communicator &)
Send the object. 
Definition: EightNodeBrick_u_p_U.cpp:399
int update(void)
Updates the element state. 
Definition: EightNodeBrick_u_p_U.cpp:523
Element * getCopy(void) const
Virtual constructor. 
Definition: EightNodeBrick_u_p_U.cpp:136
Communication parameters between processes. 
Definition: Communicator.h:66
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element. 
Definition: EightNodeBrick_u_p_U.cpp:324
Base class for hexahedral elements. 
Definition: BrickBase.h:69
Base class response objects. 
Definition: Response.h:81
void alive(void)
Reactivates the element. 
Definition: EightNodeBrick_u_p_U.cpp:252
const Matrix & getDamp(void) const
Returns the damping matrix. 
Definition: EightNodeBrick_u_p_U.cpp:179
Boris Jeremic tensor class. 
Definition: BJtensor.h:112
Base class for the finite elements. 
Definition: Element.h:112
void zeroLoad(void)
Zeroes loads on element. 
Definition: EightNodeBrick_u_p_U.cpp:264
Body forces over an element. 
Definition: BodyForces3D.h:40
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: EightNodeBrick_u_p_U.cpp:413
void Print(std::ostream &s, int flag=0) const
Print stuff. 
Definition: EightNodeBrick_u_p_U.cpp:491
int recvSelf(const Communicator &)
Receive the object. 
Definition: EightNodeBrick_u_p_U.cpp:406
Base class for 2D and 3D materials. 
Definition: NDMaterial.h:101
const Matrix & getMass(void) const
Returns the mass matrix. 
Definition: EightNodeBrick_u_p_U.cpp:223
int addLoad(ElementalLoad *theLoad, double loadFactor)
Adds to the element the load being passed as parameter. 
Definition: EightNodeBrick_u_p_U.cpp:272
Eight node hexahedral element for three-dimensional coupled problems. 
Definition: EightNodeBrick_u_p_U.h:80
Open source finite element program for structural analysis. 
Definition: ContinuaReprComponent.h:35
Matrix of floats. 
Definition: Matrix.h:111
void setDomain(Domain *theDomain)
Sets the domain for the element. 
Definition: EightNodeBrick_u_p_U.cpp:151
Base class for loads over elements. 
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model. 
Definition: Domain.h:117
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces. 
Definition: EightNodeBrick_u_p_U.cpp:358
int getNumDOF(void) const
return the number of DOF associated with the element. 
Definition: EightNodeBrick_u_p_U.cpp:147