xc
FourNodeQuadUP.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 //----------------------------------------------------------------------------
29 // Description: This file contains the class declaration for FourNodeQuadUP. //
30 // FourNodeQuadUP is a 4-node plane strain element for solid-fluid fully //
31 // coupled analysis. This implementation is a simplified u-p formulation //
32 // of Biot theory (u - solid displacement, p - fluid pressure). Each element //
33 // node has two DOFs for u and 1 DOF for p. //
34 // Written by Zhaohui Yang (May 2002) //
35 // based on FourNodeQuad element by Michael Scott //
37 
38 // $Revision: 1.1 $
39 // $Date: 2005/09/22 21:28:36 $
40 // $Source: /usr/local/cvs/OpenSees/SRC/element/UP-ucsd/FourNodeQuadUP.h,v $
41 
42 #ifndef FourNodeQuadUP_h
43 #define FourNodeQuadUP_h
44 
45 #include "domain/mesh/element/plane/QuadBase4N.h"
46 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h"
47 #include "domain/mesh/element/utils/body_forces/BodyForces2D.h"
48 
49 namespace XC {
50 class Node;
51 class NDMaterial;
52 class Response;
53 
55 //
61 class FourNodeQuadUP: public QuadBase4N<SolidMech2D>
62  {
63  private:
64 
65  static Matrix K;
66  static Vector P;
67  BodyForces2D bf;
68  Vector pressureLoad;
69 
70  double rho;
71  double kc;
72  double pressure;
73  double perm[2];
74 
75  static double shp[3][4][4];
76  static double pts[4][2];
77  static double wts[4];
78  static double dvol[4];
79  static double shpBar[3][4];
80 
81  Node *nd1Ptr(void);
82  const Node *nd1Ptr(void) const;
83  Node *nd2Ptr(void);
84  const Node *nd2Ptr(void) const;
85  Node *nd3Ptr(void);
86  const Node *nd3Ptr(void) const;
87  Node *nd4Ptr(void);
88  const Node *nd4Ptr(void) const;
89 
90  // private member functions - only objects of this class can call these
91  double mixtureRho(int ipt) const; // Mixture mass density at integration point i
92  void shapeFunction(void) const;
93  void setPressureLoadAtNodes(void);
94 
95  mutable Matrix *Ki;
96  protected:
97  int sendData(Communicator &comm);
98  int recvData(const Communicator &comm);
99  public:
100  FourNodeQuadUP(int tag, int nd1, int nd2, int nd3, int nd4, NDMaterial &m,const std::string &type, double t, double bulk, double perm1, double perm2,const BodyForces2D &bForces= BodyForces2D(), double p = 0.0);
101  Element *getCopy(void) const;
102  FourNodeQuadUP();
103  virtual ~FourNodeQuadUP();
104 
105  int getNumDOF(void) const;
106  void setDomain(Domain *theDomain);
107 
108  // public methods to set the state of the element
109  int update(void);
110 
111  // public methods to obtain stiffness, mass, damping and residual information
112  const Matrix &getTangentStiff(void) const;
113  const Matrix &getInitialStiff(void) const;
114  const Matrix &getDamp(void) const;
115  const Matrix &getMass(void) const;
116 
117  void alive(void);
118  int addLoad(ElementalLoad *theLoad, double loadFactor);
119  int addInertiaLoadToUnbalance(const Vector &accel);
120  const Vector &getResistingForce(void) const;
121  const Vector &getResistingForceIncInertia(void) const;
122 
123  // public methods for element output
124  int sendSelf(Communicator &);
125  int recvSelf(const Communicator &);
126  void Print(std::ostream &s, int flag =0) const;
127 
128  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
129  int getResponse(int responseID, Information &eleInformation);
130 
131  int setParameter(const std::vector<std::string> &argv, Parameter &param);
132  int updateParameter(int parameterID, Information &info);
133 
134 
135  // RWB; PyLiq1 & TzLiq1 need to see the excess pore pressure and initial stresses.
136  friend class PyLiq1;
137  friend class TzLiq1;
138  };
139 } // end of XC namespace
140 
141 #endif
142 
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: FourNodeQuadUP.cpp:684
Body forces over an element.
Definition: BodyForces2D.h:40
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: FourNodeQuadUP.cpp:706
Float vector abstraction.
Definition: Vector.h:94
int sendSelf(Communicator &)
Send the object.
Definition: FourNodeQuadUP.cpp:638
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
Base class response objects.
Definition: Response.h:81
void alive(void)
Reactivates the element.
Definition: FourNodeQuadUP.cpp:414
int recvSelf(const Communicator &)
Receive the object.
Definition: FourNodeQuadUP.cpp:651
int update(void)
Updates the element state.
Definition: FourNodeQuadUP.cpp:147
Uniaxial p-y material that incorporates liquefaction effects.
Definition: PyLiq1.h:61
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: FourNodeQuadUP.cpp:669
Base class for 4 node quads.
Definition: QuadBase4N.h:45
Base class for the finite elements.
Definition: Element.h:112
Element * getCopy(void) const
Virtual constructor.
Definition: FourNodeQuadUP.cpp:109
int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: FourNodeQuadUP.cpp:751
Uniaxial t-z material that incorporates liquefaction effects.
Definition: TzLiq1.h:61
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: FourNodeQuadUP.cpp:134
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
Four-node plane-strain element using bilinear isoparametric formulation.
Definition: FourNodeQuadUP.h:61
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
int sendData(Communicator &comm)
Send object members through the communicator argument.
Definition: FourNodeQuadUP.cpp:621
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: FourNodeQuadUP.cpp:541
Parameter.
Definition: Parameter.h:68
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: FourNodeQuadUP.cpp:193
Base class for loads over elements.
Definition: ElementalLoad.h:79
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: FourNodeQuadUP.cpp:480
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
Mesh node.
Definition: Node.h:111
const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: FourNodeQuadUP.cpp:290
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: FourNodeQuadUP.cpp:354
int setParameter(const std::vector< std::string > &argv, Parameter &param)
Sets the value param to the parameter argv.
Definition: FourNodeQuadUP.cpp:722
int recvData(const Communicator &comm)
Receives object members through the communicator argument.
Definition: FourNodeQuadUP.cpp:630
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: FourNodeQuadUP.cpp:137