xc
BbarBrick.h
1 // -*-c++-*-
2 //----------------------------------------------------------------------------
3 // XC program; finite element analysis code
4 // for structural analysis and design.
5 //
6 // Copyright (C) Luis C. Pérez Tato
7 //
8 // This program derives from OpenSees <http://opensees.berkeley.edu>
9 // developed by the «Pacific earthquake engineering research center».
10 //
11 // Except for the restrictions that may arise from the copyright
12 // of the original program (see copyright_opensees.txt)
13 // XC is free software: you can redistribute it and/or modify
14 // it under the terms of the GNU General Public License as published by
15 // the Free Software Foundation, either version 3 of the License, or
16 // (at your option) any later version.
17 //
18 // This software is distributed in the hope that it will be useful, but
19 // WITHOUT ANY WARRANTY; without even the implied warranty of
20 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 // GNU General Public License for more details.
22 //
23 //
24 // You should have received a copy of the GNU General Public License
25 // along with this program.
26 // If not, see <http://www.gnu.org/licenses/>.
27 //----------------------------------------------------------------------------
28 /* ****************************************************************** **
29 ** OpenSees - Open System for Earthquake Engineering Simulation **
30 ** Pacific Earthquake Engineering Research Center **
31 ** **
32 ** **
33 ** (C) Copyright 1999, The Regents of the University of California **
34 ** All Rights Reserved. **
35 ** **
36 ** Commercial use of this program without express permission of the **
37 ** University of California, Berkeley, is strictly prohibited. See **
38 ** file 'COPYRIGHT' in main directory for information on usage and **
39 ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. **
40 ** **
41 ** Developed by: **
42 ** Frank McKenna (fmckenna@ce.berkeley.edu) **
43 ** Gregory L. Fenves (fenves@ce.berkeley.edu) **
44 ** Filip C. Filippou (filippou@ce.berkeley.edu) **
45 ** **
46 ** ****************************************************************** */
47 
48 // $Revision: 1.9 $
49 // $Date: 2003/08/28 22:42:32 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/element/brick/BbarBrick.h,v $
51 
52 // Ed "C++" Love
53 //
54 // Eight node BbarBrick element
55 //
56 
57 #include <domain/mesh/element/volumetric/BrickBase.h>
58 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
59 
60 namespace XC {
68 class BbarBrick: public BrickBase
69  {
70  private :
71  //static data
72  static Matrix stiff ;
73  static Vector resid ;
74  static Matrix mass ;
75  static Matrix damping ;
76 
77  //quadrature data
78  static const double sg[2] ;
79  static const double wg[8] ;
80 
81 
82  BodyForces3D bf;
83  mutable Matrix *Ki;
84 
85  void formInertiaTerms( int tangFlag ) const;
86  void formResidAndTangent( int tang_flag ) 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(Communicator &);
91  int recvData(const Communicator &);
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 ) const;
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  void alive(void);
125  int addLoad(ElementalLoad *theLoad, double loadFactor);
126  int addInertiaLoadToUnbalance(const Vector &accel);
127 
128  //get residual and residual with inertia terms
129  const Vector &getResistingForce(void) const;
130  const Vector &getResistingForceIncInertia(void) const;
131 
132  virtual int sendSelf(Communicator &);
133  virtual int recvSelf(const Communicator &);
134 
135  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
136  int getResponse(int responseID, Information &eleInformation);
137  };
138 
139 } // end of XC namespace
int sendData(Communicator &)
Send object members through the communicator argument.
Definition: BbarBrick.cpp:918
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: BbarBrick.cpp:132
void alive(void)
Reactivates the element.
Definition: BbarBrick.cpp:338
Float vector abstraction.
Definition: Vector.h:94
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: BbarBrick.cpp:969
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
Base class for hexahedral elements.
Definition: BrickBase.h:69
Base class response objects.
Definition: Response.h:81
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: BbarBrick.cpp:992
void Print(std::ostream &s, int flag) const
Print stuff.
Definition: BbarBrick.cpp:137
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: BbarBrick.cpp:927
Base class for the finite elements.
Definition: Element.h:112
Eight-node mixed volume/pressure brick element, which uses a trilinear isoparametric formulation...
Definition: BbarBrick.h:68
Body forces over an element.
Definition: BodyForces3D.h:40
const Vector & getResistingForceIncInertia(void) const
Get residual with inertia terms.
Definition: BbarBrick.cpp:404
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: BbarBrick.cpp:328
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Element * getCopy(void) const
Virtual constructor.
Definition: BbarBrick.cpp:109
Matrix of floats.
Definition: Matrix.h:111
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: BbarBrick.cpp:125
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: BbarBrick.cpp:158
Base class for loads over elements.
Definition: ElementalLoad.h:79
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: BbarBrick.cpp:391
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
virtual int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: BbarBrick.cpp:936
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: BbarBrick.cpp:950