xc
Tri31.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 //----------------------------------------------------------------------------
28 /* ****************************************************************** **
29 ** OpenSees - Open System for Earthquake Engineering Simulation **
30 ** Pacific Earthquake Engineering Research Center **
31 ** **
32 ** **
33 ** (C) Copyright 1999, The Regents of the University of California **
34 ** All Rights Reserved. **
35 ** **
36 ** Commercial use of this program without express permission of the **
37 ** University of California, Berkeley, is strictly prohibited. See **
38 ** file 'COPYRIGHT' in main directory for information on usage and **
39 ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. **
40 ** **
41 ** Developed by: **
42 ** Frank McKenna (fmckenna@ce.berkeley.edu) **
43 ** Gregory L. Fenves (fenves@ce.berkeley.edu) **
44 ** Filip C. Filippou (filippou@ce.berkeley.edu) **
45 ** **
46 ** ****************************************************************** */
47 
48 // $Revision: 1.00 $
49 // $Date: 2010/09/08 20:01:54 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/element/triangular/Tri31.h,v $
51 
52 // Written: Roozbeh Geraili Mikola (roozbehg@berkeley.edu)
53 // Created: Sep 2010
54 // Revised: --------
55 //
56 // Description: This file contains the class definition for Tri31.
57 
58 #ifndef Tri31_h
59 #define Tri31_h
60 
61 #include "domain/mesh/element/plane/TriBase3N.h"
62 #include "domain/mesh/element/utils/physical_properties/SolidMech2D.h"
63 #include "domain/mesh/element/utils/body_forces/BodyForces2D.h"
64 
65 namespace XC {
66 class NDMaterial;
67 class Material;
68 class Response;
69 class GaussPoint;
70 
72 //
74 class Tri31: public TriBase3N<SolidMech2D>
75  {
76  private:
77  BodyForces2D bf;
78  Vector pressureLoad;
79 
80  double pressure;
81 
82  mutable Matrix *Ki;
83 
84  static double matrixData[36];
85  static Matrix K;
86  static Vector P;
87 
88  static double shp[3][3];
89 
90  // private member functions - only objects of this class can call these
91  double shapeFunction(const GaussPoint &) const;
92  void setPressureLoadAtNodes(void);
93 
94  protected:
95  int sendData(Communicator &);
96  int recvData(const Communicator &);
97 
98  public:
99  Tri31(int tag, int nd1, int nd2, int nd3,
100  NDMaterial &m, const std::string &type,
101  double t, double pressure = 0.0,
102  double rho = 0.0,
103  const BodyForces2D &bForces= BodyForces2D());
104  Tri31(int tag,const NDMaterial *ptr_mat);
105  Tri31(void);
106  Element *getCopy(void) const;
107  virtual ~Tri31(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 &getMass(void) const;
119 
120  const GaussModel &getGaussModel(void) const;
121 
122  inline double getRho(void) const
123  { return physicalProperties.getRho(); }
124  void setRho(const double &r)
125  { physicalProperties.setRho(r); }
126  double getThickness(void) const
127  { return physicalProperties.getThickness(); }
128  void setThickness(const double &t)
130 
131  void alive(void);
132  int addLoad(ElementalLoad *theLoad, double loadFactor);
133  int addInertiaLoadToUnbalance(const Vector &accel);
134 
135  const Vector &getResistingForce(void) const;
136  const Vector &getResistingForceIncInertia(void) const;
137 
138  // public methods for element output
139  int sendSelf(Communicator &);
140  int recvSelf(const Communicator &);
141  void Print(std::ostream &s, int flag =0) const;
142 
143  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
144  int getResponse(int responseID, Information &eleInformation);
145 
146  int setParameter(const std::vector<std::string> &argv, Parameter &param);
147  int updateParameter(int parameterID, Information &info);
148 
149  // RWB; PyLiq1 & TzLiq1 need to see the excess pore pressure and initial stresses.
150  friend class PyLiq1;
151  friend class TzLiq1;
152 
153  const Vector &getAvgStress(void) const;
154  const Vector &getAvgStrain(void) const;
155  };
156 
157 } //End namespace XC
158 
159 #endif
160 
int updateParameter(int parameterID, Information &info)
Updates the parameter identified by parameterID with info.
Definition: Tri31.cpp:700
int update(void)
Updates the element state.
Definition: Tri31.cpp:140
Body forces over an element.
Definition: BodyForces2D.h:40
Element * getCopy(void) const
Virtual constructor.
Definition: Tri31.cpp:97
Float vector abstraction.
Definition: Vector.h:94
const GaussModel & getGaussModel(void) const
Return the Gauss points of the element.
Definition: Tri31.cpp:339
int setParameter(const std::vector< std::string > &argv, Parameter &param)
Sets the value param to the parameter argv.
Definition: Tri31.cpp:673
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: Tri31.cpp:302
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: Tri31.cpp:184
Base class response objects.
Definition: Response.h:81
double getThickness(void) const
Return the material thickness.
Definition: SolidMech2D.h:55
int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: Tri31.cpp:522
Uniaxial p-y material that incorporates liquefaction effects.
Definition: PyLiq1.h:61
SolidMech2D physicalProperties
pointers to the material objects and physical properties.
Definition: ElemWithMaterial.h:50
3D position of Gauss points.
Definition: GaussPoint.h:38
void setRho(const double &)
Set the density for all materials.
Definition: NDMaterialPhysicalProperties.cc:197
Base class for the finite elements.
Definition: Element.h:112
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: Tri31.cpp:109
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: Tri31.cpp:649
int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: Tri31.cpp:537
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: Tri31.cpp:555
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: Tri31.cpp:112
Base class for 3 node triangles.
Definition: TriBase3N.h:46
Uniaxial t-z material that incorporates liquefaction effects.
Definition: TzLiq1.h:61
3 node triangle.
Definition: Tri31.h:74
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: Tri31.cpp:604
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: Tri31.cpp:405
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: Tri31.cpp:453
int sendData(Communicator &)
Send object members through the communicator argument.
Definition: Tri31.cpp:503
Base class for Gauss integration models.
Definition: GaussModel.h:41
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
void alive(void)
Reactivates the element.
Definition: Tri31.cpp:343
double getRho(void) const
Returns the average of the densities for each material.
Definition: NDMaterialPhysicalProperties.cc:186
void setThickness(const double &t)
Set the material thickness.
Definition: SolidMech2D.h:58
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
Parameter.
Definition: Parameter.h:68
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: Tri31.cpp:513