57 #include <domain/mesh/element/volumetric/BrickBase.h>    58 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"    59 #include "domain/mesh/element/utils/fvectors/FVectorBrick.h"    69     static const int numberGauss= 8; 
    70     static const int nShape = 4;
    78     void shape_functions_loop(
void) 
const;
    88     static double gaussPoint[numberGauss][
ndm]; 
    89     static double dvol[numberGauss]; 
    91     static double Shape[nShape][
numberNodes][numberGauss]; 
    94     static const double wg[numberGauss];
   101     void formInertiaTerms(
int tangFlag) 
const;
   103     void formResidAndTangent(
int tang_flag) 
const;
   106     const Matrix& computeB( 
int node, 
const double shp[4][8]) 
const;
   108     static size_t getVectorIndex(
const size_t &,
const size_t &);
   118     Brick( 
int tag, 
int node1,
int node2,
int node3,
int node4,
int node5,
int node6,
int node7,
int node8, 
NDMaterial &theMaterial, 
const BodyForces3D &bf= 
BodyForces3D());
   139     double getAvgStress(
const size_t &,
const size_t &) 
const;
   141     double getAvgStrain(
const size_t &,
const size_t &) 
const;
   159     void Print( std::ostream &s, 
int flag ) 
const;
 void setDomain(Domain *theDomain)
set domain 
Definition: Brick.cpp:138
Float vector abstraction. 
Definition: Vector.h:94
virtual ~Brick(void)
destructor 
Definition: Brick.cpp:127
const Matrix & getMass(void) const
Return mass matrix. 
Definition: Brick.cpp:348
int update(void)
Form residual and tangent. 
Definition: Brick.cpp:619
Vector getAvgStrain(void) const
Return the average strain in the element. 
Definition: Brick.cpp:179
Brick(void)
Default constructor. 
Definition: Brick.cpp:98
Communication parameters between processes. 
Definition: Communicator.h:66
Base class for hexahedral elements. 
Definition: BrickBase.h:69
Base class response objects. 
Definition: Response.h:81
Vector getAvgStress(void) const
Return the average stress in the element. 
Definition: Brick.cpp:148
const Matrix & getTangentStiff(void) const
Return stiffness matrix. 
Definition: Brick.cpp:192
Response * setResponse(const std::vector< std::string > &argv, Information &eleInformation)
setResponse() is a method invoked to determine if the element will respond to a request for a certain...
Definition: Brick.cpp:947
int getNumDOF(void) const
return number of dofs 
Definition: Brick.cpp:144
Eight node hexahedral element for three-dimensional problems. 
Definition: Brick.h:65
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument. 
Definition: Brick.cpp:933
Base class for the finite elements. 
Definition: Element.h:112
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis. 
Definition: Brick.cpp:970
virtual int sendSelf(Communicator &)
Sends object through the communicator argument. 
Definition: Brick.cpp:919
void alive(void)
Reactivates the element. 
Definition: Brick.cpp:358
static const int ndf
Number of DOFs per node. 
Definition: BrickBase.h:74
Ingernal forces for a brick element. 
Definition: FVectorBrick.h:40
Body forces over an element. 
Definition: BodyForces3D.h:40
const Matrix & getInitialStiff(void) const
Return initial stiffness matrix. 
Definition: Brick.cpp:265
Vector getShapeFunctionsValues(const double &r, const double &s, const double &t) const
Return the values of the shape functions. 
Definition: Brick.cpp:216
const Vector & getResistingForceIncInertia(void) const
Get residual with inertia terms. 
Definition: Brick.cpp:485
void zeroLoad(void)
Remove element loads. 
Definition: Brick.cpp:370
virtual void createInertiaLoad(const Vector &)
Creates the inertia load that corresponds to the acceleration argument. 
Definition: Brick.cpp:581
int sendData(Communicator &comm)
Send members through the communicator argument. 
Definition: Brick.cpp:899
static const int numberNodes
Number of nodes. 
Definition: BrickBase.h:72
Matrix getGaussPointsPositions(void) const
Return the coordinates of the Gauss points. 
Definition: Brick.cpp:203
void Print(std::ostream &s, int flag) const
Print out element data. 
Definition: Brick.cpp:999
Base class for 2D and 3D materials. 
Definition: NDMaterial.h:101
static const int ndm
Space dimension. 
Definition: BrickBase.h:73
int recvData(const Communicator &comm)
Receives members through the communicator argument. 
Definition: Brick.cpp:909
Open source finite element program for structural analysis. 
Definition: ContinuaReprComponent.h:35
Matrix of floats. 
Definition: Matrix.h:111
Element * getCopy(void) const
Virtual constructor. 
Definition: Brick.cpp:122
Base class for loads over elements. 
Definition: ElementalLoad.h:79
int addInertiaLoadToUnbalance(const Vector &accel)
Add inertia load due to acceleration argument. 
Definition: Brick.cpp:438
Domain (mesh and boundary conditions) of the finite element model. 
Definition: Domain.h:117
const Vector & getResistingForce(void) const
Get residual. 
Definition: Brick.cpp:469