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