xc
CorotTruss.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.7 $
48 // $Date: 2003/03/12 19:20:46 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/truss/CorotTruss.h,v $
50 
51 #ifndef CorotTruss_h
52 #define CorotTruss_h
53 
54 // Written: MHS
55 // Created: May 2001
56 
57 
58 #include "CorotTrussBase.h"
59 #include <utility/matrix/Matrix.h>
60 #include <utility/matrix/Vector.h>
61 
62 namespace XC {
63 class Channel;
64 class UniaxialMaterial;
65 
67 //
73  {
74  private:
75  UniaxialMaterial *theMaterial;
76 
77  double A;
78  protected:
79  int sendData(CommParameters &cp);
80  int recvData(const CommParameters &cp);
81 
82  public:
83  CorotTruss(int tag, int dim,int Nd1, int Nd2, UniaxialMaterial &theMaterial,double A);
84  CorotTruss(int tag,int dimension,const Material *ptr_mat);
85  CorotTruss(void);
86  CorotTruss(const CorotTruss &);
88  Element *getCopy(void) const;
89  ~CorotTruss(void);
90 
91  // public methods to obtain information about dof & connectivity
92  void setDomain(Domain *theDomain);
93 
94  // public methods to set the state of the element
95  int commitState(void);
96  int revertToLastCommit(void);
97  int revertToStart(void);
98  int update(void);
99 
100  const Material *getMaterial(void) const;
101  Material *getMaterial(void);
102  virtual double getRho(void) const;
103  inline const double &getArea(void) const
104  { return A; }
105  inline void setArea(const double &a)
106  { A= a; }
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  void zeroLoad(void);
114  int addLoad(ElementalLoad *theLoad, double loadFactor);
115  int addInertiaLoadToUnbalance(const Vector &accel);
116 
117  double getAxialForce(void) const;
118  const Vector &getResistingForce(void) const;
119  const Vector &getResistingForceIncInertia(void) const;
120 
121  // public methods for element output
122  int sendSelf(CommParameters &);
123  int recvSelf(const CommParameters &);
124  void Print(std::ostream &s, int flag =0);
125 
126  Response *setResponse(const std::vector<std::string> &argv, Information &eleInfo);
127  int getResponse(int responseID, Information &eleInformation);
128 
129  };
130 } // end of XC namespace
131 
132 #endif
133 
134 
135 
136 
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: CorotTruss.cpp:139
Float vector abstraction.
Definition: Vector.h:93
Information about an element.
Definition: Information.h:80
int update(void)
Update element state.
Definition: CorotTruss.cpp:264
Base class for corotational truss elements.
Definition: CorotTrussBase.h:40
Definition: Response.h:71
int sendData(CommParameters &cp)
Send members through the channel being passed as parameter.
Definition: CorotTruss.cpp:576
int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: CorotTruss.cpp:609
~CorotTruss(void)
Destructor.
Definition: CorotTruss.cpp:125
Base class for uniaxial materials.
Definition: UniaxialMaterial.h:92
Response * setResponse(const std::vector< std::string > &argv, Information &eleInfo)
setResponse() is a method invoked to determine if the element will respond to a request for a certain...
Definition: CorotTruss.cpp:544
Base class for materials.
Definition: Material.h:91
Truss element with corotational formulation.
Definition: CorotTruss.h:72
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: CorotTruss.cpp:562
int revertToStart(void)
Revert the element to its initial state.
Definition: CorotTruss.cpp:258
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: CorotTruss.cpp:400
void zeroLoad(void)
Zeroes loads on element.
Definition: CorotTruss.cpp:423
Base class for the finite elements.
Definition: Element.h:109
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: CorotTruss.cpp:465
int revertToLastCommit(void)
Revert to the last commited state.
Definition: CorotTruss.cpp:252
int recvData(const CommParameters &cp)
Receives members through the channel being passed as parameter.
Definition: CorotTruss.cpp:585
void Print(std::ostream &s, int flag=0)
Print stuff.
Definition: CorotTruss.cpp:527
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: CorotTruss.cpp:497
virtual double getRho(void) const
Return the material density.
Definition: CorotTruss.cpp:397
const Matrix & getTangentStiff(void) const
Return tangent stiffness matrix.
Definition: CorotTruss.cpp:299
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
CorotTruss & operator=(const CorotTruss &)
Assignment operator.
Definition: CorotTruss.cpp:107
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
int addLoad(ElementalLoad *theLoad, double loadFactor)
Adds a load.
Definition: CorotTruss.cpp:430
Base class for loads over elements.
Definition: ElementalLoad.h:77
Element * getCopy(void) const
Virtual constructor.
Definition: CorotTruss.cpp:119
int commitState(void)
Commit element state.
Definition: CorotTruss.cpp:238
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:116
int sendSelf(CommParameters &)
Sends object through the channel being passed as parameter.
Definition: CorotTruss.cpp:594
const Matrix & getInitialStiff(void) const
Return initial stiffness matrix.
Definition: CorotTruss.cpp:356