29 #ifndef ASDShellQ4CorotationalTransformation_h 30 #define ASDShellQ4CorotationalTransformation_h 33 #include "ASDShellQ4Transformation.h" 39 #define USE_POLAR_DECOMP_ALLIGN 71 std::array<QuaternionType, 4> m_QN;
72 std::array<Vector3Type, 4> m_RV;
73 std::array<QuaternionType, 4> m_QN_converged;
74 std::array<Vector3Type, 4> m_RV_converged;
84 virtual bool isLinear()
const 87 virtual void revertToStart(
void);
88 virtual void setDomain(
Domain* domain,
const ID& node_ids);
90 virtual void revertToLastCommit()
92 for (
int i = 0; i < 4; i++)
94 m_RV[i] = m_RV_converged[i];
95 m_QN[i] = m_QN_converged[i];
101 for (
int i = 0; i < 4; i++)
103 m_RV_converged[i] = m_RV[i];
104 m_QN_converged[i] = m_QN[i];
108 virtual void update(
const VectorType& globalDisplacements)
110 for (
int i = 0; i < 4; i++)
113 Vector3Type currentRotVec;
115 currentRotVec(0) = globalDisplacements(index + 3) - m_U0(index + 3);
116 currentRotVec(1) = globalDisplacements(index + 4) - m_U0(index + 4);
117 currentRotVec(2) = globalDisplacements(index + 5) - m_U0(index + 5);
120 Vector3Type incrementalRotation = currentRotVec - m_RV[i];
123 m_RV[i] = currentRotVec;
129 m_QN[i] = incrementalQuaternion * m_QN[i];
139 std::array<Vector3Type, 4> def = {
140 Vector3Type(m_nodes[0]->getCrds()),
141 Vector3Type(m_nodes[1]->getCrds()),
142 Vector3Type(m_nodes[2]->getCrds()),
143 Vector3Type(m_nodes[3]->getCrds())
145 for (
int i = 0; i < 4; i++) {
147 Vector3Type& iP = def[i];
148 iP(0) += globalDisplacements(index) - m_U0(index);
149 iP(1) += globalDisplacements(index + 1) - m_U0(index + 1);
150 iP(2) += globalDisplacements(index + 2) - m_U0(index + 2);
156 #ifndef USE_POLAR_DECOMP_ALLIGN 158 #endif // !USE_POLAR_DECOMP_ALLIGN 160 double aX1 = a.X1();
double aY1 = a.Y1();
161 double bX1 = b.X1();
double bY1 = b.Y1();
162 double aX2 = a.X2();
double aY2 = a.Y2();
163 double bX2 = b.X2();
double bY2 = b.Y2();
164 double aX3 = a.X3();
double aY3 = a.Y3();
165 double bX3 = b.X3();
double bY3 = b.Y3();
166 double aX4 = a.X4();
double aY4 = a.Y4();
167 double bX4 = b.X4();
double bY4 = b.Y4();
173 double C1 = 1.0 / (aX1 * aY2 - aX2 * aY1 - aX1 * aY4 + aX2 * aY3 - aX3 * aY2 + aX4 * aY1 + aX3 * aY4 - aX4 * aY3);
174 double C2 = bY1 / 4.0 + bY2 / 4.0 - bY3 / 4.0 - bY4 / 4.0;
175 double C3 = bY1 / 4.0 - bY2 / 4.0 - bY3 / 4.0 + bY4 / 4.0;
176 double C4 = bX1 / 4.0 + bX2 / 4.0 - bX3 / 4.0 - bX4 / 4.0;
177 double C5 = bX1 / 4.0 - bX2 / 4.0 - bX3 / 4.0 + bX4 / 4.0;
178 double C6 = aX1 + aX2 - aX3 - aX4;
179 double C7 = aX1 - aX2 - aX3 + aX4;
180 double C8 = aY1 + aY2 - aY3 - aY4;
181 double C9 = aY1 - aY2 - aY3 + aY4;
182 double f11 = 2.0 * C1 * C5 * C8 - 2.0 * C1 * C4 * C9;
183 double f12 = 2.0 * C1 * C4 * C7 - 2.0 * C1 * C5 * C6;
184 double f21 = 2.0 * C1 * C3 * C8 - 2.0 * C1 * C2 * C9;
185 double f22 = 2.0 * C1 * C2 * C7 - 2.0 * C1 * C3 * C6;
189 double alpha = std::atan2(f21 - f12, f11 + f22);
196 virtual void calculateLocalDisplacements(
203 const Vector3Type& C = LCS.Center();
205 for (
int i = 0; i < 4; i++)
210 Vector3Type X0 = Vector3Type(m_nodes[i]->getCrds());
214 Vector3Type X = X0 + Vector3Type(globalDisplacements, index);
220 Vector3Type deformationalDisplacements = X - X0;
222 localDisplacements[index] = deformationalDisplacements[0];
223 localDisplacements[index + 1] = deformationalDisplacements[1];
224 localDisplacements[index + 2] = deformationalDisplacements[2];
229 localDisplacements[index + 3],
230 localDisplacements[index + 4],
231 localDisplacements[index + 5]);
235 virtual void transformToGlobal(
251 LCS.ComputeTotalRotationMatrix(T);
262 RotationGradient(LCS, globalDisplacements, G);
263 P.addMatrixProduct(1.0, S, G, -1.0);
270 projectedLocalForces.addMatrixTransposeVector(0.0, P, RHS, 1.0);
275 RHS.addMatrixTransposeVector(0.0, T, projectedLocalForces, 1.0);
292 temp.addMatrixProduct(0.0, LHS, H, 1.0);
293 LHS.addMatrixProduct(0.0, temp, P, 1.0);
294 temp.addMatrixTransposeProduct(0.0, P, LHS, 1.0);
312 FnmT.addMatrixTranspose(0.0, Fnm, 1.0);
314 temp.addMatrixTransposeProduct(0.0, G, FnmT, 1.0);
315 LHS.addMatrixProduct(1.0, temp, P, 1.0);
327 LHS.addMatrixProduct(1.0, Fnm, G, 1.0);
332 temp.addMatrixProduct(0.0, LHS, T, 1.0);
333 LHS.addMatrixTransposeProduct(0.0, T, temp, 1.0);
336 virtual void transformToGlobal(
344 computeGlobalDisplacements(globalDisplacements);
345 calculateLocalDisplacements(LCS, globalDisplacements, localDisplacements);
346 transformToGlobal(LCS, globalDisplacements, localDisplacements, LHS, RHS, LHSrequired);
349 virtual int internalDataSize()
const 359 Vector retval(internalDataSize());
363 for (
int i = 0; i < 24; i++)
364 retval(pos++) = m_U0(i);
369 retval(pos++) = x.w();
370 retval(pos++) = x.x();
371 retval(pos++) = x.y();
372 retval(pos++) = x.z();
375 for (
int i = 0; i < 4; i++)
377 for (
int i = 0; i < 4; i++)
378 lamq(m_QN_converged[i]);
381 auto lamv = [&retval, &pos](
const Vector3Type& x)
383 retval(pos++) = x.x();
384 retval(pos++) = x.y();
385 retval(pos++) = x.z();
388 for (
int i = 0; i < 4; i++)
390 for (
int i = 0; i < 4; i++)
391 lamv(m_RV_converged[i]);
400 for (
int i = 0; i < 24; i++)
406 const double ww= v(pos++);
407 const double xx= v(pos++);
408 const double yy= v(pos++);
409 const double zz= v(pos++);
413 for(
int i = 0; i < 4; i++)
415 for(
int i = 0; i < 4; i++)
416 lamq(m_QN_converged[i]);
419 auto lamv = [&v, &pos](Vector3Type& x)
421 const double xx= v(pos++);
422 const double yy= v(pos++);
423 const double zz= v(pos++);
424 x = Vector3Type(xx, yy, zz);
427 for (
int i = 0; i < 4; i++)
429 for (
int i = 0; i < 4; i++)
430 lamv(m_RV_converged[i]);
444 inline void RotationGradient(
451 #ifdef USE_POLAR_DECOMP_ALLIGN 500 #else // !USE_POLAR_DECOMP_ALLIGN 502 const auto& P1 = LCS.P1();
503 const auto& P2 = LCS.P2();
504 const auto& P3 = LCS.P3();
505 const auto& P4 = LCS.P4();
507 double Ap = 2.0 * LCS.Area();
510 Vector3Type D12(P2 - P1);
511 Vector3Type D24(P4 - P2);
512 Vector3Type D13(P3 - P1);
525 double l12 = std::sqrt(D12(0) * D12(0) + D12(1) * D12(1));
531 G(2, 1) = -1.0 / l12;
549 #endif // USE_POLAR_DECOMP_ALLIGN 555 #endif // !ASDShellQ4CorotationalTransformation_h Float vector abstraction.
Definition: Vector.h:94
void Zero(void)
Zero's out the Matrix.
Definition: Matrix.cpp:226
static void Compute_Pt(size_t num_nodes, MatrixType &P)
Computes the Translational Projector Matrix.
Definition: ASDEICR.h:273
Vector of integers.
Definition: ID.h:95
static void Compute_H(const VectorType &displacements, MatrixType &H)
Computes the Axial Vector Jacobian.
Definition: ASDEICR.h:340
This class represent the local coordinate system of any element whose geometry is a 4-node Quadrilate...
Definition: ASDShellQ4LocalCoordinateSystem.h:41
void toRotationVector(T &rx, T &ry, T &rz) const
Extracts the Rotation Vector from this ASDQuaternion.
Definition: ASDMath.h:495
ASDQuaternion conjugate() const
Returns the Conjugate of this ASDQuaternion, which represents the opposite rotation.
Definition: ASDMath.h:456
static void Spin_AtRow(const TVec &V, TMat &S, size_t row_index)
Computes the Spin of the input vector V, and saves the result into the output matrix S...
Definition: ASDEICR.h:98
void rotateVector(const TVector3_A &a, TVector3_B &b) const
Rotates a vector using this quaternion.
Definition: ASDMath.h:556
static ASDQuaternion FromRotationVector(T rx, T ry, T rz)
Returns a ASDQuaternion from a rotation vector.
Definition: ASDMath.h:653
static void Compute_S(const NodeContainerType &nodes, MatrixType &S)
Computes the Spin Lever Matrix.
Definition: ASDEICR.h:314
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
static ASDQuaternion FromRotationMatrix(const TMatrix3x3 &m)
Returns a ASDQuaternion from a Rotation Matrix.
Definition: ASDMath.h:704
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
ASDQuaternion A simple class that implements the main features of quaternion algebra.
Definition: ASDMath.h:328