quaternionI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 :
33  w_(0),
34  v_(Zero)
35 {}
36 
37 
38 inline Foam::quaternion::quaternion(const scalar w, const vector& v)
39 :
40  w_(w),
41  v_(v)
42 {}
43 
44 
45 inline Foam::quaternion::quaternion(const vector& d, const scalar theta)
46 :
47  w_(cos(0.5*theta)),
48  v_((sin(0.5*theta)/mag(d))*d)
49 {}
50 
51 
53 (
54  const vector& d,
55  const scalar cosTheta,
56  const bool normalized
57 )
58 {
59  scalar cosHalfTheta2 = 0.5*(cosTheta + 1);
60  w_ = sqrt(cosHalfTheta2);
61 
62  if (normalized)
63  {
64  v_ = sqrt(1 - cosHalfTheta2)*d;
65  }
66  else
67  {
68  v_ = (sqrt(1 - cosHalfTheta2)/mag(d))*d;
69  }
70 }
71 
72 
73 inline Foam::quaternion::quaternion(const scalar w)
74 :
75  w_(w),
76  v_(Zero)
77 {}
78 
79 
81 :
82  w_(0),
83  v_(v)
84 {}
85 
86 
88 {
89  return quaternion(sqrt(1 - magSqr(v)), v);
90 }
91 
92 
94 (
95  const eulerOrder order,
96  const vector& angles
97 )
98 {
99  switch (order)
100  {
101  case ZYX:
102  {
103  operator=(quaternion(vector(0, 0, 1), angles.x()));
104  operator*=(quaternion(vector(0, 1, 0), angles.y()));
105  operator*=(quaternion(vector(1, 0, 0), angles.z()));
106  break;
107  }
108 
109  case ZYZ:
110  {
111  operator=(quaternion(vector(0, 0, 1), angles.x()));
112  operator*=(quaternion(vector(0, 1, 0), angles.y()));
113  operator*=(quaternion(vector(0, 0, 1), angles.z()));
114  break;
115  }
116 
117  case ZXY:
118  {
119  operator=(quaternion(vector(0, 0, 1), angles.x()));
120  operator*=(quaternion(vector(1, 0, 0), angles.y()));
121  operator*=(quaternion(vector(0, 1, 0), angles.z()));
122  break;
123  }
124 
125  case ZXZ:
126  {
127  operator=(quaternion(vector(0, 0, 1), angles.x()));
128  operator*=(quaternion(vector(1, 0, 0), angles.y()));
129  operator*=(quaternion(vector(0, 0, 1), angles.z()));
130  break;
131  }
132 
133  case YXZ:
134  {
135  operator=(quaternion(vector(0, 1, 0), angles.x()));
136  operator*=(quaternion(vector(1, 0, 0), angles.y()));
137  operator*=(quaternion(vector(0, 0, 1), angles.z()));
138  break;
139  }
140 
141  case YXY:
142  {
143  operator=(quaternion(vector(0, 1, 0), angles.x()));
144  operator*=(quaternion(vector(1, 0, 0), angles.y()));
145  operator*=(quaternion(vector(0, 1, 0), angles.z()));
146  break;
147  }
148 
149  case YZX:
150  {
151  operator=(quaternion(vector(0, 1, 0), angles.x()));
152  operator*=(quaternion(vector(0, 0, 1), angles.y()));
153  operator*=(quaternion(vector(1, 0, 0), angles.z()));
154  break;
155  }
156 
157  case YZY:
158  {
159  operator=(quaternion(vector(0, 1, 0), angles.x()));
160  operator*=(quaternion(vector(0, 0, 1), angles.y()));
161  operator*=(quaternion(vector(0, 1, 0), angles.z()));
162  break;
163  }
164 
165  case XYZ:
166  {
167  operator=(quaternion(vector(1, 0, 0), angles.x()));
168  operator*=(quaternion(vector(0, 1, 0), angles.y()));
169  operator*=(quaternion(vector(0, 0, 1), angles.z()));
170  break;
171  }
172 
173  case XYX:
174  {
175  operator=(quaternion(vector(1, 0, 0), angles.x()));
176  operator*=(quaternion(vector(0, 1, 0), angles.y()));
177  operator*=(quaternion(vector(1, 0, 0), angles.z()));
178  break;
179  }
180 
181  case XZY:
182  {
183  operator=(quaternion(vector(1, 0, 0), angles.x()));
184  operator*=(quaternion(vector(0, 0, 1), angles.y()));
185  operator*=(quaternion(vector(0, 1, 0), angles.z()));
186  break;
187  }
188 
189  case XZX:
190  {
191  operator=(quaternion(vector(1, 0, 0), angles.x()));
192  operator*=(quaternion(vector(0, 0, 1), angles.y()));
193  operator*=(quaternion(vector(1, 0, 0), angles.z()));
194  break;
195  }
196 
197  default:
199  << "Unknown euler rotation order "
200  << int(order) << abort(FatalError);
201  break;
202  }
203 }
204 
205 
207 (
208  const tensor& rotationTensor
209 )
210 {
211  scalar trace =
213  + rotationTensor.yy()
214  + rotationTensor.zz();
215 
216  if (trace > 0)
217  {
218  scalar s = 0.5/Foam::sqrt(trace + 1.0);
219 
220  w_ = 0.25/s;
221  v_[0] = (rotationTensor.zy() - rotationTensor.yz())*s;
222  v_[1] = (rotationTensor.xz() - rotationTensor.zx())*s;
223  v_[2] = (rotationTensor.yx() - rotationTensor.xy())*s;
224  }
225  else
226  {
227  if
228  (
231  )
232  {
233  const scalar s = 2.0*Foam::sqrt
234  (
235  1.0
236  + rotationTensor.xx()
237  - rotationTensor.yy()
238  - rotationTensor.zz()
239  );
240 
241  w_ = (rotationTensor.zy() - rotationTensor.yz())/s;
242  v_[0] = 0.25*s;
243  v_[1] = (rotationTensor.xy() + rotationTensor.yx())/s;
244  v_[2] = (rotationTensor.xz() + rotationTensor.zx())/s;
245  }
246  else if
247  (
249  )
250  {
251  const scalar s = 2.0*Foam::sqrt
252  (
253  1.0
254  + rotationTensor.yy()
255  - rotationTensor.xx()
256  - rotationTensor.zz()
257  );
258 
259  w_ = (rotationTensor.xz() - rotationTensor.zx())/s;
260  v_[0] = (rotationTensor.xy() + rotationTensor.yx())/s;
261  v_[1] = 0.25*s;
262  v_[2] = (rotationTensor.yz() + rotationTensor.zy())/s;
263  }
264  else
265  {
266  const scalar s = 2.0*Foam::sqrt
267  (
268  1.0
269  + rotationTensor.zz()
270  - rotationTensor.xx()
271  - rotationTensor.yy()
272  );
273 
274  w_ = (rotationTensor.yx() - rotationTensor.xy())/s;
275  v_[0] = (rotationTensor.xz() + rotationTensor.zx())/s;
276  v_[1] = (rotationTensor.yz() + rotationTensor.zy())/s;
277  v_[2] = 0.25*s;
278  }
279  }
280 }
281 
282 
283 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
284 
285 inline Foam::scalar Foam::quaternion::w() const
286 {
287  return w_;
288 }
289 
290 
291 inline const Foam::vector& Foam::quaternion::v() const
292 {
293  return v_;
294 }
295 
296 
297 inline Foam::scalar& Foam::quaternion::w()
298 {
299  return w_;
300 }
301 
302 
304 {
305  return v_;
306 }
307 
308 
310 {
311  return *this/mag(*this);
312 }
313 
314 
316 {
317  operator/=(mag(*this));
318 }
319 
320 
321 inline Foam::quaternion Foam::quaternion::mulq0v(const vector& u) const
322 {
323  return quaternion(-(v() & u), w()*u + (v() ^ u));
324 }
325 
326 
328 {
329  return (mulq0v(u)*conjugate(*this)).v();
330 }
331 
332 
334 {
335  return (conjugate(*this).mulq0v(u)*(*this)).v();
336 }
337 
338 
340 {
341  return Foam::normalize((*this)*q);
342 }
343 
344 
346 (
347  const quaternion& q
348 ) const
349 {
350  return Foam::normalize(conjugate(*this)*q);
351 }
352 
353 
355 {
356  const scalar w2 = sqr(w());
357  const scalar x2 = sqr(v().x());
358  const scalar y2 = sqr(v().y());
359  const scalar z2 = sqr(v().z());
360 
361  const scalar txy = 2*v().x()*v().y();
362  const scalar twz = 2*w()*v().z();
363  const scalar txz = 2*v().x()*v().z();
364  const scalar twy = 2*w()*v().y();
365  const scalar tyz = 2*v().y()*v().z();
366  const scalar twx = 2*w()*v().x();
367 
368  return tensor
369  (
370  w2 + x2 - y2 - z2, txy - twz, txz + twy,
371  txy + twz, w2 - x2 + y2 - z2, tyz - twx,
372  txz - twy, tyz + twx, w2 - x2 - y2 + z2
373  );
374 }
375 
376 
377 inline Foam::vector Foam::quaternion::twoAxes
378 (
379  const scalar t11,
380  const scalar t12,
381  const scalar c2,
382  const scalar t31,
383  const scalar t32
384 )
385 {
386  return vector(atan2(t11, t12), acos(c2), atan2(t31, t32));
387 }
388 
389 
390 inline Foam::vector Foam::quaternion::threeAxes
391 (
392  const scalar t11,
393  const scalar t12,
394  const scalar s2,
395  const scalar t31,
396  const scalar t32
397 )
398 {
399  return vector(atan2(t11, t12), asin(s2), atan2(t31, t32));
400 }
401 
402 
404 (
405  const eulerOrder order
406 ) const
407 {
408  const scalar w2 = sqr(w());
409  const scalar x2 = sqr(v().x());
410  const scalar y2 = sqr(v().y());
411  const scalar z2 = sqr(v().z());
412 
413  switch (order)
414  {
415  case ZYX:
416  {
417  return threeAxes
418  (
419  2*(v().x()*v().y() + w()*v().z()),
420  w2 + x2 - y2 - z2,
421  2*(w()*v().y() - v().x()*v().z()),
422  2*(v().y()*v().z() + w()*v().x()),
423  w2 - x2 - y2 + z2
424  );
425  break;
426  }
427 
428  case ZYZ:
429  {
430  return twoAxes
431  (
432  2*(v().y()*v().z() - w()*v().x()),
433  2*(v().x()*v().z() + w()*v().y()),
434  w2 - x2 - y2 + z2,
435  2*(v().y()*v().z() + w()*v().x()),
436  2*(w()*v().y() - v().x()*v().z())
437  );
438  break;
439  }
440 
441  case ZXY:
442  {
443  return threeAxes
444  (
445  2*(w()*v().z() - v().x()*v().y()),
446  w2 - x2 + y2 - z2,
447  2*(v().y()*v().z() + w()*v().x()),
448  2*(w()*v().y() - v().x()*v().z()),
449  w2 - x2 - y2 + z2
450  );
451  break;
452  }
453 
454  case ZXZ:
455  {
456  return twoAxes
457  (
458  2*(v().x()*v().z() + w()*v().y()),
459  2*(w()*v().x() - v().y()*v().z()),
460  w2 - x2 - y2 + z2,
461  2*(v().x()*v().z() - w()*v().y()),
462  2*(v().y()*v().z() + w()*v().x())
463  );
464  break;
465  }
466 
467  case YXZ:
468  {
469  return threeAxes
470  (
471  2*(v().x()*v().z() + w()*v().y()),
472  w2 - x2 - y2 + z2,
473  2*(w()*v().x() - v().y()*v().z()),
474  2*(v().x()*v().y() + w()*v().z()),
475  w2 - x2 + y2 - z2
476  );
477  break;
478  }
479 
480  case YXY:
481  {
482  return twoAxes
483  (
484  2*(v().x()*v().y() - w()*v().z()),
485  2*(v().y()*v().z() + w()*v().x()),
486  w2 - x2 + y2 - z2,
487  2*(v().x()*v().y() + w()*v().z()),
488  2*(w()*v().x() - v().y()*v().z())
489  );
490  break;
491  }
492 
493  case YZX:
494  {
495  return threeAxes
496  (
497  2*(w()*v().y() - v().x()*v().z()),
498  w2 + x2 - y2 - z2,
499  2*(v().x()*v().y() + w()*v().z()),
500  2*(w()*v().x() - v().y()*v().z()),
501  w2 - x2 + y2 - z2
502  );
503  break;
504  }
505 
506  case YZY:
507  {
508  return twoAxes
509  (
510  2*(v().y()*v().z() + w()*v().x()),
511  2*(w()*v().z() - v().x()*v().y()),
512  w2 - x2 + y2 - z2,
513  2*(v().y()*v().z() - w()*v().x()),
514  2*(v().x()*v().y() + w()*v().z())
515  );
516  break;
517  }
518 
519  case XYZ:
520  {
521  return threeAxes
522  (
523  2*(w()*v().x() - v().y()*v().z()),
524  w2 - x2 - y2 + z2,
525  2*(v().x()*v().z() + w()*v().y()),
526  2*(w()*v().z() - v().x()*v().y()),
527  w2 + x2 - y2 - z2
528  );
529  break;
530  }
531 
532  case XYX:
533  {
534  return twoAxes
535  (
536  2*(v().x()*v().y() + w()*v().z()),
537  2*(w()*v().y() - v().x()*v().z()),
538  w2 + x2 - y2 - z2,
539  2*(v().x()*v().y() - w()*v().z()),
540  2*(v().x()*v().z() + w()*v().y())
541  );
542  break;
543  }
544 
545  case XZY:
546  {
547  return threeAxes
548  (
549  2*(v().y()*v().z() + w()*v().x()),
550  w2 - x2 + y2 - z2,
551  2*(w()*v().z() - v().x()*v().y()),
552  2*(v().x()*v().z() + w()*v().y()),
553  w2 + x2 - y2 - z2
554  );
555  break;
556  }
557 
558  case XZX:
559  {
560  return twoAxes
561  (
562  2*(v().x()*v().z() - w()*v().y()),
563  2*(v().x()*v().y() + w()*v().z()),
564  w2 + x2 - y2 - z2,
565  2*(v().x()*v().z() + w()*v().y()),
566  2*(w()*v().z() - v().x()*v().y())
567  );
568  break;
569  }
570 
571  default:
573  << "Unknown euler rotation order "
574  << int(order) << abort(FatalError);
575  break;
576  }
577 
578  return Zero;
579 }
580 
581 
582 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
583 
585 {
586  w_ += q.w_;
587  v_ += q.v_;
588 }
589 
591 {
592  w_ -= q.w_;
593  v_ -= q.v_;
594 }
595 
597 {
598  scalar w0 = w();
599  w() = w()*q.w() - (v() & q.v());
600  v() = w0*q.v() + q.w()*v() + (v() ^ q.v());
601 }
602 
604 {
605  return operator*=(inv(q));
606 }
607 
608 
609 inline void Foam::quaternion::operator=(const scalar s)
610 {
611  w_ = s;
612 }
613 
614 
615 inline void Foam::quaternion::operator=(const vector& v)
616 {
617  v_ = v;
618 }
619 
620 
621 inline void Foam::quaternion::operator*=(const scalar s)
622 {
623  w_ *= s;
624  v_ *= s;
625 }
626 
627 inline void Foam::quaternion::operator/=(const scalar s)
628 {
629  w_ /= s;
630  v_ /= s;
631 }
632 
633 
634 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
635 
636 inline Foam::scalar Foam::magSqr(const quaternion& q)
637 {
638  return magSqr(q.w()) + magSqr(q.v());
639 }
640 
641 
642 inline Foam::scalar Foam::mag(const quaternion& q)
643 {
644  return sqrt(magSqr(q));
645 }
646 
647 
649 {
650  return quaternion(q.w(), -q.v());
651 }
652 
653 
655 {
656  scalar magSqrq = magSqr(q);
657  return quaternion(q.w()/magSqrq, -q.v()/magSqrq);
658 }
659 
660 
662 {
663  return q/mag(q);
664 }
665 
666 
667 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
668 
669 inline bool Foam::operator==(const quaternion& q1, const quaternion& q2)
670 {
671  return (equal(q1.w(), q2.w()) && equal(q1.v(), q2.v()));
672 }
673 
674 
675 inline bool Foam::operator!=(const quaternion& q1, const quaternion& q2)
676 {
677  return !operator==(q1, q2);
678 }
679 
680 
681 inline Foam::quaternion Foam::operator+
682 (
683  const quaternion& q1,
684  const quaternion& q2
685 )
686 {
687  return quaternion(q1.w() + q2.w(), q1.v() + q2.v());
688 }
689 
690 
692 {
693  return quaternion(-q.w(), -q.v());
694 }
695 
696 
697 inline Foam::quaternion Foam::operator-
698 (
699  const quaternion& q1,
700  const quaternion& q2
701 )
702 {
703  return quaternion(q1.w() - q2.w(), q1.v() - q2.v());
704 }
705 
706 
707 inline Foam::scalar Foam::operator&(const quaternion& q1, const quaternion& q2)
708 {
709  return q1.w()*q2.w() + (q1.v() & q2.v());
710 }
711 
712 
713 inline Foam::quaternion Foam::operator*
714 (
715  const quaternion& q1,
716  const quaternion& q2
717 )
718 {
719  return quaternion
720  (
721  q1.w()*q2.w() - (q1.v() & q2.v()),
722  q1.w()*q2.v() + q2.w()*q1.v() + (q1.v() ^ q2.v())
723  );
724 }
725 
726 
727 inline Foam::quaternion Foam::operator/
728 (
729  const quaternion& q1,
730  const quaternion& q2
731 )
732 {
733  return q1*inv(q2);
734 }
735 
736 
737 inline Foam::quaternion Foam::operator*(const scalar s, const quaternion& q)
738 {
739  return quaternion(s*q.w(), s*q.v());
740 }
741 
742 
743 inline Foam::quaternion Foam::operator*(const quaternion& q, const scalar s)
744 {
745  return quaternion(s*q.w(), s*q.v());
746 }
747 
748 
749 inline Foam::quaternion Foam::operator/(const quaternion& q, const scalar s)
750 {
751  return quaternion(q.w()/s, q.v()/s);
752 }
753 
754 
755 // ************************************************************************* //
Foam::Tensor< scalar >
Foam::quaternion::normalize
void normalize()
Normalize the quaternion by its magnitude.
Definition: quaternionI.H:315
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Foam::quaternion::operator+=
void operator+=(const quaternion &q)
Definition: quaternionI.H:584
Foam::quaternion::v
const vector & v() const
Vector part of the quaternion ( = axis of rotation)
Definition: quaternionI.H:291
Foam::Tensor::zx
const Cmpt & zx() const
Definition: TensorI.H:194
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::operator&
tmp< GeometricField< Type, fvPatchField, volMesh > > operator&(const fvMatrix< Type > &, const DimensionedField< Type, volMesh > &)
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::quaternion::operator=
quaternion & operator=(const quaternion &)=default
Copy assignment.
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::quaternion::operator/=
void operator/=(const quaternion &q)
Definition: quaternionI.H:603
Foam::atan2
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Definition: dimensionedScalar.C:312
Foam::operator-
tmp< faMatrix< Type > > operator-(const faMatrix< Type > &)
Foam::quaternion::normalized
quaternion normalized() const
Return the quaternion normalized by its magnitude.
Definition: quaternionI.H:309
Foam::quaternion::operator-=
void operator-=(const quaternion &q)
Definition: quaternionI.H:590
Foam::Tensor::yx
const Cmpt & yx() const
Definition: TensorI.H:173
Foam::quaternion::eulerAngles
vector eulerAngles(const eulerOrder order) const
Definition: quaternionI.H:404
Foam::quaternion::w
scalar w() const
Scalar part of the quaternion ( = cos(theta/2) for rotation)
Definition: quaternionI.H:285
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
Foam::quaternion::quaternion
quaternion()=default
Default construct.
Foam::Tensor::xz
const Cmpt & xz() const
Definition: TensorI.H:166
Foam::quaternion
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:56
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::Tensor::yz
const Cmpt & yz() const
Definition: TensorI.H:187
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
Foam::Tensor::zy
const Cmpt & zy() const
Definition: TensorI.H:201
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::Tensor::yy
const Cmpt & yy() const
Definition: TensorI.H:180
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::Tensor::zz
const Cmpt & zz() const
Definition: TensorI.H:208
Foam::quaternion::eulerOrder
eulerOrder
Euler-angle rotation order.
Definition: quaternion.H:102
Foam::FatalError
error FatalError
w2
#define w2
Definition: blockCreate.C:35
Foam::Tensor::xy
const Cmpt & xy() const
Definition: TensorI.H:159
Foam::quaternion::unit
static quaternion unit(const vector &v)
Definition: quaternionI.H:87
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::constant::physicoChemical::c2
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::operator/
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
Definition: dimensionedScalar.C:68
Foam::Tensor::xx
const Cmpt & xx() const
Definition: TensorI.H:152
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::quaternion::operator*=
void operator*=(const quaternion &q)
Definition: quaternionI.H:596
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::operator*
tmp< faMatrix< Type > > operator*(const areaScalarField &, const faMatrix< Type > &)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::quaternion::invTransform
vector invTransform(const vector &v) const
Rotate the given vector anti-clockwise.
Definition: quaternionI.H:333
Foam::quaternion::R
tensor R() const
The rotation tensor corresponding to the quaternion.
Definition: quaternionI.H:354
w0
#define w0
Definition: blockCreate.C:33
Foam::conjugate
quaternion conjugate(const quaternion &q)
Return the conjugate of the given quaternion.
Definition: quaternionI.H:648
Foam::rotationTensor
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:51
Foam::equal
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:46
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:61
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267
Foam::normalize
quaternion normalize(const quaternion &q)
Return the normalized (unit) quaternion of the given quaternion.
Definition: quaternionI.H:661
Foam::quaternion::transform
vector transform(const vector &v) const
Rotate the given vector.
Definition: quaternionI.H:327
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::zero
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:62
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265