63 #ifndef CorotCrdTransf3d_h 64 #define CorotCrdTransf3d_h 66 #include "CrdTransf3d.h" 67 #include <utility/matrix/Vector.h> 68 #include <utility/matrix/Matrix.h> 104 static Matrix Lr2, Lr3, A;
106 inline int computeElemtLengthAndOrient(
void)
const 108 std::cerr <<
"CorotCrdTransf3d::computeElemtLengthAndOrient; not implemented." 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;
120 const Matrix &getRotationMatrixFromQuaternion(
const Vector &q)
const;
121 const Matrix &getRotMatrixFromTangScaledPseudoVector(
const Vector &w)
const;
122 const Matrix &getSkewSymMatrix(
const Vector &theta)
const;
128 Vector &basic_to_local_element_force(
const Vector &p0)
const;
129 const Vector &local_to_global_element_force(
const Vector &)
const;
133 virtual int computeLocalAxis(
void)
const;
138 int initialize(
Node *nodeIPointer,
Node *nodeJPointer);
140 double getInitialLength(
void)
const;
141 double getDeformedLength(
void)
const;
145 int commitState(
void);
146 int revertToLastCommit(
void);
147 int revertToStart(
void);
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;
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;
164 void Print(std::ostream &s,
int flag = 0)
const;
167 const Matrix &getGlobalMatrixFromLocal(
const Matrix &local);
170 const Vector &getPointGlobalCoordFromLocal(
const Vector &)
const;
171 const Vector &getPointGlobalDisplFromBasic(
double xi,
const Vector &basicDisps)
const;
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