xc
Brick.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.11 $
49 // $Date: 2006/01/10 18:41:34 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/element/brick/Brick.h,v $
51 
52 // Ed "C++" Love
53 //
54 // Eight node Brick element
55 //
56 
57 #include <domain/mesh/element/volumetric/BrickBase.h>
58 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
59 #include "domain/mesh/element/utils/fvectors/FVectorBrick.h"
60 
61 namespace XC {
65 class Brick: public BrickBase
66  {
67  private :
68  // static constants
69  static const int numberGauss= 8;
70  static const int nShape = 4;
71 
72  BodyForces3D bf;
73  FVectorBrick p0;
74  double appliedB[ndf];
75  bool applyLoad;
76 
77  mutable Matrix *Ki;
78  void shape_functions_loop(void) const;
79 
80  //
81  // static attributes
82  //
83 
84  static Matrix stiff;
85  static Vector resid;
86  static Matrix mass;
87  static Matrix damping;
88  static double gaussPoint[numberGauss][ndm];
89  static double dvol[numberGauss];
90  static double shp[nShape][numberNodes];
91  static double Shape[nShape][numberNodes][numberGauss];
92 
93  // quadrature data
94  static const double wg[numberGauss];
95 
96  //
97  // private methods
98  //
99 
100  //inertia terms
101  void formInertiaTerms(int tangFlag) const;
102  //form residual and tangent
103  void formResidAndTangent(int tang_flag) const;
104 
105  //compute B matrix
106  const Matrix& computeB( int node, const double shp[4][8]) const;
107 
108  static size_t getVectorIndex(const size_t &,const size_t &);
109  protected:
110  int sendData(Communicator &comm);
111  int recvData(const Communicator &comm);
112  public :
113 
114  Brick(void);
115  Brick(int tag,const NDMaterial *ptr_mat);
116 
117  //full constructor
118  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());
119  Element *getCopy(void) const;
120  //destructor
121  virtual ~Brick(void);
122 
123  //set domain
124  void setDomain( Domain *theDomain );
125 
126  //return number of dofs
127  int getNumDOF(void) const;
128 
129  // update
130  int update(void);
131 
132  Matrix getGaussPointsPositions(void) const;
133  //return stiffness matrix
134  const Matrix &getTangentStiff(void) const;
135  const Matrix &getInitialStiff(void) const;
136  const Matrix &getMass(void) const;
137 
138  Vector getAvgStress(void) const;
139  double getAvgStress(const size_t &,const size_t &) const;
140  Vector getAvgStrain(void) const;
141  double getAvgStrain(const size_t &,const size_t &) const;
142 
143  void alive(void);
144  void zeroLoad(void);
145  int addLoad(ElementalLoad *theLoad, double loadFactor);
146  int addInertiaLoadToUnbalance(const Vector &accel);
147 
148  Vector getShapeFunctionsValues(const double &r, const double &s, const double &t) const;
149 
150  const Vector &getResistingForce(void) const; //get residual
151  const Vector &getResistingForceIncInertia(void) const; //get residual with inertia terms
152  virtual void createInertiaLoad(const Vector &);
153 
154  // public methods for element output
155  virtual int sendSelf(Communicator &);
156  virtual int recvSelf(const Communicator &);
157 
158  //print out element data
159  void Print( std::ostream &s, int flag ) const;
160 
161  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
162  int getResponse(int responseID, Information &eleInformation);
163  };
164 
165 } // end of XC namespace
void setDomain(Domain *theDomain)
set domain
Definition: Brick.cpp:138
Float vector abstraction.
Definition: Vector.h:94
virtual ~Brick(void)
destructor
Definition: Brick.cpp:127
const Matrix & getMass(void) const
Return mass matrix.
Definition: Brick.cpp:348
int update(void)
Form residual and tangent.
Definition: Brick.cpp:619
Vector getAvgStrain(void) const
Return the average strain in the element.
Definition: Brick.cpp:179
Information about an element.
Definition: Information.h:81
Brick(void)
Default constructor.
Definition: Brick.cpp:98
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
Vector getAvgStress(void) const
Return the average stress in the element.
Definition: Brick.cpp:148
const Matrix & getTangentStiff(void) const
Return stiffness matrix.
Definition: Brick.cpp:192
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: Brick.cpp:947
int getNumDOF(void) const
return number of dofs
Definition: Brick.cpp:144
Eight node hexahedral element for three-dimensional problems.
Definition: Brick.h:65
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: Brick.cpp:933
Base class for the finite elements.
Definition: Element.h:112
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Brick.cpp:970
virtual int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: Brick.cpp:919
void alive(void)
Reactivates the element.
Definition: Brick.cpp:358
static const int ndf
Number of DOFs per node.
Definition: BrickBase.h:74
Ingernal forces for a brick element.
Definition: FVectorBrick.h:40
Body forces over an element.
Definition: BodyForces3D.h:40
const Matrix & getInitialStiff(void) const
Return initial stiffness matrix.
Definition: Brick.cpp:265
Vector getShapeFunctionsValues(const double &r, const double &s, const double &t) const
Return the values of the shape functions.
Definition: Brick.cpp:216
const Vector & getResistingForceIncInertia(void) const
Get residual with inertia terms.
Definition: Brick.cpp:485
void zeroLoad(void)
Remove element loads.
Definition: Brick.cpp:370
virtual void createInertiaLoad(const Vector &)
Creates the inertia load that corresponds to the acceleration argument.
Definition: Brick.cpp:581
int sendData(Communicator &comm)
Send members through the communicator argument.
Definition: Brick.cpp:899
static const int numberNodes
Number of nodes.
Definition: BrickBase.h:72
Matrix getGaussPointsPositions(void) const
Return the coordinates of the Gauss points.
Definition: Brick.cpp:203
void Print(std::ostream &s, int flag) const
Print out element data.
Definition: Brick.cpp:999
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
static const int ndm
Space dimension.
Definition: BrickBase.h:73
int recvData(const Communicator &comm)
Receives members through the communicator argument.
Definition: Brick.cpp:909
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
Element * getCopy(void) const
Virtual constructor.
Definition: Brick.cpp:122
Base class for loads over elements.
Definition: ElementalLoad.h:79
int addInertiaLoadToUnbalance(const Vector &accel)
Add inertia load due to acceleration argument.
Definition: Brick.cpp:438
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
const Vector & getResistingForce(void) const
Get residual.
Definition: Brick.cpp:469