xc
Nine_Four_Node_QuadUP.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 //
30 // NineFourNodeQuadUP, a 9-4-node (9 node for solid and 4 node for fluid) //
31 // plane strain element for solid-fluid fully coupled analysis. This //
32 // implementation is a simplified u-p formulation of Biot theory //
33 // (u - solid displacement, p - fluid pressure). Each element node has two //
34 // DOFs for u and 1 DOF for p. //
35 // //
36 // Written by Zhaohui Yang (March 2004) //
37 // //
39 
40 // $Revision: 1.5 $
41 // $Date: 2007-02-02 01:44:56 $
42 // $Source: /usr/local/cvs/OpenSees/SRC/element/UP-ucsd/Nine_Four_Node_QuadUP.h,v $
43 
44 #ifndef NineFourNodeQuadUP_h
45 #define NineFourNodeQuadUP_h
46 
47 #include <domain/mesh/element/ElemWithMaterial.h>
48 #include <utility/matrix/Matrix.h>
49 #include <utility/matrix/Vector.h>
50 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h"
51 #include "domain/mesh/element/utils/body_forces/BodyForces2D.h"
52 
53 namespace XC {
54 class Node;
55 class NDMaterial;
56 class Response;
57 
58 class NineFourNodeQuadUP: public ElemWithMaterial<9,SolidMech2D>
59  {
60  private:
61  BodyForces2D bf;
62  mutable Matrix *Ki;
63 
64  double kc;
65  double fluidRho;
66  double perm[2];
67  double appliedB[2]; // Body forces applied with load pattern, C.McGann, U.Washington
68 
69  int applyLoad; // flag for body force in load, C.McGann, U.Washington
70 
71  static Matrix K;
72  static Vector P; //"< Element resisting force vector
73  static const int nintu;
74  static const int nintp;
75  static const int nenu;
76  static const int nenp;
77 
78 
79  static double shgu[3][9][9];
80  static double shgp[3][4][4];
81  static double shgq[3][9][4];
82  static double shlu[3][9][9];
83  static double shlp[3][4][4];
84  static double shlq[3][9][4];
85  static double wu[9];
86  static double wp[4];
87  static double dvolu[9];
88  static double dvolp[4];
89  static double dvolq[4]; // Stores detJacobian (overwritten)
90 
91  // private member functions - only objects of this class can call these
92  double mixtureRho(int ipt) const; // Mixture mass density at integration point i
93  void shapeFunction(double *w, int nint, int nen, int mode) const;
94  void globalShapeFunction(double *dvol, double *w, int nint, int nen, int mode) const;
95 
96  double *initNodeDispl;
97  protected:
98  int sendData(Communicator &);
99  int recvData(const Communicator &);
100  public:
101  NineFourNodeQuadUP(int tag= 0,const NDMaterial *ptr_mat= nullptr);
102  NineFourNodeQuadUP(int tag, int nd1, int nd2, int nd3, int nd4,
103  int nd5, int nd6, int nd7, int nd8, int nd9,
104  NDMaterial &m, const char *type,
105  double t, double bulk, double frho, double perm1, double perm2,
106  const BodyForces2D &bForces= BodyForces2D());
107 
108  Element *getCopy(void) const;
109  virtual ~NineFourNodeQuadUP(void);
110 
111  inline double getCombinedBulkModulus(void) const
112  { return kc; }
113  void setCombinedBulkModulus(const double &d)
114  { kc= d; }
115  inline double getHorizontalPermeability(void) const
116  { return perm[0]; }
117  void setHorizontalPermeability(const double &d)
118  { perm[0]= d; }
119  inline double getVerticalPermeability(void) const
120  { return perm[1]; }
121  void setVerticalPermeability(const double &d)
122  { perm[1]= d; }
123  inline double getFluidRho(void) const
124  { return fluidRho; }
125  void setFluidRho(const double &d)
126  { fluidRho= d; }
127  inline const BodyForces2D &getBodyForces(void) const
128  { return bf; }
129  void setBodyForces(const BodyForces2D &f)
130  { bf= f; }
131 
132  int getNumDOF(void) const;
133  void setDomain(Domain *theDomain);
134 
135  // public methods to set the state of the element
136  int update(void);
137 
138  // public methods to obtain stiffness, mass, damping and residual information
139  const Matrix &getTangentStiff(void) const;
140  const Matrix &getInitialStiff(void) const;
141  const Matrix &getDamp(void) const;
142  const Matrix &getMass(void) const;
143 
144  void alive(void);
145  void zeroLoad(void);
146  int addLoad(ElementalLoad *theLoad, double loadFactor);
147  int addInertiaLoadToUnbalance(const Vector &);
148  const Vector &getResistingForce(void) const;
149  const Vector &getResistingForceIncInertia(void) const;
150 
151  // public methods for element output
152  int sendSelf(Communicator &);
153  int recvSelf(const Communicator &);
154  void Print(std::ostream &s, int flag =0) const;
155 
156  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
157  int getResponse(int responseID, Information &eleInformation);
158 
159  int setParameter(const std::vector<std::string> &argv, Parameter &param);
160  int updateParameter(int parameterID, Information &info);
161 
162  // RWB; PyLiq1 & TzLiq1 need to see the excess pore pressure and
163  // initial stresses.
164  friend class PyLiq1;
165  friend class TzLiq1;
166  };
167 } // end of XC namespace
168 
169 #endif
Body forces over an element.
Definition: BodyForces2D.h:40
Float vector abstraction.
Definition: Vector.h:94
int update(void)
Updates the element state.
Definition: Nine_Four_Node_QuadUP.cpp:177
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
Base class response objects.
Definition: Response.h:81
int sendData(Communicator &)
Send object members through the communicator argument.
Definition: Nine_Four_Node_QuadUP.cpp:858
const Matrix & getMass(void) const
stores the mass matrix in the K member.
Definition: Nine_Four_Node_QuadUP.cpp:501
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Nine_Four_Node_QuadUP.cpp:991
Element * getCopy(void) const
Virtual constructor.
Definition: Nine_Four_Node_QuadUP.cpp:124
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: Nine_Four_Node_QuadUP.cpp:933
Uniaxial p-y material that incorporates liquefaction effects.
Definition: PyLiq1.h:61
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: Nine_Four_Node_QuadUP.cpp:227
void setDomain(Domain *theDomain)
Set the domain for the element.
Definition: Nine_Four_Node_QuadUP.cpp:142
Base class for the finite elements.
Definition: Element.h:112
Uniaxial t-z material that incorporates liquefaction effects.
Definition: TzLiq1.h:61
virtual ~NineFourNodeQuadUP(void)
Destructor.
Definition: Nine_Four_Node_QuadUP.cpp:128
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: Nine_Four_Node_QuadUP.cpp:868
int setParameter(const std::vector< std::string > &argv, Parameter &param)
Sets the value param to the parameter argv.
Definition: Nine_Four_Node_QuadUP.cpp:1030
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
void alive(void)
Reactivates the element.
Definition: Nine_Four_Node_QuadUP.cpp:575
const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: Nine_Four_Node_QuadUP.cpp:369
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: Nine_Four_Node_QuadUP.cpp:773
int recvSelf(const Communicator &)
Receive the object.
Definition: Nine_Four_Node_QuadUP.cpp:893
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
Parameter.
Definition: Parameter.h:68
void zeroLoad(void)
Makes zero the loads on the element.
Definition: Nine_Four_Node_QuadUP.cpp:565
Base class for loads over elements.
Definition: ElementalLoad.h:79
Definition: Nine_Four_Node_QuadUP.h:58
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: Nine_Four_Node_QuadUP.cpp:688
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: Nine_Four_Node_QuadUP.cpp:918
int sendSelf(Communicator &)
Send the object.
Definition: Nine_Four_Node_QuadUP.cpp:878
int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Nine_Four_Node_QuadUP.cpp:1080
int getNumDOF(void) const
Return the number of degrees of freedom.
Definition: Nine_Four_Node_QuadUP.cpp:138