xc
CorotTrussSection.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.6 $
48 // $Date: 2003/03/12 19:20:46 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/element/truss/CorotTrussSection.h,v $
50 
51 #ifndef CorotTrussSection_h
52 #define CorotTrussSection_h
53 
54 // Written: MHS
55 // Created: May 2001
56 
57 #include "CorotTrussBase.h"
58 #include <utility/matrix/Matrix.h>
59 #include <utility/matrix/Vector.h>
60 
61 namespace XC {
62 class Channel;
63 class SectionForceDeformation;
64 
66 //
73  {
74  private:
75  SectionForceDeformation *theSection; // pointer to a material
76  public:
77  CorotTrussSection(int tag, int dim,int Nd1, int Nd2, SectionForceDeformation &theMaterial);
78  CorotTrussSection(int tag,int dimension,const Material *ptr_mat);
80  CorotTrussSection(void);
82  Element *getCopy(void) const;
83  ~CorotTrussSection(void);
84 
85  void setDomain(Domain *theDomain);
86 
87  // public methods to set the state of the element
88  int commitState(void);
89  int revertToLastCommit(void);
90  int revertToStart(void);
91  int update(void);
92 
93  const Material *getMaterial(void) const;
94  Material *getMaterial(void);
95  virtual double getRho(void) const;
96 
97  // public methods to obtain stiffness, mass, damping and residual information
98  const Matrix &getTangentStiff(void) const;
99  const Matrix &getInitialStiff(void) const;
100  const Matrix &getMass(void) const;
101 
102  int addLoad(ElementalLoad *theLoad, double loadFactor);
103  int addInertiaLoadToUnbalance(const Vector &accel);
104 
105  const Vector &getResistingForce(void) const;
106  const Vector &getResistingForceIncInertia(void) const;
107 
108  // public methods for element output
109  int sendSelf(CommParameters &);
110  int recvSelf(const CommParameters &);
111  void Print(std::ostream &s, int flag =0);
112 
113  Response *setResponse(const std::vector<std::string> &argv, Information &eleInfo);
114  int getResponse(int responseID, Information &eleInformation);
115  };
116 } // end of XC namespace
117 
118 #endif
119 
120 
121 
122 
Base class for force deformation section models.
Definition: SectionForceDeformation.h:86
Float vector abstraction.
Definition: Vector.h:93
Truss element with corotatinal formulation and material of type SectionForceDeformation.
Definition: CorotTrussSection.h:72
Information about an element.
Definition: Information.h:80
int recvSelf(const CommParameters &)
Receive the object.
Definition: CorotTrussSection.cpp:531
Base class for corotational truss elements.
Definition: CorotTrussBase.h:40
Definition: Response.h:71
int update(void)
Update element state.
Definition: CorotTrussSection.cpp:267
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: CorotTrussSection.cpp:460
Base class for materials.
Definition: Material.h:91
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: CorotTrussSection.cpp:549
int revertToLastCommit(void)
Revert the element to its last commited state.
Definition: CorotTrussSection.cpp:253
Base class for the finite elements.
Definition: Element.h:109
Element * getCopy(void) const
Virtual constructor.
Definition: CorotTrussSection.cpp:124
virtual double getRho(void) const
Return the density of the section.
Definition: CorotTrussSection.cpp:424
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: CorotTrussSection.cpp:558
CorotTrussSection & operator=(const CorotTrussSection &)
Assignment operator.
Definition: CorotTrussSection.cpp:113
int revertToStart(void)
Revert the element to its initial state.
Definition: CorotTrussSection.cpp:260
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: CorotTrussSection.cpp:500
void Print(std::ostream &s, int flag=0)
Print stuff.
Definition: CorotTrussSection.cpp:534
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: CorotTrussSection.cpp:309
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: CorotTrussSection.cpp:144
Base class for loads over elements.
Definition: ElementalLoad.h:77
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:116
int sendSelf(CommParameters &)
Send the object.
Definition: CorotTrussSection.cpp:528
int commitState(void)
Commit element state.
Definition: CorotTrussSection.cpp:239
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: CorotTrussSection.cpp:427