xc
CorotTruss.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.7 $
49 // $Date: 2003/03/12 19:20:46 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/element/truss/CorotTruss.h,v $
51 
52 #ifndef CorotTruss_h
53 #define CorotTruss_h
54 
55 // Written: MHS
56 // Created: May 2001
57 
58 
59 #include "CorotTrussBase.h"
60 #include <utility/matrix/Matrix.h>
61 #include <utility/matrix/Vector.h>
62 
63 namespace XC {
64 class Channel;
65 class UniaxialMaterial;
66 
68 //
74  {
75  private:
76  UniaxialMaterial *theMaterial;
77 
78  double A;
79  double persistentInitialDeformation;
80 
81  protected:
82  void free_material(void);
83  void free_mem(void);
84  void set_material(const UniaxialMaterial &);
85 
86  int sendData(Communicator &comm);
87  int recvData(const Communicator &comm);
88 
89  public:
90  CorotTruss(int tag= 0);
91  CorotTruss(int tag, int dim,int Nd1, int Nd2, const UniaxialMaterial &theMaterial,double A);
92  CorotTruss(int tag,int dimension,const Material *ptr_mat);
93  CorotTruss(const CorotTruss &);
95  Element *getCopy(void) const;
96  ~CorotTruss(void);
97 
98  // public methods to obtain information about dof & connectivity
99  void setDomain(Domain *theDomain);
100 
101  // Element birth and death stuff.
102  const double &getPersistentInitialSectionDeformation(void) const;
104 
105  // public methods to set the state of the element
106  int commitState(void);
107  int revertToLastCommit(void);
108  int revertToStart(void);
109  int update(void);
110 
111  const Material *getMaterial(void) const;
112  Material *getMaterial(void);
113  virtual double getRho(void) const;
114  double getLinearRho(void) const;
115  inline const double &getSectionArea(void) const
116  { return A; }
117  inline void setSectionArea(const double &a)
118  { A= a; }
119 
120  // public methods to obtain stiffness, mass, damping and residual information
121  const Matrix &getTangentStiff(void) const;
122  const Matrix &getInitialStiff(void) const;
123  const Matrix &getMass(void) const;
124 
125  void alive(void);
126  void zeroLoad(void);
127  int addLoad(ElementalLoad *, double loadFactor);
128  int addInertiaLoadToUnbalance(const Vector &accel);
129 
130  double getAxialForce(void) const;
131  const Vector &getResistingForce(void) const;
132  const Vector &getResistingForceIncInertia(void) const;
133  double getInitialStrain(void) const;
134 
135  // public methods for element output
136  int sendSelf(Communicator &);
137  int recvSelf(const Communicator &);
138  void Print(std::ostream &s, int flag =0) const;
139 
140  Response *setResponse(const std::vector<std::string> &argv, Information &eleInfo);
141  int getResponse(int responseID, Information &eleInformation);
142 
143  };
144 } // end of XC namespace
145 
146 #endif
147 
148 
149 
150 
void setDomain(Domain *theDomain)
setDomain() to set a link to the enclosing XC::Domain and to set the node pointers.
Definition: CorotTruss.cpp:174
void alive(void)
Reactivates the element.
Definition: CorotTruss.cpp:470
void incrementPersistentInitialDeformationWithCurrentDeformation(void)
Increments the persistent (does not get wiped out by zeroLoad) initial deformation of the section...
Definition: CorotTruss.cpp:155
int sendData(Communicator &comm)
Send members through the communicator argument.
Definition: CorotTruss.cpp:642
Float vector abstraction.
Definition: Vector.h:94
int addLoad(ElementalLoad *, double loadFactor)
Adds a load.
Definition: CorotTruss.cpp:488
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
int update(void)
Update element state.
Definition: CorotTruss.cpp:305
int recvData(const Communicator &comm)
Receives members through the communicator argument.
Definition: CorotTruss.cpp:651
Base class for corotational truss elements.
Definition: CorotTrussBase.h:41
Base class response objects.
Definition: Response.h:81
const double & getPersistentInitialSectionDeformation(void) const
Returns the value of the persistent (does not get wiped out by zeroLoad) initial deformation of the s...
Definition: CorotTruss.cpp:148
~CorotTruss(void)
Destructor.
Definition: CorotTruss.cpp:161
Base class for uniaxial materials.
Definition: UniaxialMaterial.h:93
CorotTruss(int tag=0)
Constructor.
Definition: CorotTruss.cpp:96
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:613
Base class for materials.
Definition: Material.h:93
Truss element with corotational formulation.
Definition: CorotTruss.h:73
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: CorotTruss.cpp:628
double getInitialStrain(void) const
Return the element initial strain.
Definition: CorotTruss.cpp:529
int revertToStart(void)
Revert the element to its initial state.
Definition: CorotTruss.cpp:297
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: CorotTruss.cpp:447
void zeroLoad(void)
Zeroes loads on element.
Definition: CorotTruss.cpp:481
Base class for the finite elements.
Definition: Element.h:112
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: CorotTruss.cpp:534
int revertToLastCommit(void)
Revert to the last committed state.
Definition: CorotTruss.cpp:287
int recvSelf(const Communicator &)
Receives object through the communicator argument.
Definition: CorotTruss.cpp:675
void set_material(const UniaxialMaterial &)
Assign the material.
Definition: CorotTruss.cpp:80
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: CorotTruss.cpp:566
virtual double getRho(void) const
Return the material density.
Definition: CorotTruss.cpp:440
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: CorotTruss.cpp:596
const Matrix & getTangentStiff(void) const
Return tangent stiffness matrix.
Definition: CorotTruss.cpp:341
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
CorotTruss & operator=(const CorotTruss &)
Assignment operator.
Definition: CorotTruss.cpp:132
Matrix of floats.
Definition: Matrix.h:111
Base class for loads over elements.
Definition: ElementalLoad.h:79
Element * getCopy(void) const
Virtual constructor.
Definition: CorotTruss.cpp:143
int commitState(void)
Commit element state.
Definition: CorotTruss.cpp:273
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
double getLinearRho(void) const
Returns the element mass per unit length.
Definition: CorotTruss.cpp:444
const Matrix & getInitialStiff(void) const
Return initial stiffness matrix.
Definition: CorotTruss.cpp:398
int sendSelf(Communicator &)
Sends object through the communicator argument.
Definition: CorotTruss.cpp:660
void free_material(void)
Free the material pointer.
Definition: CorotTruss.cpp:70