xc
LinearCrdTransf3d.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.9 $
49 // $Date: 2005/12/15 00:30:38 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/LinearCrdTransf3d.h,v $
51 
52 // Written: Remo Magalhaes de Souza (rmsouza@ce.berkeley.edu)
53 // Created: 04/2000
54 // Revision: A
55 //
56 // Description: This file contains the class definition for
57 // LinearCrdTransf3d.h. LinearCrdTransf3d provides the
58 // abstraction of a linear transformation for a spatial frame
59 // between the global and basic coordinate systems
60 
61 // What: "@(#) LinearCrdTransf3d.h, revA"
62 
63 #ifndef LinearCrdTransf3d_h
64 #define LinearCrdTransf3d_h
65 
66 #include "SmallDispCrdTransf3d.h"
67 
68 namespace XC {
69 
71 //
76  {
77  public:
78  LinearCrdTransf3d(int tag);
79  LinearCrdTransf3d(int tag, const Vector &vecInLocXZPlane);
80  LinearCrdTransf3d(int tag, const Vector &vecInLocXZPlane, const Vector &rigJntOffsetI, const Vector &rigJntOffsetJ);
81 
82  LinearCrdTransf3d(void);
83 
84  int update(void);
85 
86  int commitState(void);
87  int revertToLastCommit(void);
88  int revertToStart(void);
89 
90  const Vector &getGlobalResistingForce(const Vector &basicForce, const Vector &p0) const;
91  const Matrix &getGlobalStiffMatrix(const Matrix &basicStiff, const Vector &basicForce) const;
92 
93  CrdTransf3d *getCopy(void) const;
94 
95  void Print(std::ostream &s, int flag = 0) const;
96 
97  // functions used in post-processing only
98  const Vector &getPointGlobalCoordFromLocal(const Vector &) const;
99  const Vector &getPointGlobalDisplFromBasic(double xi, const Vector &) const;
100  };
101 } // end of XC namespace
102 
103 #endif
104 
Float vector abstraction.
Definition: Vector.h:94
CrdTransf3d * getCopy(void) const
Virtual constructor.
Definition: LinearCrdTransf3d.cpp:132
const Matrix & getGlobalStiffMatrix(const Matrix &basicStiff, const Vector &basicForce) const
Returns the stiffness matrix expressed in global coordinates.
Definition: LinearCrdTransf3d.cpp:125
int commitState(void)
Commits state of the coordinate transformation.
Definition: LinearCrdTransf3d.cpp:97
Linear coordinate transformation.
Definition: LinearCrdTransf3d.h:75
int revertToLastCommit(void)
Returns to the last committed state.
Definition: LinearCrdTransf3d.cpp:102
const Vector & getGlobalResistingForce(const Vector &basicForce, const Vector &p0) const
Returns the load vector expressed in global coordinates.
Definition: LinearCrdTransf3d.cpp:117
Base class for small displacements 3D coordinate transformations.
Definition: SmallDispCrdTransf3d.h:40
LinearCrdTransf3d(void)
Constructor; invoked by a FEM_ObjectBroker, recvSelf() needs to be invoked on this object...
Definition: LinearCrdTransf3d.cpp:92
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: LinearCrdTransf3d.cpp:206
int update(void)
Updates the transformation.
Definition: LinearCrdTransf3d.cpp:113
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
int revertToStart(void)
Returns to the initial state.
Definition: LinearCrdTransf3d.cpp:106
Base class for 3D coordinate transformation.
Definition: CrdTransf3d.h:81