xc
BrickUP.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 
49 // Description: This file contains the class declaration for //
50 // BrickUP, an 8-node cubic element for solid-fluid fully coupled analysis. //
51 // This implementation is a simplified u-p formulation of Biot theory //
52 // (u - solid displacement, p - fluid pressure). Each element node has three //
53 // DOFs for u and 1 DOF for p. //
54 // //
55 // Written by Zhaohui Yang (March 2004) //
56 // //
58 
59 // $Revision: 1.1 $
60 // $Date: 2005/09/22 21:28:36 $
61 // $Source: /usr/local/cvs/OpenSees/SRC/element/UP-ucsd/BrickUP.h,v $
62 
63 // by Zhaohui Yang (Modified based on Ed "C++" Love's Brick element)
64 //
65 // Eight node BrickUP element
66 //
67 
68 #include <domain/mesh/element/volumetric/BrickBase.h>
69 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
70 
71 namespace XC{
72 class NDMaterial;
76 class BrickUP: public BrickBase
77  {
78  private :
79  BodyForces3D bf;
80  double rho_f;
81  double kc;
82  double perm[3];
83  mutable Matrix *Ki;
84 
85  //static data
86  static Matrix stiff ;
87  static Vector resid ;
88  static Matrix mass ;
89  static Matrix damp ;
90 
91  //quadrature data
92  static const double sg[2] ;
93  static const double wg[8] ;
94 
95 
96  void formInertiaTerms(int tangFlag) const;
97  void formDampingTerms(int tangFlag) const;
98  void formResidAndTangent(int tang_flag) const;
99  double mixtureRho(int ipt) const;
100  const Matrix& computeBbar(int node,const double shp[4][8],const double shpBar[4][8] ) ;
101  const Matrix& computeB( int node, const double shp[4][8] ) const;
102  Matrix transpose( int dim1, int dim2, const Matrix &M ) ;
103  protected:
104  int sendData(Communicator &);
105  int recvData(const Communicator &);
106  public :
107  //null constructor
108  BrickUP(void);
109 
110  //full constructor
111  BrickUP( int tag, int node1,int node2,int node3,int node4,int node5,int node6,int node7,int node8,NDMaterial &theMaterial,double bulk, double rhof, double perm1, double perm2, double perm3, const BodyForces3D &bForces= BodyForces3D()) ;
112  Element *getCopy(void) const;
113  //destructor
114  virtual ~BrickUP(void) ;
115 
116  //set domain
117  void setDomain( Domain *theDomain );
118 
119  //return number of dofs
120  int getNumDOF(void) const;
121 
122  //print out element data
123  void Print(std::ostream &s, int flag) const;
124 
125  //return stiffness matrix
126  const Matrix &getTangentStiff(void) const;
127  const Matrix &getInitialStiff(void) const;
128  const Matrix &getDamp(void) const;
129  const Matrix &getMass(void) const;
130 
131  void alive(void);
132  int addLoad(ElementalLoad *theLoad, double loadFactor);
133  int addInertiaLoadToUnbalance(const Vector &accel);
134 
135  //get residual
136  const Vector &getResistingForce(void) const;
137  const Vector &getResistingForceIncInertia(void) const;
138 
139  virtual int sendSelf(Communicator &);
140  virtual int recvSelf(const Communicator &);
141 
142  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
143  int getResponse(int responseID, Information &eleInformation);
144 
145 
146  // RWB; PyLiq1 & TzLiq1 need to see the excess pore pressure and initial stresses.
147  friend class PyLiq1;
148  friend class TzLiq1;
149 } ;
150 
151 } // end of XC namespace
Float vector abstraction.
Definition: Vector.h:94
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: BrickUP.cpp:1069
Information about an element.
Definition: Information.h:81
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
const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: BrickUP.cpp:378
Uniaxial p-y material that incorporates liquefaction effects.
Definition: PyLiq1.h:61
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: BrickUP.cpp:566
Eight node exahedron.
Definition: BrickUP.h:76
Base class for the finite elements.
Definition: Element.h:112
int sendData(Communicator &)
Send object members through the communicator argument.
Definition: BrickUP.cpp:991
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: BrickUP.cpp:1001
Body forces over an element.
Definition: BodyForces3D.h:40
Uniaxial t-z material that incorporates liquefaction effects.
Definition: TzLiq1.h:61
Element * getCopy(void) const
Virtual constructor.
Definition: BrickUP.cpp:125
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: BrickUP.cpp:211
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
void alive(void)
Reactivates the element.
Definition: BrickUP.cpp:514
virtual int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: BrickUP.cpp:1011
int addLoad(ElementalLoad *theLoad, double loadFactor)
Adds to the element the load being passed as parameter.
Definition: BrickUP.cpp:525
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
void Print(std::ostream &s, int flag) const
Print stuff.
Definition: BrickUP.cpp:150
Matrix of floats.
Definition: Matrix.h:111
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: BrickUP.cpp:1044
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
virtual int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: BrickUP.cpp:1025
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: BrickUP.cpp:137
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: BrickUP.cpp:367
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: BrickUP.cpp:145
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: BrickUP.cpp:580