xc
CrdTransf3d.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.3 $
48 // $Date: 2005/12/15 00:30:38 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/CrdTransf3d.h,v $
50 
51 
52 // File: ~/CrdTransf/CrdTransf3d.h
53 //
54 // Written: Remo Magalhaes de Souza (rmsouza@ce.berkeley.edu)
55 // Created: 04/2000
56 // Revision: A
57 //
58 // Description: This file contains the class definition for
59 // CrdTransf3d. CrdTransf3d provides the abstraction of spatial
60 // coordinate transformation for a 3d frame.
61 
62 //
63 // What: "@(#) CrdTransf3d.h, revA"
64 
65 #ifndef CrdTransf3d_h
66 #define CrdTransf3d_h
67 
68 #include "CrdTransf.h"
69 #include "utility/matrix/Matrix.h"
70 
71 class Pos3d;
72 class Ref3d3d;
73 
74 namespace XC {
75 class Vector;
76 
78 //
80 class CrdTransf3d: public CrdTransf
81  {
82  protected:
83  mutable Matrix R;
84 
85  void set_rigid_joint_offsetI(const Vector &rigJntOffsetI);
86  void set_rigid_joint_offsetJ(const Vector &rigJntOffsetJ);
87  static void inic_ug(const Vector &d1,const Vector &d2,double *ug);
88  void modif_ug_init_disp(double *ug) const;
89  void global_to_local(const double *ug,double *ul) const;
90  void calc_Wu(const double *ug,double *ul,double *Wu) const;
91  const Vector &calc_ub(const double *ul,Vector &) const;
92 
93  static Vector vectorI;
94  static Vector vectorJ;
95  static Vector vectorK;
96  static Vector vectorCoo;
97  virtual int computeElemtLengthAndOrient(void) const= 0;
98  virtual int computeLocalAxis(void) const= 0;
99 
100  int sendData(CommParameters &cp);
101  int recvData(const CommParameters &cp);
102  public:
103  CrdTransf3d(int tag, int classTag);
104  CrdTransf3d(int tag, int classTag, const Vector &vecInLocXZPlane);
105 
106  virtual CrdTransf3d *getCopy(void) const= 0;
107  int initialize(Node *node1Pointer, Node *node2Pointer);
108  virtual void set_xz_vector(const Vector &vecInLocXZPlane);
109  Vector get_xz_vector(void) const;
110  const Vector &getI(void) const;
111  const Vector &getJ(void) const;
112  const Vector &getK(void) const;
113  Matrix getLocalAxes(bool) const;
114  int getLocalAxes(Vector &xAxis, Vector &yAxis, Vector &zAxis) const;
115  Pos3d getPosNodeI(void) const;
116  Pos3d getPosNodeJ(void) const;
117  Ref3d3d getLocalReference(void) const;
118  Vector getPointLocalCoordFromGlobal(const Vector &xg) const;
119  const Vector &getPointGlobalCoordFromBasic(const double &xi) const;
120  const Matrix &getPointsGlobalCoordFromBasic(const Vector &) const;
121  const Vector &getVectorGlobalCoordFromLocal(const Vector &localCoords) const;
122  const Matrix &getVectorGlobalCoordFromLocal(const Matrix &localCoords) const;
123  const Vector &getVectorLocalCoordFromGlobal(const Vector &globalCoords) const;
124 
125  const Matrix &getCooNodes(void) const;
126  const Matrix &getCooPoints(const size_t &ndiv) const;
127  const Vector &getCooPoint(const double &xrel) const;
128 
129  void gira(const double &);
130  };
131 } // end of XC namespace
132 
133 #endif
const Vector & getCooPoint(const double &xrel) const
Return the point that correspond to the relative coordinate 0<=xrel<=1.
Definition: CrdTransf3d.cpp:478
const Vector & getI(void) const
Returns the ${i}$ unit vector of the local axis expressed in global coordinates for the current geome...
Definition: CrdTransf3d.cpp:286
void modif_ug_init_disp(double *ug) const
brief Modifical el vector de displacement globales of the nodes de acuerdo con los displacements inic...
Definition: CrdTransf3d.cpp:112
Float vector abstraction.
Definition: Vector.h:93
void set_rigid_joint_offsetJ(const Vector &rigJntOffsetJ)
check rigid joint offset for node J
Definition: CrdTransf3d.cpp:223
void gira(const double &)
Counterclockwise rotate the coordinate transformation by the specified angle.
Definition: CrdTransf3d.cpp:493
const Vector & calc_ub(const double *ul, Vector &) const
Sean dx1,dy1,dz1,gx1,gy1,gz1 los displacements y giros of the node dorsal y dx2,dy2,dz2,gx2,gy2,gz2 los of the node frontal expressed in local coordinates.
Definition: CrdTransf3d.cpp:175
Vector getPointLocalCoordFromGlobal(const Vector &xg) const
Returns the local coordinates del point a partir of the globales.
Definition: CrdTransf3d.cpp:366
CrdTransf provides the abstraction of a frame coordinate transformation.
Definition: CrdTransf.h:87
const Matrix & getPointsGlobalCoordFromBasic(const Vector &) const
Returns the points expressed in global coordinates.
Definition: CrdTransf3d.cpp:386
Pos3d getPosNodeI(void) const
Returns the position of node I.
Definition: CrdTransf3d.cpp:338
const Matrix & getCooNodes(void) const
Returns the coordinates of the nodes.
Definition: CrdTransf3d.cpp:442
Vector get_xz_vector(void) const
Returns the vector that defines the local XZ plane.
Definition: CrdTransf3d.cpp:91
Matrix R
Transformation matrix.
Definition: CrdTransf3d.h:83
Ref3d3d getLocalReference(void) const
Returns the local reference system.
Definition: CrdTransf3d.cpp:358
const Vector & getK(void) const
Returns the ${k}$ unit vector of the local axis expressed in global coordinates for the current geome...
Definition: CrdTransf3d.cpp:302
Matrix getLocalAxes(bool) const
Returs a matrix with the axes of the element as matrix rows [[x1,y1,z1],[x2,y2,z2],...·].
Definition: CrdTransf3d.cpp:323
static void inic_ug(const Vector &d1, const Vector &d2, double *ug)
brief Rellena el vector de displacement globales of the nodes.
Definition: CrdTransf3d.cpp:101
const Matrix & getCooPoints(const size_t &ndiv) const
Return points distributed between the nodes as a matrix with the coordinates as rows.
Definition: CrdTransf3d.cpp:459
const Vector & getVectorLocalCoordFromGlobal(const Vector &globalCoords) const
Returns the vector expresado en local coordinates.
Definition: CrdTransf3d.cpp:432
const Vector & getVectorGlobalCoordFromLocal(const Vector &localCoords) const
Returns the vector expressed in global coordinates.
Definition: CrdTransf3d.cpp:402
void set_rigid_joint_offsetI(const Vector &rigJntOffsetI)
check rigid joint offset for node I
Definition: CrdTransf3d.cpp:205
int sendData(CommParameters &cp)
Sends object members through the channel being passed as parameter.
Definition: CrdTransf3d.cpp:497
CrdTransf3d(int tag, int classTag)
Default constructor.
Definition: CrdTransf3d.cpp:191
virtual void set_xz_vector(const Vector &vecInLocXZPlane)
Set the vector that defines the local XZ plane.
Definition: CrdTransf3d.cpp:80
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
const Vector & getJ(void) const
Returns the ${j}$ unit vector of the local axis expressed in global coordinates for the current geome...
Definition: CrdTransf3d.cpp:294
int recvData(const CommParameters &cp)
Receives object members through the channel being passed as parameter.
Definition: CrdTransf3d.cpp:507
const Vector & getPointGlobalCoordFromBasic(const double &xi) const
Returns the point expresado en global coordinates.
Definition: CrdTransf3d.cpp:376
Communication parameters between processes.
Definition: CommParameters.h:65
Matrix of floats.
Definition: Matrix.h:108
Mesh node.
Definition: Node.h:110
virtual CrdTransf3d * getCopy(void) const =0
Virtual constructor.
Pos3d getPosNodeJ(void) const
Return the position of node J.
Definition: CrdTransf3d.cpp:348
void global_to_local(const double *ug, double *ul) const
Returns node displacements expressed in local coordinates.
Definition: CrdTransf3d.cpp:128
Base class for 3D coordinate transformation.
Definition: CrdTransf3d.h:80