54 #ifndef TOTALLAGRANGIANFD20BRICK_H    55 #define TOTALLAGRANGIANFD20BRICK_H    57 #include <domain/mesh/element/ElemWithMaterial.h>    58 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h"    59 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"    78     static const double pts[3]; 
    79     static const double wts[3]; 
    82     double det_of_Jacobian;
    86     static const int  NumIntegrationPts; 
    87     static const int  NumTotalGaussPts;
    88     static const int  NumNodes; 
    89     static const int  NumDof; 
    90     static const int  NumElemDof;
    92     static BJtensor shapeFunction(
double , 
double , 
double );
    93     static BJtensor shapeFunctionDerivative(
double , 
double , 
double );
   101     BJtensor getStiffnessTensor(
void) 
const;
   104     BJtensor getSurfaceForce(
void) 
const;
   108     int node_numb_1,  
int node_numb_2,  
int node_numb_3,  
int node_numb_4,
   109     int node_numb_5,  
int node_numb_6,  
int node_numb_7,  
int node_numb_8,
   110     int node_numb_9,  
int node_numb_10, 
int node_numb_11, 
int node_numb_12,
   111     int node_numb_13, 
int node_numb_14, 
int node_numb_15, 
int node_numb_16,
   112     int node_numb_17, 
int node_numb_18, 
int node_numb_19, 
int node_numb_20,
   124     const Matrix &getInitialStiff(
void) 
const;
   129     int addInertiaLoadToUnbalance(
const Vector &accel);
   137     void Print(std::ostream &s, 
int flag =0) 
const;
 const Matrix & getMass(void) const
Returns the mass matrix. 
Definition: TotalLagrangianFD20NodeBrick.cpp:470
Float vector abstraction. 
Definition: Vector.h:94
Total lagrangian formulation twenty node hexahedral element for three-dimensional problems...
Definition: TotalLagrangianFD20NodeBrick.h:71
Communication parameters between processes. 
Definition: Communicator.h:66
Base class response objects. 
Definition: Response.h:81
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix. 
Definition: TotalLagrangianFD20NodeBrick.cpp:429
Boris Jeremic tensor class. 
Definition: BJtensor.h:112
int getNumDOF(void) const
return the number of DOF associated with the element. 
Definition: TotalLagrangianFD20NodeBrick.cpp:123
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: TotalLagrangianFD20NodeBrick.cpp:659
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element. 
Definition: TotalLagrangianFD20NodeBrick.cpp:564
void setDomain(Domain *theDomain)
Sets the domain for the element. 
Definition: TotalLagrangianFD20NodeBrick.cpp:127
Base class for the finite elements. 
Definition: Element.h:112
virtual int recvSelf(const Communicator &)
Receive the object. 
Definition: TotalLagrangianFD20NodeBrick.cpp:622
Body forces over an element. 
Definition: BodyForces3D.h:40
Element * getCopy(void) const
Virtual constructor. 
Definition: TotalLagrangianFD20NodeBrick.cpp:111
void Print(std::ostream &s, int flag=0) const
Print stuff. 
Definition: TotalLagrangianFD20NodeBrick.cpp:630
Base class for 2D and 3D materials. 
Definition: NDMaterial.h:101
void alive(void)
Reactivates the element. 
Definition: TotalLagrangianFD20NodeBrick.cpp:512
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces. 
Definition: TotalLagrangianFD20NodeBrick.cpp:584
int update()
Updates the element state. 
Definition: TotalLagrangianFD20NodeBrick.cpp:135
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
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis. 
Definition: TotalLagrangianFD20NodeBrick.cpp:685
virtual int sendSelf(Communicator &)
Send the object. 
Definition: TotalLagrangianFD20NodeBrick.cpp:615