53 #ifndef TWENTYSEVENNODEBRICK_H 54 #define TWENTYSEVENNODEBRICK_H 56 #include <domain/mesh/element/ElementBase.h> 57 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h" 73 double determinant_of_Jacobian;
75 int r_integration_order;
76 int s_integration_order;
77 int t_integration_order;
82 std::vector<MatPoint3D> matpoint;
107 int node_numb_1,
int node_numb_2,
int node_numb_3,
int node_numb_4,
108 int node_numb_5,
int node_numb_6,
int node_numb_7,
int node_numb_8,
109 int node_numb_9,
int node_numb_10,
int node_numb_11,
int node_numb_12,
110 int node_numb_13,
int node_numb_14,
int node_numb_15,
int node_numb_16,
111 int node_numb_17,
int node_numb_18,
int node_numb_19,
int node_numb_20,
112 int node_numb_21,
int node_numb_22,
int node_numb_23,
int node_numb_24,
113 int node_numb_25,
int node_numb_26,
int node_numb_27,
135 const Matrix &getInitialStiff(
void)
const;
138 const Matrix &getConsMass(
void)
const;
139 const Matrix &getLumpedMass(
void)
const;
144 int addInertiaLoadToUnbalance(
const Vector &accel);
145 const Vector FormEquiBodyForce(
void);
152 void Print(std::ostream &s,
int flag =0)
const;
160 void incremental_Update(
void);
163 static BJtensor H_3D(
double r1,
double r2,
double r3);
164 BJtensor interp_poli_at(
double r,
double s,
double t);
165 static BJtensor dh_drst_at(
double r,
double s,
double t);
182 BJtensor total_disp(FILE *fp,
double * u);
188 int get_global_number_of_node(
int local_node_number);
189 int get_Brick_Number(
void);
196 static double get_Gauss_p_c(
short order,
short point_numb);
197 static double get_Gauss_p_w(
short order,
short point_numb);
202 BJtensor iterative_nodal_forces(
void)
const;
207 BJtensor linearized_nodal_forces(
void)
const;
213 void reportshort(
char *);
214 void reportPAK(
char *);
215 void reportpqtheta(
int);
217 void computeGaussPoint(
void);
218 void reportCIPIC(
char *);
219 void reportTensorF(FILE *);
220 Vector getWeightofGP(
void);
int revertToStart(void)
Reverts the element to its initial state.
Definition: TwentySevenNodeBrick.cpp:3042
BJtensor mass_matrix(const BJtensor &)
Returns the mass matrix.
Definition: TwentySevenNodeBrick.cpp:1328
Float vector abstraction.
Definition: Vector.h:94
Twenty seven node hexahedral element for three-dimensional problems.
Definition: TwentySevenNodeBrick.h:70
int revertToLastCommit(void)
Revert to the last committed state.
Definition: TwentySevenNodeBrick.cpp:3017
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: TwentySevenNodeBrick.cpp:2880
Stress tensor.
Definition: stresst.h:70
Communication parameters between processes.
Definition: Communicator.h:66
virtual int sendSelf(Communicator &)
Send the object.
Definition: TwentySevenNodeBrick.cpp:3813
Base class response objects.
Definition: Response.h:81
Boris Jeremic tensor class.
Definition: BJtensor.h:112
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: TwentySevenNodeBrick.cpp:2876
BJtensor Nodal_Coordinates(void) const
Returns the coordinates of the nodes.
Definition: TwentySevenNodeBrick.cpp:1362
BJtensor getMassTensor(void) const
Returns the tensor de masas.
Definition: TwentySevenNodeBrick.cpp:1197
Base class for the finite elements.
Definition: Element.h:112
void alive(void)
Reactivates the element.
Definition: TwentySevenNodeBrick.cpp:3181
void set_strain_stress_tensor(FILE *fp, double *u)
Asigna el stiffness tensor.
Definition: TwentySevenNodeBrick.cpp:1108
Body forces over an element.
Definition: BodyForces3D.h:40
int commitState(void)
Commit the current element state.
Definition: TwentySevenNodeBrick.cpp:2887
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: TwentySevenNodeBrick.cpp:3668
BJtensor getStiffnessTensor(void) const
Returns the stiffness tensor.
Definition: TwentySevenNodeBrick.cpp:909
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TwentySevenNodeBrick.cpp:3144
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: TwentySevenNodeBrick.cpp:3067
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
Element * getCopy(void) const
Virtual constructor.
Definition: TwentySevenNodeBrick.cpp:195
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: TwentySevenNodeBrick.cpp:3958
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: TwentySevenNodeBrick.cpp:3828
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: TwentySevenNodeBrick.cpp:3643
BJtensor stiffness_matrix(const BJtensor &)
Returns the stiffness matrix.
Definition: TwentySevenNodeBrick.cpp:1298
virtual int recvSelf(const Communicator &)
Receive the object.
Definition: TwentySevenNodeBrick.cpp:3820
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Base class for finite element with pointer to nodes container.
Definition: ElementBase.h:47
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: TwentySevenNodeBrick.cpp:3896
int update(void)
Updates the element state.
Definition: TwentySevenNodeBrick.cpp:4940
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
TwentySevenNodeBrick(void)
Constructor.
Definition: TwentySevenNodeBrick.cpp:187
~TwentySevenNodeBrick()
Destructor.
Definition: TwentySevenNodeBrick.cpp:199