XC Open source finite element analysis program
EightNodeBrick.h
1 //----------------------------------------------------------------------------
2 // XC program; finite element analysis code
3 // for structural analysis and design.
4 //
5 // Copyright (C) Luis Claudio Pérez Tato
6 //
7 // This program derives from OpenSees <http://opensees.berkeley.edu>
8 // developed by the «Pacific earthquake engineering research center».
9 //
10 // Except for the restrictions that may arise from the copyright
11 // of the original program (see copyright_opensees.txt)
12 // XC is free software: you can redistribute it and/or modify
13 // it under the terms of the GNU General Public License as published by
14 // the Free Software Foundation, either version 3 of the License, or
15 // (at your option) any later version.
16 //
17 // This software is distributed in the hope that it will be useful, but
18 // WITHOUT ANY WARRANTY; without even the implied warranty of
19 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 // GNU General Public License for more details.
21 //
22 //
23 // You should have received a copy of the GNU General Public License
24 // along with this program.
25 // If not, see <http://www.gnu.org/licenses/>.
26 //----------------------------------------------------------------------------
28 //
29 // COPYLEFT (C): :-))
30 //``This source code is Copyrighted in U.S., by the The Regents of the University
31 //of California, for an indefinite period, and anybody caught using it without our
32 //permission, will be mighty good friends of ourn, cause we don't give a darn.
33 //Hack it. Compile it. Debug it. Run it. Yodel it. Enjoy it. We wrote it, that's
34 //all we wanted to do.'' bj
35 // PROJECT: Object Oriented Finite Element Program
36 // FILE: EightNodeBrick.h
37 // CLASS: EightNodeBrick
38 // MEMBER FUNCTIONS:
39 //
40 // MEMBER VARIABLES
41 //
42 // PURPOSE: Finite Element Class
43 // RETURN:
44 // VERSION:
45 // LANGUAGE: C++.ver >= 3.0
46 // TARGET OS: DOS || UNIX || . . .
47 // DESIGNER: Boris Jeremic, Zhaohui Yang and Xiaoyan Wu
48 // PROGRAMMER: Boris Jeremic, Zhaohui Yang and Xiaoyan Wu
49 // DATE: Aug. 2000
50 // UPDATE HISTORY: Modified from Brick3D and FourNodeQuad.hh 07/06/00
51 // Sept. - Oct 2000 connected to OpenSees by Zhaohui
52 // Sept 2001 optimized to some extent (static tensors...)
53 //
55 //
56 
57 
58 #ifndef EIGHTNODEBRICK_H
59 #define EIGHTNODEBRICK_H
60 
61 // Commented by Xiaoyan. Use ~/fem/node.hh 08/04/00
62 // Released Node.h now. Wu use Opensees's Node.09/27/00
63 
64 
65 #include <domain/mesh/element/ElementBase.h>
66 #include <utility/matrix/Matrix.h>
67 #include <utility/matrix/Vector.h>
68 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
69 
70 namespace XC {
71 class Node;
72 class NDMaterial;
73 class MatPoint3D;
74 class stresstensor;
75 class Information;
76 class BJtensor;
77 //class QuadRule1d;
78 
80 //
82 class EightNodeBrick: public ElementBase<8>
83  {
84  private:
85  int numDOF;
86 
87  mutable Matrix *Ki;
88 
89  static Matrix K;
90  static Matrix C;
91  static Matrix M;
92  static Vector P;
93  BodyForces3D bf;
94 
95  //Vector q0; //!< hold element load affects q0 and p0 in one vector
96 
97  double rho;
98  double pressure;
99  int order;
100 
101  //int dir; //!< Direction of vertial coord.
102  //double surflevel; //!< free surface level above or below this element
103 
104  //QuadRule1d *theQuadRule; //!< Integration rule
105 
106 // Matrix J; //!< Jacobian of transformation
107 // Matrix L; //!< Inverse of J
108 // Matrix B; //!< Strain interpolation matrix
109 
110 
111  // // static data - single copy for all objects of the class
112  // static G3Matrix N; //!< Displacement interpolation matrix
113 
114  // // private member functions - only objects of this class can call these
115  // void setJacobian (double r, double s, double t); //!< Xiaoyan changed
116  // double formDetJ (double r, double s, double t); //!< xi, eta to r and s
117  // void formBMatrix (double r, double s, double t); //!< and added t
118  // static void formNMatrix (double r, double s, double t); //!< 07/06/00
119 
120  private:
121  // element number (tag)
122  //unsigned int elem_numb;
123 
124  double determinant_of_Jacobian;
125  //int G_N_numbs[8]; //!< Global node numbers for this element Xiaoyan changed from 20 to 8
126 
127  NDMaterial *mmodel;
128 
129  int r_integration_order;
130  int s_integration_order;
131  int t_integration_order;
132 
133  // Now I want 3D array of Material points!
134  // MatPoint3D[r_integration_order][s_integration_order][t_integration_order]
135  // 3D array of Material points
136  std::vector<MatPoint3D> matpoint;
137 
138  // 3D array of material models for each Material points
139  // NDMaterial *GPmmodel; //!< pointer to array of material models for Material Points
140  // Do we need this one?
141 
142  //..NDMaterial *MatPoint; //!< Zhaohui 10-01-2000
143 
144 
145  int LM[24];
146  public:
147  EightNodeBrick(int element_number,
148  int node_numb_1, int node_numb_2, int node_numb_3, int node_numb_4,
149  int node_numb_5, int node_numb_6, int node_numb_7, int node_numb_8,
150  NDMaterial * Globalmmodel, const BodyForces3D &bForces,
151  double r, double p);
152  // int dir, double surflevel);
153  //, EPState *InitEPS); const std::string &type,
154 
155  EightNodeBrick(void);
156  Element *getCopy(void) const;
157  ~EightNodeBrick(void);
158 
159  //Not needed Zhaohui
160  //int Initialize(int element_number,
161  // int node_numb_1, int node_numb_2, int node_numb_3, int node_numb_4,
162  // int node_numb_5, int node_numb_6, int node_numb_7, int node_numb_8,
163  // NDMaterial * Globalmmodel, const BodyForces3D &=,
164  // double p, double r);
165  // //, EPState * InitEPS);const std::string &type,
166 
167  // update, Guanzhou added Apr. 2004 to update incremental strain in the domain
168  int update(void);
169 
170  int getNumDOF(void) const;
171  void setDomain(Domain *theDomain);
172 
173  // public methods to set the state of the element
174  int commitState();
175  int revertToLastCommit ();
176  int revertToStart();
177 
178  // public methods to obtain stiffness, mass, damping and residual information
179  // We haven't build the following functions.
180  // All the value of K M Dmp and F are nothing. just
181  // want to test the program. Xiaoyan 08/16/00
182  const Matrix &getTangentStiff(void) const;
183  const Matrix &getInitialStiff(void) const;
184  const Matrix &getMass(void) const;
185 
186  const Matrix &getConsMass(void);
187 
188  int addLoad(ElementalLoad *theLoad, double loadFactor);
189  int addInertiaLoadToUnbalance(const Vector &accel);
190 
191  const Vector FormEquiBodyForce(void);
192  const Vector &getResistingForce(void) const;
193  const Vector &getResistingForceIncInertia(void) const;
194 
195  int sendSelf(CommParameters &);
196  int recvSelf(const CommParameters &);
197 
198  void Print(std::ostream &s, int flag =0);
199  // Do nothing with void Print (std::ostream &s, int flag =0);
200  // use Brick3D report. 08/16/00
201  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
202  int getResponse (int responseID, Information &eleInformation);
203 
204  void incremental_Update(void);
205  //void iterative_Update(void);
206 
207 
208  BJtensor H_3D(double r1, double r2, double r3) const;
209  BJtensor interp_poli_at(double r, double s, double t);
210  BJtensor dh_drst_at(double r, double s, double t) const;
211 
212 
213  //CE Dynamic Allocation for for brick3d s.
214  //Finite_Element * new_el( int total );
215  EightNodeBrick &operator[](int subscript);
216  //Finite_Element & operator[](int subscript);
217  //Finite_Element & operator[](int subscript);
218 
219  BJtensor getStiffnessTensor(void) const;
220  //matrix stiffness_tensor(void);
221 
222  void set_strain_stress_tensor(FILE *fp, double * u);
223  BJtensor getMassTensor(void) const;
224 //out19Jan2001 double Potential_Energy(void);
225 
226  BJtensor Jacobian_3D(BJtensor dh) const;
227  BJtensor Jacobian_3Dinv(BJtensor dh) const;
228  BJtensor Nodal_Coordinates(void) const;
229 
230  BJtensor incr_disp(void) const;
231  BJtensor total_disp(void) const;
232 
233  BJtensor total_disp(FILE *fp, double * u);
234 
235  BJtensor stiffness_matrix(const BJtensor & K);
236  BJtensor mass_matrix(const BJtensor & M);
237 
238 
239  int get_global_number_of_node(int local_node_number);
240  int get_Brick_Number(void);
241 
242  int * get_LM(void);
243  //void set_LM(Node * node); //!< commented out temporarily 09-27-2000 Zhaohui
244 
245  //these two files are originally in fe.h
246  double get_Gauss_p_c(short order, short point_numb) const;
247  double get_Gauss_p_w(short order, short point_numb) const;
248 
249  // returns nodal forces for given stress field in an element
250  BJtensor nodal_forces(void) const;
251  // returns nodal forces for ITERATIVE stress field in an element
252  BJtensor iterative_nodal_forces(void) const;
253  // returns nodal forces for given constant stress field in the element
254  BJtensor nodal_forces_from_stress(stresstensor &) const;
255  // returns nodal forces for given incremental strain field in an element
256  // by using the linearized constitutive BJtensor from the begining of the step !
257  BJtensor linearized_nodal_forces(void) const;
258 
259  // updates Material point stresses and strains from given displacements
260  BJtensor update_stress_strain(BJtensor & disp);
261 
262  void report(char *);
263  void reportshort(char *);
264  void reportPAK(char *);
265  void reportpqtheta(int);
266  //void reportLM(char *);
267  void computeGaussPoint(void);
268  void reportCIPIC(char *);
269  void reportTensorF(FILE *);
270 
271  // Setting initial E according to the initial pressure
272  //void setInitE(void);
273  //void reportStressTensorF(FILE *);
274  };
275 } // end of XC namespace
276 
277 
278 #endif
279 
BJtensor getMassTensor(void) const
Returns the tensor de masas.
Definition: EightNodeBrick.cpp:1075
Definition: BJtensor.h:110
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
Base class for 2D and 3D materials.
Definition: NDMaterial.h:91
Definition: Vector.h:82
void Print(std::ostream &s, int flag=0)
Imprime el objeto.
Definition: EightNodeBrick.cpp:3062
Hexaedro de ocho nodos.
Definition: EightNodeBrick.h:82
Information about an element.
Definition: Information.h:80
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
Element * getCopy(void) const
Virtual constructor.
Definition: EightNodeBrick.cpp:408
Definition: stresst.h:68
BJtensor dh_drst_at(double r, double s, double t) const
Definition: EightNodeBrick.cpp:690
Definition: Matrix.h:82
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: EightNodeBrick.cpp:2406
void incremental_Update(void)
Definition: EightNodeBrick.cpp:422
Base class for finite element with pointer to nodes container.
Definition: ElementBase.h:46
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: EightNodeBrick.cpp:2976
Body forces over an element.
Definition: BodyForces3D.h:39
int update(void)
Actualiza el estado of the element.
Definition: EightNodeBrick.cpp:4029
int commitState()
Consuma el estado of the element.
Definition: EightNodeBrick.cpp:2413
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: EightNodeBrick.cpp:2711
Communication parameters between processes.
Definition: CommParameters.h:65
================================================================================
Definition: ContinuaReprComponent.h:34
Definition: Response.h:71