xc
CorotTrussSection.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.6 $
49 // $Date: 2003/03/12 19:20:46 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/element/truss/CorotTrussSection.h,v $
51 
52 #ifndef CorotTrussSection_h
53 #define CorotTrussSection_h
54 
55 // Written: MHS
56 // Created: May 2001
57 
58 #include "CorotTrussBase.h"
59 #include <utility/matrix/Matrix.h>
60 #include <utility/matrix/Vector.h>
61 #include "domain/mesh/element/utils/physical_properties/SectionFDPhysicalProperties.h"
62 
63 namespace XC {
64 
66 //
73  {
74  private:
75  SectionFDPhysicalProperties physicalProperties;
76  double persistentInitialDeformation;
77 
78  double computeCurrentStrain(void) const;
79  public:
80  CorotTrussSection(int tag, int dim,int Nd1, int Nd2, SectionForceDeformation &theMaterial);
81  CorotTrussSection(int tag,int dimension,const Material *ptr_mat);
82  CorotTrussSection(void);
83  Element *getCopy(void) const;
84  inline virtual ~CorotTrussSection(void) {}
85 
86  void setDomain(Domain *theDomain);
87 
88  // Element birth and death stuff.
89  const double &getPersistentInitialSectionDeformation(void) const;
91 
92  // public methods to set the state of the element
93  int commitState(void);
94  int revertToLastCommit(void);
95  int revertToStart(void);
96  int update(void);
97 
98  const Material *getMaterial(void) const;
99  Material *getMaterial(void);
100  virtual double getRho(void) const;
101  double getLinearRho(void) const;
102 
103  // public methods to obtain stiffness, mass, damping and residual information
104  const Matrix &getTangentStiff(void) const;
105  const Matrix &getInitialStiff(void) const;
106  const Matrix &getMass(void) const;
107 
108  void alive(void);
109  void zeroLoad(void);
110  int addLoad(ElementalLoad *, double loadFactor);
111  int addInertiaLoadToUnbalance(const Vector &accel);
112 
113  double getAxialForce(void) const;
114  const Vector &getResistingForce(void) const;
115  const Vector &getResistingForceIncInertia(void) const;
116  double getInitialStrain(void) const;
117 
118  // public methods for element output
119  int sendSelf(Communicator &);
120  int recvSelf(const Communicator &);
121  void Print(std::ostream &s, int flag =0) const;
122 
123  Response *setResponse(const std::vector<std::string> &argv, Information &eleInfo);
124  int getResponse(int responseID, Information &eleInformation);
125  };
126 } // end of XC namespace
127 
128 #endif
129 
130 
131 
132 
double getInitialStrain(void) const
Return the element initial strain.
Definition: CorotTrussSection.cpp:512
Base class for force deformation section models.
Definition: SectionForceDeformation.h:88
Float vector abstraction.
Definition: Vector.h:94
Truss element with corotational formulation and material of type SectionForceDeformation.
Definition: CorotTrussSection.h:72
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
double getAxialForce(void) const
Return the axial internal force.
Definition: CorotTrussSection.cpp:494
Base class for corotational truss elements.
Definition: CorotTrussBase.h:41
Base class response objects.
Definition: Response.h:81
int sendSelf(Communicator &)
Send the object.
Definition: CorotTrussSection.cpp:585
int update(void)
Update element state.
Definition: CorotTrussSection.cpp:271
const Vector & getResistingForce(void) const
Return the element resisting force.
Definition: CorotTrussSection.cpp:528
Base class for materials.
Definition: Material.h:93
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:607
void incrementPersistentInitialDeformationWithCurrentDeformation(void)
Increments the persistent (does not get wiped out by zeroLoad) initial deformation of the section...
Definition: CorotTrussSection.cpp:105
int addInertiaLoadToUnbalance(const Vector &accel)
Add inertia load to the out of balance vector.
Definition: CorotTrussSection.cpp:488
int revertToLastCommit(void)
Revert the element to its last committed state.
Definition: CorotTrussSection.cpp:223
void zeroLoad(void)
Make loads zero.
Definition: CorotTrussSection.cpp:454
Base class for the finite elements.
Definition: Element.h:112
Element * getCopy(void) const
Virtual constructor.
Definition: CorotTrussSection.cpp:92
virtual double getRho(void) const
Return the linear density of the section.
Definition: CorotTrussSection.cpp:411
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: CorotTrussSection.cpp:616
int revertToStart(void)
Revert the element to its initial state.
Definition: CorotTrussSection.cpp:234
CorotTrussSection(void)
Constructor:
Definition: CorotTrussSection.cpp:86
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: CorotTrussSection.cpp:557
double getLinearRho(void) const
Returns the mass per unit length.
Definition: CorotTrussSection.cpp:415
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: CorotTrussSection.cpp:591
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: CorotTrussSection.cpp:294
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Physical properties for shells.
Definition: SectionFDPhysicalProperties.h:41
Matrix of floats.
Definition: Matrix.h:111
const double & getPersistentInitialSectionDeformation(void) const
Returns the value of the persistent (does not get wiped out by zeroLoad) initial deformation of the s...
Definition: CorotTrussSection.cpp:98
void setDomain(Domain *theDomain)
setDomain() to set a link to the enclosing XC::Domain and to set the node pointers.
Definition: CorotTrussSection.cpp:114
int recvSelf(const Communicator &)
Receive the object.
Definition: CorotTrussSection.cpp:588
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
void alive(void)
Reactivates the element.
Definition: CorotTrussSection.cpp:443
int commitState(void)
Commit element state.
Definition: CorotTrussSection.cpp:209
const Matrix & getMass(void) const
Return the mass matrix.
Definition: CorotTrussSection.cpp:419