XC Open source finite element analysis program
TwentyEightNodeBrickUP.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 // 20-8 Noded brick element: TwentyEightNodeBrickUP
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 namespace XC {
58 //
60 class TwentyEightNodeBrickUP : public ElemWithMaterial<20,NDMaterialPhysicalProperties>
61  {
62  private:
63  BodyForces3D bf;
64  double rho; // Fluid mass per unit volume
65  double kc; // combined bulk modulus
66  double perm[3]; // permeability
67  mutable Matrix *Ki;
68 
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 nintp;
78  static const int nenu;
79  static const int nenp;
80 
81 
82  //local nodal coordinates, three coordinates for each of twenty nodes
83  // static double xl[3][20] ;
84  static double xl[3][20] ;
85 
86  static double shgu[4][20][27]; // Stores shape functions and derivatives (overwritten)
87  static double shgp[4][8][8]; // Stores shape functions and derivatives (overwritten)
88  static double shgq[4][20][8]; // Stores shape functions and derivatives (overwritten)
89  static double shlu[4][20][27]; // Stores shape functions and derivatives
90  static double shlp[4][8][8]; // Stores shape functions and derivatives
91  static double shlq[4][20][8]; // Stores shape functions and derivatives
92  static double wu[27]; // Stores quadrature weights
93  static double wp[8]; // Stores quadrature weights
94  static double dvolu[27]; // Stores detJacobian (overwritten)
95  static double dvolp[8]; // Stores detJacobian (overwritten)
96  static double dvolq[8]; // Stores detJacobian (overwritten)
97 
98  //inertia terms
99  void formInertiaTerms( int tangFlag ) const;
100  //damping terms
101  void formDampingTerms( int tangFlag ) const;
102 
103  // Mixture mass density at integration point i
104  double mixtureRho(int ipt) const;
105 
106  //compute coordinate system
107  void computeBasis(void) const;
108 
109  // compute local shape functions
110  void compuLocalShapeFunction();
111  void Jacobian3d(int gaussPoint, double& xsj, int mode) const;
112  const Matrix &getStiff( int flag) const;
113  protected:
114  int sendData(CommParameters &);
115  int recvData(const CommParameters &);
116  public :
117  //null constructor
119 
120  //full constructor
121  TwentyEightNodeBrickUP( int tag,int node1,
122  int node2,
123  int node3,
124  int node4,
125  int node5,
126  int node6,
127  int node7,
128  int node8,
129  int node9,
130  int node10,
131  int node11,
132  int node12,
133  int node13,
134  int node14,
135  int node15,
136  int node16,
137  int node17,
138  int node18,
139  int node19,
140  int node20,
141  NDMaterial &theMaterial,double bulk, double rhof,
142  double perm1, double perm2, double perm3,
143  const BodyForces3D &bForces= BodyForces3D());
144 
145  //destructor
146  virtual ~TwentyEightNodeBrickUP(void) ;
147  Element *getCopy(void) const;
148  //set domain
149  void setDomain( Domain *theDomain ) ;
150 
151  //return number of dofs
152  int getNumDOF(void) const;
153 
154  //print out element data
155  void Print( std::ostream &s, int flag ) ;
156 
157  int update(void);
158 
159  //return stiffness matrix
160  const Matrix &getTangentStiff(void) const;
161  const Matrix &getInitialStiff(void) const;
162  const Matrix &getDamp(void) const;
163  const Matrix &getMass(void) const;
164 
165  int addLoad(ElementalLoad *theLoad, double loadFactor);
166  int addInertiaLoadToUnbalance(const Vector &accel);
167 
168  //get residual
169  const Vector &getResistingForce(void) const;
170  const Vector &getResistingForceIncInertia(void) const;
171 
172  virtual int sendSelf(CommParameters &);
173  virtual int recvSelf(const CommParameters &);
174 
175  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
176  int getResponse(int responseID, Information &eleInformation);
177 
178  // RWB; PyLiq1 & TzLiq1 need to see the excess pore pressure and initial stresses.
179  friend class PyLiq1;
180  friend class TzLiq1;
181  };
182 
183 } // end of XC namespace
??.
Definition: TzLiq1.h:60
Element with material.
Definition: ElemWithMaterial.h:40
Element * getCopy(void) const
Virtual constructor.
Definition: TwentyEightNodeBrickUP.cpp:146
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
Base class for 2D and 3D materials.
Definition: NDMaterial.h:91
Definition: Vector.h:82
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TwentyEightNodeBrickUP.cpp:489
void Print(std::ostream &s, int flag)
Imprime el objeto.
Definition: TwentyEightNodeBrickUP.cpp:183
Information about an element.
Definition: Information.h:80
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: TwentyEightNodeBrickUP.cpp:759
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: TwentyEightNodeBrickUP.cpp:930
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
Definition: Matrix.h:82
??.
Definition: PyLiq1.h:62
Body forces over an element.
Definition: BodyForces3D.h:39
virtual int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: TwentyEightNodeBrickUP.cpp:964
int update(void)
Actualiza el estado of the element.
Definition: TwentyEightNodeBrickUP.cpp:256
Communication parameters between processes.
Definition: CommParameters.h:65
================================================================================
Definition: ContinuaReprComponent.h:34
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: TwentyEightNodeBrickUP.cpp:940
const Matrix & getDamp(void) const
Returns the matriz de amortiguamiento.
Definition: TwentyEightNodeBrickUP.cpp:500
Definition: Response.h:71
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: TwentyEightNodeBrickUP.cpp:161
int getNumDOF(void) const
Return the number of element DOFs.
Definition: TwentyEightNodeBrickUP.cpp:179
virtual int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: TwentyEightNodeBrickUP.cpp:950
Hexaedro de veintiocho nodos.
Definition: TwentyEightNodeBrickUP.h:60