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