XC Open source finite element analysis program
Tri31.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.00 $
48 // $Date: 2010/09/08 20:01:54 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/triangular/Tri31.h,v $
50 
51 // Written: Roozbeh Geraili Mikola (roozbehg@berkeley.edu)
52 // Created: Sep 2010
53 // Revised: --------
54 //
55 // Description: This file contains the class definition for Tri31.
56 
57 #ifndef Tri31_h
58 #define Tri31_h
59 
60 #include "domain/mesh/element/plane/TriBase3N.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 Tri31: public TriBase3N<SolidMech2D>
74  {
75  private:
76  BodyForces2D bf;
77  Vector pressureLoad;
78 
79  double pressure;
80 
81  mutable Matrix *Ki;
82 
83  static double matrixData[36];
84  static Matrix K;
85  static Vector P;
86 
87  static double shp[3][3];
88 
89  // private member functions - only objects of this class can call these
90  double shapeFunction(const GaussPoint &) const;
91  void setPressureLoadAtNodes(void);
92 
93  protected:
94  int sendData(CommParameters &);
95  int recvData(const CommParameters &);
96 
97  public:
98  Tri31(int tag, int nd1, int nd2, int nd3,
99  NDMaterial &m, const std::string &type,
100  double t, double pressure = 0.0,
101  double rho = 0.0,
102  const BodyForces2D &bForces= BodyForces2D());
103  Tri31(int tag,const NDMaterial *ptr_mat);
104  Tri31(void);
105  Element *getCopy(void) const;
106  virtual ~Tri31(void);
107 
108  int getNumDOF(void) const;
109  void setDomain(Domain *theDomain);
110 
111  // public methods to set the state of the element
112  int update(void);
113 
114  // public methods to obtain stiffness, mass, damping and residual information
115  const Matrix &getTangentStiff(void) const;
116  const Matrix &getInitialStiff(void) const;
117  const Matrix &getMass(void) const;
118 
119  const GaussModel &getGaussModel(void) const;
120 
121  inline double getRho(void) const
122  { return physicalProperties.getRho(); }
123  void setRho(const double &r)
124  { physicalProperties.setRho(r); }
125  double getThickness(void) const
126  { return physicalProperties.getThickness(); }
127  void setThickness(const double &t)
128  { physicalProperties.setThickness(t); }
129 
130  int addLoad(ElementalLoad *theLoad, double loadFactor);
131  int addInertiaLoadToUnbalance(const Vector &accel);
132 
133  const Vector &getResistingForce(void) const;
134  const Vector &getResistingForceIncInertia(void) const;
135 
136  // public methods for element output
137  int sendSelf(CommParameters &);
138  int recvSelf(const CommParameters &);
139  void Print(std::ostream &s, int flag =0);
140 
141  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
142  int getResponse(int responseID, Information &eleInformation);
143 
144  int setParameter(const std::vector<std::string> &argv, Parameter &param);
145  int updateParameter(int parameterID, Information &info);
146 
147  // RWB; PyLiq1 & TzLiq1 need to see the excess pore pressure and initial stresses.
148  friend class PyLiq1;
149  friend class TzLiq1;
150 
151  const Vector &getAvgStress(void) const;
152  const Vector &getAvgStrain(void) const;
153  };
154 
155 } //End namespace XC
156 
157 #endif
158 
??.
Definition: TzLiq1.h:60
int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: Tri31.cpp:520
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
Base class for 3 node triangles.
Definition: TriBase3N.h:46
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: Tri31.cpp:496
3D position of Gauss points.
Definition: GaussPoint.h:37
void Print(std::ostream &s, int flag=0)
Imprime el objeto.
Definition: Tri31.cpp:538
Information about an element.
Definition: Information.h:80
int sendData(CommParameters &)
Send object members through the channel being passed as parameter.
Definition: Tri31.cpp:486
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 the domain for the element.
Definition: Tri31.cpp:112
3 node triangle.
Definition: Tri31.h:73
int update(void)
Actualiza el estado of the element.
Definition: Tri31.cpp:140
Definition: Matrix.h:82
??.
Definition: PyLiq1.h:62
Body forces over an element.
Definition: BodyForces2D.h:39
Element * getCopy(void) const
Virtual constructor.
Definition: Tri31.cpp:97
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Tri31.cpp:302
Definition: Parameter.h:65
Communication parameters between processes.
Definition: CommParameters.h:65
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: Tri31.cpp:436
int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: Tri31.cpp:505
================================================================================
Definition: ContinuaReprComponent.h:34
const GaussModel & getGaussModel(void) const
Returns the puntos de Gauss of the element.
Definition: Tri31.cpp:339
Definition: Response.h:71