XC Open source finite element analysis program
CrdTransf.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.4 $
48 // $Date: 2005/12/15 00:30:38 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/coordTransformation/CrdTransf.h,v $
50 
51 
52 // File: ~/crdTransf/CrdTransf.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 //
60 // What: "@(#) CrdTransf.h, revA"
61 
62 #ifndef CrdTransf_h
63 #define CrdTransf_h
64 
65 #include "utility/actor/actor/MovableObject.h"
66 #include <utility/tagged/TaggedObject.h>
67 #include <utility/matrix/Vector.h>
68 
69 class TritrizPos3d;
70 
71 namespace XC {
72 class Matrix;
73 class Node;
74 class TransfCooLoader;
75 
76 
78 //
80 //
82 //
87 class CrdTransf: public TaggedObject, public MovableObject
88  {
89  protected:
90  Node *nodeIPtr, *nodeJPtr;
91  mutable double L;
92  Vector nodeIOffset,nodeJOffset;
93  std::vector<double> nodeIInitialDisp;
94  std::vector<double> nodeJInitialDisp;
95  bool initialDispChecked;
96 
97  int set_node_ptrs(Node *nodeIPointer, Node *nodeJPointer);
98 
99  int sendData(CommParameters &cp);
100  int recvData(const CommParameters &cp);
101  virtual void set_rigid_joint_offsetI(const Vector &rigJntOffsetI)= 0;
102  virtual void set_rigid_joint_offsetJ(const Vector &rigJntOffsetJ)= 0;
103  public:
104  CrdTransf(int tag, int classTag, int dim_joint_offset);
105  CrdTransf(void);
106  virtual ~CrdTransf(void);
107 
108  const TransfCooLoader *GetTransfCooLoader(void) const;
110  std::string getName(void) const;
111 
112  virtual int initialize(Node *node1Pointer, Node *node2Pointer) = 0;
113  virtual int update(void) = 0;
114  virtual double getInitialLength(void) const= 0;
115  virtual double getDeformedLength(void) const= 0;
116  double getLength(bool initialGeometry= true) const;
117 
118  virtual int commitState(void) = 0;
119  virtual int revertToLastCommit(void) = 0;
120  virtual int revertToStart(void) = 0;
121 
122  virtual const Vector &getBasicTrialDisp(void) const= 0;
123  virtual const Vector &getBasicIncrDisp(void) const= 0;
124  virtual const Vector &getBasicIncrDeltaDisp(void) const= 0;
125  virtual const Vector &getBasicTrialVel(void) const= 0;
126  virtual const Vector &getBasicTrialAccel(void) const= 0;
127 
128  // AddingSensitivity:BEGIN //////////////////////////////////
129  virtual const Vector &getBasicDisplSensitivity(int gradNumber);
130  virtual const Vector &getGlobalResistingForceShapeSensitivity(const Vector &basicForce, const Vector &uniformLoad);
131  virtual const Vector &getBasicTrialDispShapeSensitivity(void);
132  // AddingSensitivity:END //////////////////////////////////
133 
134  virtual const Vector &getGlobalResistingForce(const Vector &basicForce, const Vector &uniformLoad) const= 0;
135  virtual const Matrix &getGlobalStiffMatrix(const Matrix &basicStiff, const Vector &basicForce) const= 0;
136  virtual const Matrix &getInitialGlobalStiffMatrix(const Matrix &basicStiff) const= 0;
137 
138  virtual const Vector &getI(void) const= 0;
139  virtual const Vector &getJ(void) const= 0;
140  virtual Matrix getLocalAxes(bool) const= 0;
141 
142  virtual const Vector &getPointGlobalCoordFromLocal(const Vector &localCoords) const= 0;
143  virtual const Vector &getPointGlobalCoordFromBasic(const double &xi) const= 0;
144  virtual Vector getPointLocalCoordFromGlobal(const Vector &globalCoords) const= 0;
145  double getPointBasicCoordFromGlobal(const Vector &globalCoords) const;
146  const Matrix &getPointsGlobalCoordFromLocal(const Matrix &localCoords) const;
147  virtual const Matrix &getPointsGlobalCoordFromBasic(const Vector &) const= 0;
148  virtual const Vector &getPointGlobalDisplFromBasic(double xi, const Vector &basicDisps) const= 0;
149 
150  virtual const Vector &getVectorGlobalCoordFromLocal(const Vector &localCoords) const= 0;
151  virtual const Matrix &getVectorGlobalCoordFromLocal(const Matrix &localCoords) const= 0;
152  virtual const Vector &getVectorLocalCoordFromGlobal(const Vector &globalCoords) const= 0;
153 
154  virtual const Matrix &getCooNodos(void) const= 0;
155  virtual const Matrix &getCooPuntos(const size_t &ndiv) const= 0;
156  virtual const Vector &getCooPunto(const double &xrel) const= 0;
157 
158  };
159 } // end of XC namespace
160 
161 #endif
CrdTransf provides the abstraction of a frame coordinate transformation. It is an abstract base class...
Definition: CrdTransf.h:87
int sendData(CommParameters &cp)
Sends object members through the channel being passed as parameter.
Definition: CrdTransf.cpp:189
Node * nodeJPtr
pointers to los nodos extremos of the element.
Definition: CrdTransf.h:90
Definition: Vector.h:82
std::string getName(void) const
Returns the name of the material.
Definition: CrdTransf.cpp:105
std::vector< double > nodeIInitialDisp
Initial displacement for I node.
Definition: CrdTransf.h:93
Mesh node.
Definition: Node.h:99
Object that can move between processes.
Definition: MovableObject.h:91
std::vector< double > nodeJInitialDisp
Initial displacement for J node.
Definition: CrdTransf.h:94
int set_node_ptrs(Node *nodeIPointer, Node *nodeJPointer)
Asigna los pointers to node dorsal y frontal.
Definition: CrdTransf.cpp:113
Definition: Matrix.h:82
int recvData(const CommParameters &cp)
Receives object members through the channel being passed as parameter.
Definition: CrdTransf.cpp:204
Object idenfied by an integer (tag).
Definition: TaggedObject.h:82
double L
element length
Definition: CrdTransf.h:91
virtual ~CrdTransf(void)
Destructor virtual.
Definition: CrdTransf.cpp:83
Communication parameters between processes.
Definition: CommParameters.h:65
================================================================================
Definition: ContinuaReprComponent.h:34
const TransfCooLoader * GetTransfCooLoader(void) const
Returns (if possible) a pointer to the coordinate transformation handler (owner). ...
Definition: CrdTransf.cpp:87
Vector nodeJOffset
rigid joint offsets
Definition: CrdTransf.h:92
Manager for the creation/deletion of coordinate transformations.
Definition: TransfCooLoader.h:49