XC Open source finite element analysis program
FourNodeQuad.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 // $Revision: 1.10 $
48 // $Date: 2003/10/07 21:18:50 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/fourNodeQuad/FourNodeQuad.h,v $
50 
51 // Written: MHS
52 // Created: Feb 2000
53 // Revised: Dec 2000 for efficiency
54 //
55 // Description: This file contains the class definition for FourNodeQuad.
56 
57 #ifndef FourNodeQuad_h
58 #define FourNodeQuad_h
59 
60 #include "domain/mesh/element/plane/QuadBase4N.h"
61 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h"
62 #include "domain/mesh/element/utils/body_forces/BodyForces2D.h"
63 
64 namespace XC {
65 class NDMaterial;
66 class Material;
67 class Response;
68 class GaussPoint;
69 
71 //
73 class FourNodeQuad: public QuadBase4N<SolidMech2D>
74  {
75  private:
76  BodyForces2D bf;
77  Vector pressureLoad;
78 
79  double pressure;
80  mutable Matrix *Ki;
81 
82  static double matrixData[64];
83  static Matrix K;
84  static Vector P;
85  static double shp[3][4];
86 
87  // private member functions - only objects of this class can call these
88  double shapeFunction(const GaussPoint &gp) const;
89  void setPressureLoadAtNodes(void);
90 
91  protected:
92  int sendData(CommParameters &);
93  int recvData(const CommParameters &);
94  bool check_material_type(const std::string &type) const;
95  public:
96  FourNodeQuad(int tag, int nd1, int nd2, int nd3, int nd4, NDMaterial &m, const std::string &type, double t, double pressure = 0.0, double rho = 0.0, const BodyForces2D &bForces= BodyForces2D());
97  FourNodeQuad(int tag,const NDMaterial *ptr_mat);
98  FourNodeQuad(void);
99  Element *getCopy(void) const;
100  virtual ~FourNodeQuad(void);
101 
102  int getNumDOF(void) const;
103  void setDomain(Domain *theDomain);
104 
105  // public methods to set the state of the element
106  int update(void);
107 
108  // public methods to obtain stiffness, mass, damping and residual information
109  const Matrix &getTangentStiff(void) const;
110  const Matrix &getInitialStiff(void) const;
111  const Matrix &getMass(void) const;
112 
113  const GaussModel &getGaussModel(void) const;
114 
115  inline double getRho(void) const
116  { return physicalProperties.getRho(); }
117  void setRho(const double &r)
118  { physicalProperties.setRho(r); }
119  double getThickness(void) const
120  { return physicalProperties.getThickness(); }
121  void setThickness(const double &t)
122  { physicalProperties.setThickness(t); }
123 
124  int addLoad(ElementalLoad *theLoad, double loadFactor);
125  int addInertiaLoadToUnbalance(const Vector &accel);
126 
127  const Vector &getResistingForce(void) const;
128  const Vector &getResistingForceIncInertia(void) const;
129 
130  // public methods for element output
131  int sendSelf(CommParameters &);
132  int recvSelf(const CommParameters &);
133  void Print(std::ostream &s, int flag =0);
134 
135  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
136  int getResponse(int responseID, Information &eleInformation);
137 
138  int setParameter(const std::vector<std::string> &argv, Parameter &param);
139  int updateParameter(int parameterID, Information &info);
140 
141  // RWB; PyLiq1 & TzLiq1 need to see the excess pore pressure and initial stresses.
142  friend class PyLiq1;
143  friend class TzLiq1;
144 
145  double detJ(const double &,const double &) const;
146 
147  };
148 } // end of XC namespace
149 
150 #endif
151 
FourNodeQuad(void)
Constructor.
Definition: FourNodeQuad.cpp:94
??.
Definition: TzLiq1.h:60
int addLoad(ElementalLoad *theLoad, double loadFactor)
Adds a load over element.
Definition: FourNodeQuad.cpp:336
int getNumDOF(void) const
Return the number of element DOFs.
Definition: FourNodeQuad.cpp:118
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: FourNodeQuad.cpp:501
const Vector & getResistingForceIncInertia(void) const
Return the fuerza de respuesta of the element incluyendo la debida a la inercia.
Definition: FourNodeQuad.cpp:438
int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: FourNodeQuad.cpp:511
int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: FourNodeQuad.cpp:526
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:98
Base class for Gauss integration models.
Definition: GaussModel.h:40
Base class for 2D and 3D materials.
Definition: NDMaterial.h:91
Definition: Vector.h:82
int update(void)
Actualiza los valores de las variables de estado.
Definition: FourNodeQuad.cpp:134
3D position of Gauss points.
Definition: GaussPoint.h:37
Four node quad.
Definition: FourNodeQuad.h:73
const Matrix & getMass(void) const
Return the matriz de masas.
Definition: FourNodeQuad.cpp:292
virtual ~FourNodeQuad(void)
Destructor.
Definition: FourNodeQuad.cpp:104
Element * getCopy(void) const
Virtual constructor.
Definition: FourNodeQuad.cpp:100
void Print(std::ostream &s, int flag=0)
Prints element information.
Definition: FourNodeQuad.cpp:546
Information about an element.
Definition: Information.h:80
bool check_material_type(const std::string &type) const
Checks the material type.
Definition: FourNodeQuad.cpp:114
int addInertiaLoadToUnbalance(const Vector &accel)
Adds inertia loads.
Definition: FourNodeQuad.cpp:343
const Matrix & getInitialStiff(void) const
Return the initial tangent stiffness matrix.
Definition: FourNodeQuad.cpp:238
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
SolidMech2D physicalProperties
pointers to the material objects and physical properties.
Definition: ElemWithMaterial.h:43
void setDomain(Domain *theDomain)
Sets domain pointer and computes the consistent load vector due to pressure.
Definition: FourNodeQuad.cpp:122
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: FourNodeQuad.cpp:180
Base class for 4 node quads.
Definition: QuadBase4N.h:45
Definition: Matrix.h:82
??.
Definition: PyLiq1.h:62
Body forces over an element.
Definition: BodyForces2D.h:39
Definition: Parameter.h:65
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: FourNodeQuad.cpp:491
Communication parameters between processes.
Definition: CommParameters.h:65
================================================================================
Definition: ContinuaReprComponent.h:34
const GaussModel & getGaussModel(void) const
Returns the puntos de Gauss of the element.
Definition: FourNodeQuad.cpp:332
const Vector & getResistingForce(void) const
Return the fuerza de respuesta of the element.
Definition: FourNodeQuad.cpp:388
Definition: Response.h:71