XC Open source finite element analysis program
Brick.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 //----------------------------------------------------------------------------
27 /* ****************************************************************** **
28 ** OpenSees - Open System for Earthquake Engineering Simulation **
29 ** Pacific Earthquake Engineering Research Center **
30 ** **
31 ** **
32 ** (C) Copyright 1999, The Regents of the University of California **
33 ** All Rights Reserved. **
34 ** **
35 ** Commercial use of this program without express permission of the **
36 ** University of California, Berkeley, is strictly prohibited. See **
37 ** file 'COPYRIGHT' in main directory for information on usage and **
38 ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. **
39 ** **
40 ** Developed by: **
41 ** Frank McKenna (fmckenna@ce.berkeley.edu) **
42 ** Gregory L. Fenves (fenves@ce.berkeley.edu) **
43 ** Filip C. Filippou (filippou@ce.berkeley.edu) **
44 ** **
45 ** ****************************************************************** */
46 
47 // $Revision: 1.11 $
48 // $Date: 2006/01/10 18:41:34 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/brick/Brick.h,v $
50 
51 // Ed "C++" Love
52 //
53 // Eight node Brick element
54 //
55 
56 #include <domain/mesh/element/volumen/BrickBase.h>
57 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
58 
59 namespace XC {
61 //
63 class Brick : public BrickBase
64  {
65  private :
66  BodyForces3D bf;
67  mutable Matrix *Ki;
68 
69  //
70  // static attributes
71  //
72 
73  static Matrix stiff;
74  static Vector resid;
75  static Matrix mass;
76  static Matrix damping;
77 
78  //quadrature data
79  static const double sg[2];
80  static const double wg[8];
81 
82  //local nodal coordinates, three coordinates for each of eight nodes
83  static double xl[3][8];
84 
85  //
86  // private methods
87  //
88 
89  //inertia terms
90  void formInertiaTerms(int tangFlag) const;
91  //form residual and tangent
92  void formResidAndTangent(int tang_flag) const;
93 
94  //compute coordinate system
95  void computeBasis(void) const;
96 
97  //compute B matrix
98  const Matrix& computeB( int node, const double shp[4][8]) const;
99 
100  //Matrix transpose
101  Matrix transpose( int dim1, int dim2, const Matrix &M );
102  static size_t getVectorIndex(const size_t &,const size_t &);
103  protected:
104  int sendData(CommParameters &cp);
105  int recvData(const CommParameters &cp);
106  public :
107 
108  //null constructor
109  Brick(void);
110  Brick(int tag,const NDMaterial *ptr_mat);
111 
112  //full constructor
113  Brick( int tag, int node1,int node2,int node3,int node4,int node5,int node6,int node7,int node8, NDMaterial &theMaterial, const BodyForces3D &bf= BodyForces3D());
114  Element *getCopy(void) const;
115  //destructor
116  virtual ~Brick(void);
117 
118  //set domain
119  void setDomain( Domain *theDomain );
120 
121  //return number of dofs
122  int getNumDOF(void) const;
123 
124  // update
125  int update(void);
126 
127  //return stiffness matrix
128  const Matrix &getTangentStiff() const;
129  const Matrix &getInitialStiff() const;
130  const Matrix &getMass() const;
131 
132  Vector getAvgStress(void) const;
133  double getAvgStress(const size_t &,const size_t &) const;
134  Vector getAvgStrain(void) const;
135  double getAvgStrain(const size_t &,const size_t &) const;
136 
137  int addLoad(ElementalLoad *theLoad, double loadFactor);
138  int addInertiaLoadToUnbalance(const Vector &accel);
139 
140  //get residual
141  const Vector &getResistingForce(void) const;
142 
143  //get residual with inertia terms
144  const Vector &getResistingForceIncInertia(void) const;
145 
146  // public methods for element output
147  virtual int sendSelf(CommParameters &);
148  virtual int recvSelf(const CommParameters &);
149 
150  //print out element data
151  void Print( std::ostream &s, int flag );
152 
153  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
154  int getResponse(int responseID, Information &eleInformation);
155  };
156 
157 } // end of XC namespace
virtual int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: Brick.cpp:1121
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: Brick.cpp:471
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
Hexaedro de ocho nodos.
Definition: Brick.h:63
Base class for 2D and 3D materials.
Definition: NDMaterial.h:91
Definition: Vector.h:82
int sendData(CommParameters &cp)
Send members through the channel being passed as parameter.
Definition: Brick.cpp:1090
const Matrix & getMass() const
Returns the mass matrix.
Definition: Brick.cpp:404
Brick(void)
null constructor
Definition: Brick.cpp:92
void Print(std::ostream &s, int flag)
Imprime el objeto.
Definition: Brick.cpp:185
Information about an element.
Definition: Information.h:80
Element * getCopy(void) const
Virtual constructor.
Definition: Brick.cpp:117
Vector getAvgStrain(void) const
Return the average strain in the element.
Definition: Brick.cpp:174
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
void setDomain(Domain *theDomain)
set domain
Definition: Brick.cpp:133
virtual ~Brick(void)
destructor
Definition: Brick.cpp:122
int getNumDOF(void) const
return number of dofs
Definition: Brick.cpp:140
Definition: Matrix.h:82
int update(void)
Actualiza el estado of the element.
Definition: Brick.cpp:637
Vector getAvgStress(void) const
Return the tensión media in the element.
Definition: Brick.cpp:144
Body forces over an element.
Definition: BodyForces3D.h:39
int recvData(const CommParameters &cp)
Receives members through the channel being passed as parameter.
Definition: Brick.cpp:1099
virtual int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: Brick.cpp:1108
Communication parameters between processes.
Definition: CommParameters.h:65
Base class for hexahedra.
Definition: BrickBase.h:44
================================================================================
Definition: ContinuaReprComponent.h:34
Definition: Response.h:71