xc
Twenty_Node_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 // by Jinchi Lu and Zhaohui Yang (May 2004)
49 //
50 // 20NodeBrick element
51 //
52 
53 #include <domain/mesh/element/ElemWithMaterial.h>
54 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h"
55 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
56 
57 
58 namespace XC {
59  class NDMaterial;
63 class Twenty_Node_Brick: public ElemWithMaterial<20,NDMaterialPhysicalProperties>
64  {
65  private:
66  BodyForces3D bf;
67  mutable Matrix *Ki;
68 
70  //static data
71  static Matrix stiff ;
72  static Vector resid ;
73  static Matrix mass ;
74  static Matrix damp ;
75 
76  //quadrature data
77  static const int nintu;
78  static const int nenu;
79 
80  //local nodal coordinates, three coordinates for each of twenty nodes
81  // static double xl[3][20] ;
82  static double xl[3][20] ;
83 
84  static double shgu[4][20][27]; // Stores shape functions and derivatives (overwritten)
85  static double shlu[4][20][27]; // Stores shape functions and derivatives
86  static double wu[27]; // Stores quadrature weights
87  static double dvolu[27]; // Stores detJacobian (overwritten)
88 
89  void formInertiaTerms( int tangFlag ) const;
90  void formDampingTerms( int tangFlag ) const;
91  double mixtureRho(int ipt) const;
92  void computeBasis(void) const;
93 
94 
95  // compute local shape functions
96  void compuLocalShapeFunction(void);
97  static void Jacobian3d(int gaussPoint, double& xsj, int mode);
98  const Matrix &getStiff( int flag ) const;
99 
100  protected:
101  int sendData(Communicator &);
102  int recvData(const Communicator &);
103  public :
104  //null constructor
105  Twenty_Node_Brick( ) ;
106 
107  //full constructor
108  Twenty_Node_Brick( int tag,
109  int node1,
110  int node2,
111  int node3,
112  int node4,
113  int node5,
114  int node6,
115  int node7,
116  int node8,
117  int node9,
118  int node10,
119  int node11,
120  int node12,
121  int node13,
122  int node14,
123  int node15,
124  int node16,
125  int node17,
126  int node18,
127  int node19,
128  int node20,
129  NDMaterial &theMaterial,
130  const BodyForces3D &bForces= BodyForces3D()) ;
131  Element *getCopy(void) const;
132  //destructor
133  virtual ~Twenty_Node_Brick(void) ;
134 
135  //set domain
136  void setDomain( Domain *theDomain ) ;
137 
138  //return number of dofs
139  int getNumDOF(void) const;
140 
141  //print out element data
142  void Print(std::ostream &s, int flag) const;
143 
144  int update(void);
145 
146  //return stiffness matrix
147  const Matrix &getTangentStiff(void) const;
148  const Matrix &getInitialStiff(void) const;
149  const Matrix &getDamp(void) const;
150  const Matrix &getMass(void) const;
151 
152  void alive(void);
153  int addLoad(ElementalLoad *theLoad, double loadFactor);
154  int addInertiaLoadToUnbalance(const Vector &accel);
155 
156  //get residual
157  const Vector &getResistingForce(void) const;
158  const Vector &getResistingForceIncInertia(void) const;
159 
160  // public methods for element output
161  int sendSelf(Communicator &);
162  int recvSelf(const Communicator &);
163 
164  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
165  int getResponse(int responseID, Information &eleInformation);
166 
167  //plotting
168  int displaySelf(Renderer &theViewer, int displayMode, float fact);
169  };
170 
171 } // end of XC namespace
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: Twenty_Node_Brick.cpp:767
Float vector abstraction.
Definition: Vector.h:94
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: Twenty_Node_Brick.cpp:305
Base class response objects.
Definition: Response.h:81
int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: Twenty_Node_Brick.cpp:776
int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: Twenty_Node_Brick.cpp:790
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Twenty_Node_Brick.cpp:436
void Print(std::ostream &s, int flag) const
Print stuff.
Definition: Twenty_Node_Brick.cpp:144
Twenty node exahedron.
Definition: Twenty_Node_Brick.h:63
const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: Twenty_Node_Brick.cpp:447
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: Twenty_Node_Brick.cpp:139
int sendData(Communicator &)
Send object members through the communicator argument.
Definition: Twenty_Node_Brick.cpp:758
Base class for the finite elements.
Definition: Element.h:112
Element * getCopy(void) const
Virtual constructor.
Definition: Twenty_Node_Brick.cpp:115
Body forces over an element.
Definition: BodyForces3D.h:40
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: Twenty_Node_Brick.cpp:809
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: Twenty_Node_Brick.cpp:541
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: Twenty_Node_Brick.cpp:623
void alive(void)
Reactivates the element.
Definition: Twenty_Node_Brick.cpp:484
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
int update(void)
Updates the element state.
Definition: Twenty_Node_Brick.cpp:216
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Twenty_Node_Brick.cpp:835
Element with material.
Definition: ElemWithMaterial.h:45
Matrix of floats.
Definition: Matrix.h:111
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: Twenty_Node_Brick.cpp:131