XC Open source finite element analysis program
TwentyNodeBrick.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 // COPYRIGHT (C): :-))
30 // PROJECT: Object Oriented Finite Element Program
31 // FILE: TwentyNodeBrick.h
32 // CLASS: TwentyNodeBrick
33 // MEMBER FUNCTIONS:
34 //
35 // MEMBER VARIABLES
36 //
37 // PURPOSE: Finite Element Class
38 // RETURN:
39 // VERSION:
40 // LANGUAGE: C++.ver >= 3.0
41 // TARGET OS: DOS || UNIX || . . .
42 // DESIGNER: Boris Jeremic, Zhaohui Yang and Xiaoyan Wu
43 // PROGRAMMER: Boris Jeremic, Zhaohui Yang and Xiaoyan Wu
44 // DATE: Aug. 2001
45 // UPDATE HISTORY:
46 //
48 //
49 
50 
51 
52 #ifndef TWENTYNODEBRICK_H
53 #define TWENTYNODEBRICK_H
54 
55 #include <domain/mesh/element/ElementBase.h>
56 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
57 
58 
59 
60 namespace XC {
61 class Node;
62  class MatPoint3D;
63  class BJtensor;
64  class NDMaterial;
65 
67 //
69 class TwentyNodeBrick: public ElementBase<20>
70  {
71  private:
72  // private attributes - a copy for each object of the class
73 
74  int numDOF; // Number of element DOF
75 
76  Matrix *Ki;
77 
78  static Matrix K; // Element stiffness Matrix
79  static Matrix C; // Element damping matrix
80  static Matrix M; // Element mass matrix
81  static Vector P; // Element resisting force vector
82  BodyForces3D bf;
83 
84  // double thickness; // Element thickness
85  double rho; // Mass per unit volume DO WE GET THIS ONE OUT!!!
86  double pressure; // Normal surface traction (pressure) over entire element
87  int order; // Order of the quadrature rule
88 
89  //Matrix J; // Jacobian of transformation
90  //Matrix L; // Inverse of J
91  //Matrix B; // Strain interpolation matrix
92 
93 
94  private:
95 
96  double determinant_of_Jacobian;
97 
98  NDMaterial * mmodel; // pointer to GLOBAL material models
99 
100  int r_integration_order; // Gauss-Legendre integration order in r direction
101  int s_integration_order; // Gauss-Legendre integration order in s direction
102  int t_integration_order; // Gauss-Legendre integration order in t direction
103 
104  // Now I want 3D array of Material points!
105  // MatPoint3D[r_integration_order][s_integration_order][t_integration_order]
106  // 3D array of Material points
107  MatPoint3D ** matpoint; // pointer to array of Material Points
108 
109  // this is LM array. This array holds DOFs for this element
110  //int LM[60]; // for 20noded x 3 = 60
111  public:
112 
113  void incremental_Update(void);
114  //void iterative_Update(void);
115 
116  static BJtensor H_3D(double r1, double r2, double r3);
117  BJtensor interp_poli_at(double r, double s, double t);
118  static BJtensor dh_drst_at(double r, double s, double t);
119 
120 
121  TwentyNodeBrick & operator[](int subscript);
122 
123  BJtensor getStiffnessTensor(void) const;
124 
125  void set_strain_stress_tensor(FILE *fp, double * u);
126  BJtensor getMassTensor(void) const;
127 
128  BJtensor Jacobian_3D(BJtensor dh) const;
129  BJtensor Jacobian_3Dinv(BJtensor dh) const;
130  BJtensor Nodal_Coordinates(void) const;
131 
132  BJtensor incr_disp(void) const;
133  BJtensor total_disp(void) const;
134 
135  BJtensor total_disp(FILE *fp, double * u);
136 
137  BJtensor stiffness_matrix(const BJtensor & K);
138  BJtensor mass_matrix(const BJtensor & M);
139 
140 
141  int get_global_number_of_node(int local_node_number);
142  int get_Brick_Number(void);
143 
144 
145  //int * get_LM(void);
146  //void set_LM(Node * node); // commented out temporarily 09-27-2000 Zhaohui
147 
148  //these two files are originally in fe.h
149  static double get_Gauss_p_c(short order, short point_numb);
150  static double get_Gauss_p_w(short order, short point_numb);
151 
152  // returns nodal forces for given stress field in an element
153  BJtensor nodal_forces(void) const;
154  // returns nodal forces for ITERATIVE stress field in an element
155  BJtensor iterative_nodal_forces(void) const;
156  // returns nodal forces for given constant stress field in the element
157  BJtensor nodal_forces_from_stress(stresstensor &) const;
158  // returns nodal forces for given incremental strain field in an element
159  // by using the linearized constitutive tensor from the begining of the step !
160  BJtensor linearized_nodal_forces(void) const;
161 
162  // updates Material point stresses and strains from given displacements
163  BJtensor update_stress_strain(BJtensor & disp);
164 
165  void report(const std::string &);
166  void reportshort(const std::string &);
167  void reportPAK(const std::string &);
168  void reportpqtheta(int);
169  //void reportLM(const std::string &);
170  void computeGaussPoint(void);
171  void reportCIPIC(const std::string &);
172  void reportTensorF(FILE *);
173  Vector getWeightofGP(void);
174 
175 
176  // Setting initial E according to the initial pressure
177  //void setInitE(void);
178  //void reportStressTensorF(FILE *);
179  public:
180  TwentyNodeBrick(int element_number,
181  int node_numb_1, int node_numb_2, int node_numb_3, int node_numb_4,
182  int node_numb_5, int node_numb_6, int node_numb_7, int node_numb_8,
183  int node_numb_9, int node_numb_10, int node_numb_11, int node_numb_12,
184  int node_numb_13, int node_numb_14, int node_numb_15, int node_numb_16,
185  int node_numb_17, int node_numb_18, int node_numb_19, int node_numb_20,
186  NDMaterial * Globalmmodel, const BodyForces3D &bForces, double r, double p);
187 
188  TwentyNodeBrick(void);
189  Element *getCopy(void) const;
190  ~TwentyNodeBrick(void);
191 
192  int getNumDOF(void) const;
193  void setDomain(Domain *theDomain);
194 
195  // public methods to set the state of the element
196  int commitState ();
197  int revertToLastCommit ();
198  int revertToStart ();
199 
200  // update, Guanzhou added Apr. 2004 to update incremental strain in the domain
201  int update(void);
202 
203  // public methods to obtain stiffness, mass, damping and residual information
204  // We haven't build the following functions.
205  // All the value of K M Dmp and F are nothing.
206  const Matrix &getTangentStiff(void) const;
207  const Matrix &getInitialStiff(void) const;
208  const Matrix &getMass(void) const;
209 
210  const Matrix &getConsMass(void) const;
211  const Matrix &getLumpedMass(void) const;
212 
213  int addLoad(ElementalLoad *theLoad, double loadFactor);
214  //int addLoad(const Vector &addP);
215  int addInertiaLoadToUnbalance(const Vector &accel);
216  const Vector FormEquiBodyForce(void);
217  const Vector &getResistingForce(void) const;
218  const Vector &getResistingForceIncInertia(void) const;
219 
220  virtual int sendSelf(CommParameters &);
221  virtual int recvSelf(const CommParameters &);
222 
223  void Print(std::ostream &s, int flag =0);
224  // Do nothing with void Print (std::ostream &s, int flag =0);
225  // use Brick3D report. 08/16/00
226  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
227  int getResponse(int responseID, Information &eleInformation);
228 };
229 } // end of XC namespace
230 
231 
232 #endif
233 
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
Hexaedro de veinte nodos.
Definition: TwentyNodeBrick.h:69
Information about an element.
Definition: Information.h:80
static BJtensor dh_drst_at(double r, double s, double t)
Definition: TwentyNodeBrick.cpp:510
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
int commitState()
Consuma el estado of the element.
Definition: TwentyNodeBrick.cpp:2509
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: TwentyNodeBrick.cpp:3199
Definition: stresst.h:68
Definition: Matrix.h:82
int update(void)
Actualiza el estado of the element.
Definition: TwentyNodeBrick.cpp:4442
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: TwentyNodeBrick.cpp:2502
Base class for finite element with pointer to nodes container.
Definition: ElementBase.h:46
Body forces over an element.
Definition: BodyForces3D.h:39
Element * getCopy(void) const
Virtual constructor.
Definition: TwentyNodeBrick.cpp:208
Integration point on three-dimensional space.
Definition: MatPoint3D.h:66
void Print(std::ostream &s, int flag=0)
Imprime el objeto.
Definition: TwentyNodeBrick.cpp:3331
Communication parameters between processes.
Definition: CommParameters.h:65
================================================================================
Definition: ContinuaReprComponent.h:34
Definition: Response.h:71
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TwentyNodeBrick.cpp:2764