28 #ifndef TOTALLAGRANGIANFD8BRICK_H 29 #define TOTALLAGRANGIANFD8BRICK_H 31 #include <domain/mesh/element/ElemWithMaterial.h> 32 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h" 33 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h" 54 static const double pts[2];
55 static const double wts[2];
61 double det_of_Jacobian;
64 static const int NumIntegrationPts;
65 static const int NumTotalGaussPts;
66 static const int NumNodes;
67 static const int NumDof;
68 static const int NumElemDof;
72 static BJtensor shapeFunction(
double ,
double ,
double );
73 static BJtensor shapeFunctionDerivative(
double ,
double ,
double );
76 BJtensor Jacobian_3D(
double ,
double ,
double)
const;
77 BJtensor Jacobian_3Dinv(
double ,
double ,
double)
const;
78 BJtensor dh_Global(
double ,
double ,
double)
const;
82 BJtensor getStiffnessTensor(
void)
const;
85 BJtensor getSurfaceForce(
void)
const;
90 int node_numb_1,
int node_numb_2,
int node_numb_3,
int node_numb_4,
92 int node_numb_5,
int node_numb_6,
int node_numb_7,
int node_numb_8,
113 const Matrix &getInitialStiff(
void)
const;
120 int addInertiaLoadToUnbalance(
const Vector &accel);
134 void Print(std::ostream &s,
int flag =0)
const;
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: TotalLagrangianFD8NodeBrick.cpp:1093
Float vector abstraction.
Definition: Vector.h:94
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: TotalLagrangianFD8NodeBrick.cpp:1020
Communication parameters between processes.
Definition: Communicator.h:66
Base class response objects.
Definition: Response.h:81
Element * getCopy(void) const
Virtual constructor.
Definition: TotalLagrangianFD8NodeBrick.cpp:91
virtual int sendSelf(Communicator &)
Send the object.
Definition: TotalLagrangianFD8NodeBrick.cpp:1003
Boris Jeremic tensor class.
Definition: BJtensor.h:112
int commitState()
Commit the current element state.
Definition: TotalLagrangianFD8NodeBrick.cpp:113
int update()
Updates the element state.
Definition: TotalLagrangianFD8NodeBrick.cpp:190
virtual int recvSelf(const Communicator &)
Receive the object.
Definition: TotalLagrangianFD8NodeBrick.cpp:1012
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: TotalLagrangianFD8NodeBrick.cpp:660
Base class for the finite elements.
Definition: Element.h:112
int revertToLastCommit()
Revert to the last committed state.
Definition: TotalLagrangianFD8NodeBrick.cpp:142
Body forces over an element.
Definition: BodyForces3D.h:40
void zeroLoad()
Zeroes the loads over the element.
Definition: TotalLagrangianFD8NodeBrick.cpp:810
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TotalLagrangianFD8NodeBrick.cpp:734
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: TotalLagrangianFD8NodeBrick.cpp:911
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
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: TotalLagrangianFD8NodeBrick.cpp:1063
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: TotalLagrangianFD8NodeBrick.cpp:100
int revertToStart()
Reverts the element to its initial state.
Definition: TotalLagrangianFD8NodeBrick.cpp:165
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
Total lagrangian formulation eight node hexahedral element for three-dimensional problems.
Definition: TotalLagrangianFD8NodeBrick.h:44
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: TotalLagrangianFD8NodeBrick.cpp:104
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: TotalLagrangianFD8NodeBrick.cpp:946