XC Open source finite element analysis program
BbarBrick.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.9 $
48 // $Date: 2003/08/28 22:42:32 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/brick/BbarBrick.h,v $
50 
51 // Ed "C++" Love
52 //
53 // Eight node BbarBrick 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 BbarBrick : public BrickBase
64  {
65  private :
66  //static data
67  static Matrix stiff ;
68  static Vector resid ;
69  static Matrix mass ;
70  static Matrix damping ;
71 
72  //quadrature data
73  static const double sg[2] ;
74  static const double wg[8] ;
75 
76 
77  //local nodal coordinates, three coordinates for each of four nodes
78  // static double xl[3][8] ;
79  static double xl[][8] ;
80 
81  BodyForces3D bf;
82  mutable Matrix *Ki;
83 
84  void formInertiaTerms( int tangFlag ) const;
85  void formResidAndTangent( int tang_flag ) const;
86  void computeBasis(void) const;
87  const Matrix& computeBbar(int node, const double shp[4][8], const double shpBar[4][8]) const;
88  Matrix transpose( int dim1, int dim2, const Matrix &M ) ;
89  protected:
90  int sendData(CommParameters &);
91  int recvData(const CommParameters &);
92  public :
93  //null constructor
94  BbarBrick( ) ;
95  //full constructor
96  BbarBrick( int tag, int node1,
97  int node2,
98  int node3,
99  int node4,
100  int node5,
101  int node6,
102  int node7,
103  int node8,
104  NDMaterial &theMaterial,
105  const BodyForces3D &bForces= BodyForces3D()) ;
106  Element *getCopy(void) const;
107  //destructor
108  virtual ~BbarBrick( ) ;
109 
110  //set domain
111  void setDomain( Domain *theDomain ) ;
112 
113  //return number of dofs
114  int getNumDOF(void) const;
115 
116  //print out element data
117  void Print( std::ostream &s, int flag ) ;
118 
119  //return stiffness matrix
120  const Matrix &getTangentStiff(void) const;
121  const Matrix &getInitialStiff(void) const;
122  const Matrix &getMass(void) const;
123 
124  int addLoad(ElementalLoad *theLoad, double loadFactor);
125  int addInertiaLoadToUnbalance(const Vector &accel);
126 
127  //get residual and residual with inertia terms
128  const Vector &getResistingForce(void) const;
129  const Vector &getResistingForceIncInertia(void) const;
130 
131  virtual int sendSelf(CommParameters &);
132  virtual int recvSelf(const CommParameters &);
133 
134  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
135  int getResponse(int responseID, Information &eleInformation);
136  };
137 
138 } // end of XC namespace
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
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: BbarBrick.cpp:332
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: BbarBrick.cpp:927
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
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: BbarBrick.cpp:397
Definition: Matrix.h:82
Body forces over an element.
Definition: BodyForces3D.h:39
void Print(std::ostream &s, int flag)
Imprime el objeto.
Definition: BbarBrick.cpp:141
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: BbarBrick.cpp:129
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: BbarBrick.cpp:936
Hexaedro.
Definition: BbarBrick.h:63
virtual int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: BbarBrick.cpp:945
Communication parameters between processes.
Definition: CommParameters.h:65
Element * getCopy(void) const
Virtual constructor.
Definition: BbarBrick.cpp:113
virtual int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: BbarBrick.cpp:959
Base class for hexahedra.
Definition: BrickBase.h:44
================================================================================
Definition: ContinuaReprComponent.h:34
Definition: Response.h:71