xc
TotalLagrangianFD8NodeBrick.h
1 //===============================================================================
2 //# COPYRIGHT (C): Woody's license (by BJ):
3 // ``This source code is Copyrighted in
4 // U.S., for an indefinite period, and anybody
5 // caught using it without our permission, will be
6 // mighty good friends of ourn, cause we don't give
7 // a darn. Hack it. Compile it. Debug it. Run it.
8 // Yodel it. Enjoy it. We wrote it, that's all we
9 // wanted to do.''
10 //
11 //# PROJECT: Object Oriented Finite Element Program
12 //# PURPOSE: Finite Deformation Hyper-Elastic classes
13 //# CLASS:
14 //#
15 //# VERSION: 0.6_(1803398874989) (golden section)
16 //# LANGUAGE: C++
17 //# TARGET OS: all...
18 //# DESIGN: Zhao Cheng, Boris Jeremic (jeremic@ucdavis.edu)
19 //# PROGRAMMER(S): Zhao Cheng, Boris Jeremic
20 //#
21 //#
22 //# DATE: Sept2005
23 //# UPDATE HISTORY:
24 //#
25 //#
26 //===============================================================================
27 
28 #ifndef TOTALLAGRANGIANFD8BRICK_H
29 #define TOTALLAGRANGIANFD8BRICK_H
30 
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"
34 
35 
36 namespace XC {
37 
44 class TotalLagrangianFD8NodeBrick: public ElemWithMaterial<8,NDMaterialPhysicalProperties>
45  {
46  private:
47 
48 
49  static Matrix K; // Element stiffness Matrix
50 // static Matrix C; // Element damping matrix
51 
52  static Matrix M; // Element mass matrix
53  static Vector P; // Element resisting force vector
54  static const double pts[2]; // Stores quadrature points
55  static const double wts[2]; // Stores quadrature weights
56 
57  Vector *Q; // Applied nodal loads
58 
59  BodyForces3D bf; // Body forces
60 
61  double det_of_Jacobian;
62  mutable Matrix *Ki;
63 
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;
69 
70 
71 
72  static BJtensor shapeFunction(double , double , double );
73  static BJtensor shapeFunctionDerivative(double , double , double );
74 
75 
76  BJtensor Jacobian_3D(double , double , double) const;
77  BJtensor Jacobian_3Dinv(double , double , double) const;
78  BJtensor dh_Global(double , double , double) const;
79  BJtensor getNodesCrds(void) const;
80  BJtensor getNodesDisp(void) const;
81 
82  BJtensor getStiffnessTensor(void) const;
83  BJtensor getRtensor(void) const;
84  BJtensor getBodyForce(void) const;
85  BJtensor getSurfaceForce(void) const;
86  BJtensor getForces(void) const;
87  public:
89 
90  int node_numb_1, int node_numb_2, int node_numb_3, int node_numb_4,
91 
92  int node_numb_5, int node_numb_6, int node_numb_7, int node_numb_8,
93 
94  NDMaterial &m, const BodyForces3D &bForces= BodyForces3D());
95 
96 
97 
99  Element *getCopy(void) const;
101 
102 
103  int getNumDOF(void) const;
104 
105  void setDomain(Domain *theDomain);
106 
107  int commitState();
108  int revertToLastCommit ();
109  int revertToStart ();
110  int update();
111 
112  const Matrix &getTangentStiff(void) const;
113  const Matrix &getInitialStiff(void) const;
114  const Matrix &getMass(void) const;
115 
116 
117 
118  void zeroLoad ();
119  int addLoad(ElementalLoad *theLoad, double loadFactor);
120  int addInertiaLoadToUnbalance(const Vector &accel);
121 
122 
123 
124  const Vector &getResistingForce(void) const;
125  const Vector &getResistingForceIncInertia(void) const;
126 
127 
128 
129  // public methods for element output
130 
131  virtual int sendSelf(Communicator &);
132  virtual int recvSelf(const Communicator &);
133 
134  void Print(std::ostream &s, int flag =0) const;
135 
136 
137 
138  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
139  int getResponse(int responseID, Information &eleInformation);
140 
141 
142 // int setParameter(const char **argv, int argc, Information &info);
143 // int updateParameter(int parameterID, Information &info);
144  };
145 
146 } // end of XC namespace
147 
148 #endif
149 
150 
151 
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: TotalLagrangianFD8NodeBrick.cpp:1093
Float vector abstraction.
Definition: Vector.h:94
Information about an element.
Definition: Information.h:81
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