xc
TwentyEightNodeBrickUP.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 // 20-8 Noded brick element: TwentyEightNodeBrickUP
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 namespace XC {
72 class TwentyEightNodeBrickUP: public ElemWithMaterial<20,NDMaterialPhysicalProperties>
73  {
74  private:
75  BodyForces3D bf;
76  double rho_f;
77  double kc;
78  double perm[3];
79  mutable Matrix *Ki;
80 
81  //static data
82  static Matrix stiff ;
83  static Vector resid ;
84  static Matrix mass ;
85  static Matrix damp ;
86 
87  //quadrature data
88  static const int nintu;
89  static const int nintp;
90  static const int nenu;
91  static const int nenp;
92 
93 
94  //local nodal coordinates, three coordinates for each of twenty nodes
95  // static double xl[3][20] ;
96  static double xl[3][20] ;
97 
98  static double shgu[4][20][27]; // Stores shape functions and derivatives (overwritten)
99  static double shgp[4][8][8]; // Stores shape functions and derivatives (overwritten)
100  static double shgq[4][20][8]; // Stores shape functions and derivatives (overwritten)
101  static double shlu[4][20][27]; // Stores shape functions and derivatives
102  static double shlp[4][8][8]; // Stores shape functions and derivatives
103  static double shlq[4][20][8]; // Stores shape functions and derivatives
104  static double wu[27]; // Stores quadrature weights
105  static double wp[8]; // Stores quadrature weights
106  static double dvolu[27]; // Stores detJacobian (overwritten)
107  static double dvolp[8]; // Stores detJacobian (overwritten)
108  static double dvolq[8]; // Stores detJacobian (overwritten)
109 
110  //inertia terms
111  void formInertiaTerms( int tangFlag ) const;
112  //damping terms
113  void formDampingTerms( int tangFlag ) const;
114 
115  // Mixture mass density at integration point i
116  double mixtureRho(int ipt) const;
117 
118  //compute coordinate system
119  void computeBasis(void) const;
120 
121  // compute local shape functions
122  void compuLocalShapeFunction();
123  void Jacobian3d(int gaussPoint, double& xsj, int mode) const;
124  const Matrix &getStiff( int flag) const;
125  protected:
126  int sendData(Communicator &);
127  int recvData(const Communicator &);
128  public :
129  //null constructor
131 
132  //full constructor
133  TwentyEightNodeBrickUP( int tag,int node1,
134  int node2,
135  int node3,
136  int node4,
137  int node5,
138  int node6,
139  int node7,
140  int node8,
141  int node9,
142  int node10,
143  int node11,
144  int node12,
145  int node13,
146  int node14,
147  int node15,
148  int node16,
149  int node17,
150  int node18,
151  int node19,
152  int node20,
153  NDMaterial &theMaterial,double bulk, double rhof,
154  double perm1, double perm2, double perm3,
155  const BodyForces3D &bForces= BodyForces3D());
156 
157  //destructor
158  virtual ~TwentyEightNodeBrickUP(void) ;
159  Element *getCopy(void) const;
160  //set domain
161  void setDomain( Domain *theDomain ) ;
162 
163  //return number of dofs
164  int getNumDOF(void) const;
165 
166  //print out element data
167  void Print( std::ostream &s, int flag ) const;
168 
169  int update(void);
170 
171  //return stiffness matrix
172  const Matrix &getTangentStiff(void) const;
173  const Matrix &getInitialStiff(void) const;
174  const Matrix &getDamp(void) const;
175  const Matrix &getMass(void) const;
176 
177  void alive(void);
178  int addLoad(ElementalLoad *theLoad, double loadFactor);
179  int addInertiaLoadToUnbalance(const Vector &accel);
180 
181  //get residual
182  const Vector &getResistingForce(void) const;
183  const Vector &getResistingForceIncInertia(void) const;
184 
185  virtual int sendSelf(Communicator &);
186  virtual int recvSelf(const Communicator &);
187 
188  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
189  int getResponse(int responseID, Information &eleInformation);
190 
191  // RWB; PyLiq1 & TzLiq1 need to see the excess pore pressure and initial stresses.
192  friend class PyLiq1;
193  friend class TzLiq1;
194  };
195 
196 } // end of XC namespace
Float vector abstraction.
Definition: Vector.h:94
Information about an element.
Definition: Information.h:81
Twenty eight node hexahedral element for three-dimensional coupled problems.
Definition: TwentyEightNodeBrickUP.h:72
Communication parameters between processes.
Definition: Communicator.h:66
Base class response objects.
Definition: Response.h:81
Uniaxial p-y material that incorporates liquefaction effects.
Definition: PyLiq1.h:61
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TwentyEightNodeBrickUP.cpp:488
void Print(std::ostream &s, int flag) const
Print stuff.
Definition: TwentyEightNodeBrickUP.cpp:182
Base class for the finite elements.
Definition: Element.h:112
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: TwentyEightNodeBrickUP.cpp:344
Body forces over an element.
Definition: BodyForces3D.h:40
Uniaxial t-z material that incorporates liquefaction effects.
Definition: TzLiq1.h:61
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: TwentyEightNodeBrickUP.cpp:975
int sendData(Communicator &)
Send object members through the communicator argument.
Definition: TwentyEightNodeBrickUP.cpp:941
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: TwentyEightNodeBrickUP.cpp:1020
int getNumDOF(void) const
Return the number of element DOFs.
Definition: TwentyEightNodeBrickUP.cpp:178
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: TwentyEightNodeBrickUP.cpp:951
virtual int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: TwentyEightNodeBrickUP.cpp:961
Element * getCopy(void) const
Virtual constructor.
Definition: TwentyEightNodeBrickUP.cpp:145
int update(void)
Updates the element state.
Definition: TwentyEightNodeBrickUP.cpp:255
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: TwentyEightNodeBrickUP.cpp:770
void alive(void)
Reactivates the element.
Definition: TwentyEightNodeBrickUP.cpp:605
const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: TwentyEightNodeBrickUP.cpp:499
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: TwentyEightNodeBrickUP.cpp:160
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: TwentyEightNodeBrickUP.cpp:666
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: TwentyEightNodeBrickUP.cpp:994
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
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