56 #include <domain/mesh/element/volumetric/BrickBase.h> 57 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h" 73 static const double sg[2] ;
74 static const double wg[8] ;
79 static double xl[][8] ;
84 void formInertiaTerms(
int tangFlag )
const;
85 void formResidAndTangent(
int tang_flag )
const;
86 void computeBasis(
void)
const;
87 const Matrix& computeBbar(
int node,
const double shp[4][8],
const double shpBar[4][8])
const;
88 Matrix transpose(
int dim1,
int dim2,
const Matrix &M ) ;
117 void Print( std::ostream &s,
int flag ) ;
121 const Matrix &getInitialStiff(
void)
const;
125 int addInertiaLoadToUnbalance(
const Vector &accel);
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: BbarBrick.cpp:135
Float vector abstraction.
Definition: Vector.h:93
virtual int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: BbarBrick.cpp:958
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: BbarBrick.cpp:977
Base class for hexahedra.
Definition: BrickBase.h:44
Definition: Response.h:71
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: BbarBrick.cpp:926
void Print(std::ostream &s, int flag)
Print stuff.
Definition: BbarBrick.cpp:140
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: BbarBrick.cpp:935
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: BbarBrick.cpp:1000
Base class for the finite elements.
Definition: Element.h:109
Hexaedro.
Definition: BbarBrick.h:63
Body forces over an element.
Definition: BodyForces3D.h:39
const Vector & getResistingForceIncInertia(void) const
Get residual with inertia terms.
Definition: BbarBrick.cpp:396
Base class for 2D and 3D materials.
Definition: NDMaterial.h:97
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: BbarBrick.cpp:331
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
Element * getCopy(void) const
Virtual constructor.
Definition: BbarBrick.cpp:112
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: BbarBrick.cpp:128
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: BbarBrick.cpp:161
Base class for loads over elements.
Definition: ElementalLoad.h:77
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: BbarBrick.cpp:383
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:116
virtual int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: BbarBrick.cpp:944