XC Open source finite element analysis program
CorotCrdTransf2d.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.8 $
48 // $Date: 2006/01/13 01:07:48 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/CorotCrdTransf2d.h,v $
50 
51 // Written: Remo Magalhaes de Souza (rmsouza@ce.berkeley.edu)
52 // Created: 04/2000
53 // Revision: A
54 //
55 // Description: This file contains the class definition for
56 // CorotCrdTransf2d.h. CorotCrdTransf2d provides the
57 // abstraction of a corotation transformation for a planar frame element
58 
59 // What: "@(#) CorotCrdTransf2d.h, revA"
60 
61 #ifndef CorotCrdTransf2d_h
62 #define CorotCrdTransf2d_h
63 
64 #include "CrdTransf2d.h"
65 #include <utility/matrix/Matrix.h>
66 
67 namespace XC {
69 //
72  {
73  private:
74  double cosAlpha, sinAlpha; // direction cossines of deformed element wrt to local system
75  double Ln; // deformed element length
76  double Lx, Ly; // components of the deformed member
77  mutable double Lxdot, Lydot;
78  mutable double Lxdotdot, Lydotdot;
79 
80  Vector ub; // basic displacements
81  Vector ubcommit; // commited basic displacements
82  Vector ubpr; // previous basic displacements
83 
84  bool nodeOffsets;
85 
86  static Matrix Tlg; // matrix that transforms from global to local coordinates
87  static Matrix Tbl; // matrix that transforms from local to basic coordinates
88  static Matrix kg;
89  static Vector uxg;
90  static Vector pg;
91  static Vector dub;
92  static Vector Dub;
93 
94  int compElemtLengthAndOrient(void);
95  int compElemtLengthAndOrientWRTLocalSystem(const Vector &ul);
96  void transfLocalDisplsToBasic(const Vector &ul);
97  void getTransfMatrixLocalGlobal(Matrix &Tlg) const;
98  void getTransfMatrixBasicLocal(Matrix &Tbl) const;
99  const Matrix &getGeomStiffMatrix(const Vector &pb) const;
100  protected:
101  int sendData(CommParameters &);
102  int recvData(const CommParameters &);
103  public:
104  CorotCrdTransf2d(int tag= 0);
105  CorotCrdTransf2d(int tag, const Vector &rigJntOffsetI, const Vector &rigJntOffsetJ);
106 
107  int initialize(Node *nodeIPointer, Node *nodeJPointer);
108  int update(void);
109  double getInitialLength(void) const;
110  double getDeformedLength(void) const;
111 
112  int commitState(void);
113  int revertToLastCommit(void);
114  int revertToStart(void);
115 
116  const Vector &getBasicTrialDisp(void) const;
117  const Vector &getBasicIncrDisp(void) const;
118  const Vector &getBasicIncrDeltaDisp(void) const;
119  const Vector &getBasicTrialVel(void) const;
120  const Vector &getBasicTrialAccel(void) const;
121 
122  const Vector &getGlobalResistingForce(const Vector &basicForce, const Vector &uniformLoad) const;
123  const Matrix &getGlobalStiffMatrix(const Matrix &basicStiff, const Vector &basicForce) const;
124  const Matrix &getInitialGlobalStiffMatrix(const Matrix &basicStiff) const;
125 
126 // AddingSensitivity:BEGIN //////////////////////////////////
127  const Vector &getBasicDisplSensitivity(int gradNumber);
128  const Vector &getGlobalResistingForceShapeSensitivity(const Vector &q,const Vector &p0,int gradNumber);
129  const Vector &getBasicTrialDispShapeSensitivity(void);
130  bool isShapeSensitivity(void);
131  double getdLdh(void);
132  double getd1overLdh(void);
133 // AddingSensitivity:END //////////////////////////////////
134 
135  CrdTransf2d *getCopy(void) const;
136 
137  int sendSelf(CommParameters &);
138  int recvSelf(const CommParameters &);
139 
140  void Print(std::ostream &s, int flag = 0);
141 
142  // functions used in post-processing only
143  const Vector &getPointGlobalCoordFromLocal(const Vector &) const;
144  const Vector &getPointGlobalDisplFromBasic(double xi, const Vector &basicDisps) const;
145  const Vector &getI(void);
146  const Vector &getJ(void);
147  int getLocalAxes(Vector &xAxis, Vector &yAxis) const;
148 
149  };
150 } // end of XC namespace
151 #endif
int getLocalAxes(Vector &xAxis, Vector &yAxis) const
element position.
Definition: CorotCrdTransf2d.cpp:1027
const Vector & getJ(void)
Returns the ${j}$ unit vector of the local axis expressed in global coordinates for the current geome...
Definition: CorotCrdTransf2d.cpp:1017
Coordinate transformation corrotacional en 3d.
Definition: CorotCrdTransf2d.h:71
int recvData(const CommParameters &)
Receives object members through the channel being passed as parameter.
Definition: CorotCrdTransf2d.cpp:947
Definition: Vector.h:82
const Vector & getGlobalResistingForce(const Vector &basicForce, const Vector &uniformLoad) const
Definition: CorotCrdTransf2d.cpp:500
Mesh node.
Definition: Node.h:99
CrdTransf2d * getCopy(void) const
Virtual constructor.
Definition: CorotCrdTransf2d.cpp:928
int sendSelf(CommParameters &)
Envia el objeto through the channel being passed as parameter.
Definition: CorotCrdTransf2d.cpp:962
int sendData(CommParameters &)
Sends object members through the channel being passed as parameter.
Definition: CorotCrdTransf2d.cpp:932
Definition: Matrix.h:82
int recvSelf(const CommParameters &)
Receives object through the channel being passed as parameter.
Definition: CorotCrdTransf2d.cpp:976
CorotCrdTransf2d(int tag=0)
Constructor.
Definition: CorotCrdTransf2d.cpp:95
Communication parameters between processes.
Definition: CommParameters.h:65
Base class for 2D coordinate transformation.
Definition: CrdTransf2d.h:77
================================================================================
Definition: ContinuaReprComponent.h:34
const Vector & getI(void)
Returns the ${i}$ unit vector of the local axis expressed in global coordinates for the current geome...
Definition: CorotCrdTransf2d.cpp:1007
void Print(std::ostream &s, int flag=0)
Printing.
Definition: CorotCrdTransf2d.cpp:1035