xc
CorotCrdTransf3d.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: 2006/01/13 01:07:48 $
50 // $Source: /usr/local/cvs/OpenSees/SRC/CorotCrdTransf3d.h,v $
51 
52 
53 // Written: Remo Magalhaes de Souza (rmsouza@ce.berkeley.edu)
54 // Created: 04/2000
55 // Revision: A
56 //
57 // Description: This file contains the class definition for
58 // CorotCrdTransf3d.h. CorotCrdTransf3d provides the
59 // abstraction of a corotation transformation for a spatial frame element
60 
61 // What: "@(#) CorotCrdTransf3d.h, revA"
62 
63 #ifndef CorotCrdTransf3d_h
64 #define CorotCrdTransf3d_h
65 
66 #include "CrdTransf3d.h"
67 #include <utility/matrix/Vector.h>
68 #include <utility/matrix/Matrix.h>
69 
70 namespace XC {
72 //
75  {
76  private:
77  Vector vAxis;
78 
79  mutable Vector xAxis;
80  mutable double Ln;
81 
82  mutable Vector alphaIq;
83  mutable Vector alphaJq;
84 
85  Vector alphaIqcommit;
86  Vector alphaJqcommit;
87  mutable Vector alphaI;
88  mutable Vector alphaJ;
89 
90  Vector ul;
91  Vector ulcommit;
92  Vector ulpr;
93 
94  static Matrix RI;
95  static Matrix RJ;
96  static Matrix Rbar;
97  static Matrix e;
98  static Matrix Tp;
99  static Matrix T;
100  static Matrix Tlg;
101  static Matrix TlgInv;
102  static Matrix Tbl;
103  static Matrix kg;
104  static Matrix Lr2, Lr3, A;
105 
106  inline int computeElemtLengthAndOrient(void) const
107  {
108  std::cerr << "CorotCrdTransf3d::computeElemtLengthAndOrient; not implemented."
109  << std::endl;
110  return 0;
111  }
112  void compTransfMatrixBasicGlobal(void);
113  void compTransfMatrixBasicGlobalNew(void);
114  void compTransfMatrixLocalGlobal(Matrix &Tlg) const;
115  void compTransfMatrixBasicLocal(Matrix &Tbl) const;
116  const Vector &getQuaternionFromRotMatrix(const Matrix &RotMatrix) const;
117  const Vector &getQuaternionFromPseudoRotVector(const Vector &theta) const;
118  const Vector &getTangScaledPseudoVectorFromQuaternion(const Vector &theta) const;
119  const Vector &quaternionProduct(const Vector &q1, const Vector &q2) const;
120  const Matrix &getRotationMatrixFromQuaternion(const Vector &q) const;
121  const Matrix &getRotMatrixFromTangScaledPseudoVector(const Vector &w) const;
122  const Matrix &getSkewSymMatrix(const Vector &theta) const;
123  const Matrix &getLMatrix(const Vector &ri) const;
124  const Matrix &getKs2Matrix(const Vector &ri, const Vector &z) const;
125 
126 
127 
128  Vector &basic_to_local_element_force(const Vector &p0) const;
129  const Vector &local_to_global_element_force(const Vector &) const;
130  protected:
131  int sendData(Communicator &);
132  int recvData(const Communicator &);
133  virtual int computeLocalAxis(void) const;
134  public:
135  CorotCrdTransf3d(int tag= 0);
136  CorotCrdTransf3d(int tag, const Vector &vecInLocXZPlane,const Vector &rigJntOffsetI, const Vector &rigJntOffsetJ);
137 
138  int initialize(Node *nodeIPointer, Node *nodeJPointer);
139  int update(void);
140  double getInitialLength(void) const;
141  double getDeformedLength(void) const;
142 
143  virtual void set_xz_vector(const Vector &vecInLocXZPlane);
144 
145  int commitState(void);
146  int revertToLastCommit(void);
147  int revertToStart(void);
148 
149  const Vector &getBasicTrialDisp(void) const;
150  const Vector &getBasicIncrDisp(void) const;
151  const Vector &getBasicIncrDeltaDisp(void) const;
152  const Vector &getBasicTrialVel(void) const;
153  const Vector &getBasicTrialAccel(void) const;
154 
155  const Vector &getGlobalResistingForce(const Vector &basicForce, const Vector &uniformLoad) const;
156  const Matrix &getGlobalStiffMatrix(const Matrix &basicStiff, const Vector &basicForce) const;
157  const Matrix &getInitialGlobalStiffMatrix(const Matrix &basicStiff) const;
158 
159  CrdTransf3d *getCopy(void) const;
160 
161  int sendSelf(Communicator &);
162  int recvSelf(const Communicator &);
163 
164  void Print(std::ostream &s, int flag = 0) const;
165 
166  // method used to rotate consistent mass matrix
167  const Matrix &getGlobalMatrixFromLocal(const Matrix &local);
168 
169  // functions used in post-processing only
170  const Vector &getPointGlobalCoordFromLocal(const Vector &) const;
171  const Vector &getPointGlobalDisplFromBasic(double xi, const Vector &basicDisps) const;
172  const Vector &getVectorGlobalCoordFromLocal(const Vector &localCoords) const;
173  const Matrix &getVectorGlobalCoordFromLocal(const Matrix &localCoords) const;
174  const Vector &getVectorLocalCoordFromGlobal(const Vector &globalCoords) const;
175  };
176 } // end of XC namespace
177 #endif
Float vector abstraction.
Definition: Vector.h:94
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: CorotCrdTransf3d.cpp:2094
Communication parameters between processes.
Definition: Communicator.h:66
int recvSelf(const Communicator &)
Receive the object.
Definition: CorotCrdTransf3d.cpp:2056
const Vector & getVectorLocalCoordFromGlobal(const Vector &globalCoords) const
Returns the local coordinates of the vector.
Definition: CorotCrdTransf3d.cpp:1581
int recvData(const Communicator &)
Receives object members through the communicator argument.
Definition: CorotCrdTransf3d.cpp:2021
CrdTransf3d * getCopy(void) const
Virtual constructor.
Definition: CorotCrdTransf3d.cpp:1994
int sendData(Communicator &)
Sends object members through the communicator argument.
Definition: CorotCrdTransf3d.cpp:1998
int sendSelf(Communicator &)
Send the object.
Definition: CorotCrdTransf3d.cpp:2042
virtual void set_xz_vector(const Vector &vecInLocXZPlane)
Set the vector that defines the local XZ plane.
Definition: CorotCrdTransf3d.cpp:203
const Vector & getVectorGlobalCoordFromLocal(const Vector &localCoords) const
Returns the global coordinates of the vector.
Definition: CorotCrdTransf3d.cpp:1556
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
Mesh node.
Definition: Node.h:111
Coordinate transformation corrotacional en 3d.
Definition: CorotCrdTransf3d.h:74
Base class for 3D coordinate transformation.
Definition: CrdTransf3d.h:81