All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
vectorTools.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::vectorTools
26 
27 Description
28  Functions for analysing the relationships between vectors
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef vectorTools_H
33 #define vectorTools_H
34 
35 #include "vector.H"
36 #include "unitConversion.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 /*---------------------------------------------------------------------------*\
44  Namespace vectorTools Declaration
45 \*---------------------------------------------------------------------------*/
46 
47 //- Collection of functions for testing relationships between two vectors.
48 namespace vectorTools
49 {
50  //- Test if a and b are parallel: a^b = 0
51  // Uses the cross product, so the tolerance is proportional to
52  // the sine of the angle between a and b in radians
53  template<class T>
55  (
56  const Vector<T>& a,
57  const Vector<T>& b,
58  const T& tolerance = small
59  )
60  {
61  return (mag(a ^ b) < tolerance) ? true : false;
62 // return ( mag( mag(a & b)/(mag(a)*mag(b)) - 1.0 ) < tolerance )
63 // ? true
64 // : false;
65  }
66 
67  //- Test if a and b are orthogonal: a.b = 0
68  // Uses the dot product, so the tolerance is proportional to
69  // the cosine of the angle between a and b in radians
70  template<class T>
72  (
73  const Vector<T>& a,
74  const Vector<T>& b,
75  const T& tolerance = small
76  )
77  {
78  return (mag(a & b) < tolerance) ? true : false;
79  }
80 
81  //- Test if angle between a and b is acute: a.b > 0
82  template<class T>
83  bool areAcute
84  (
85  const Vector<T>& a,
86  const Vector<T>& b
87  )
88  {
89  return ((a & b) > 0) ? true : false;
90  }
91 
92  //- Test if angle between a and b is obtuse: a.b < 0
93  template<class T>
94  bool areObtuse
95  (
96  const Vector<T>& a,
97  const Vector<T>& b
98  )
99  {
100  return ((a & b) < 0) ? true : false;
101  }
102 
103  //- Calculate angle between a and b in radians
104  template<class T>
106  (
107  const Vector<T>& a,
108  const Vector<T>& b,
109  const T& tolerance = small
110  )
111  {
112  scalar cosPhi = (a & b)/(mag(a)*mag(b) + tolerance);
113 
114  // Enforce bounding between -1 and 1
115  return min(max(cosPhi, -1), 1);
116  }
117 
118  //- Calculate angle between a and b in radians
119  template<class T>
121  (
122  const Vector<T>& a,
123  const Vector<T>& b,
124  const T& tolerance = small
125  )
126  {
127  scalar cosPhi = (a & b)/(mag(a)*mag(b) + tolerance);
128 
129  // Enforce bounding between -1 and 1
130  return acos( min(max(cosPhi, -1), 1) );
131  }
132 
133  //- Calculate angle between a and b in degrees
134  template<class T>
136  (
137  const Vector<T>& a,
138  const Vector<T>& b,
139  const T& tolerance = small
140  )
141  {
142  return radToDeg(radAngleBetween(a, b, tolerance));
143  }
144 
145 } // End namespace vectorTools
146 
147 
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 
150 } // End namespace Foam
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 #endif
155 
156 // ************************************************************************* //
dimensionedScalar acos(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
bool areOrthogonal(const Vector< T > &a, const Vector< T > &b, const T &tolerance=small)
Test if a and b are orthogonal: a.b = 0.
Definition: vectorTools.H:71
Unit conversion functions.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
T cosPhi(const Vector< T > &a, const Vector< T > &b, const T &tolerance=small)
Calculate angle between a and b in radians.
Definition: vectorTools.H:105
bool areObtuse(const Vector< T > &a, const Vector< T > &b)
Test if angle between a and b is obtuse: a.b < 0.
Definition: vectorTools.H:94
T degAngleBetween(const Vector< T > &a, const Vector< T > &b, const T &tolerance=small)
Calculate angle between a and b in degrees.
Definition: vectorTools.H:135
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
T radAngleBetween(const Vector< T > &a, const Vector< T > &b, const T &tolerance=small)
Calculate angle between a and b in radians.
Definition: vectorTools.H:120
Functions for analysing the relationships between vectors.
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross product operators.
Definition: Vector.H:57
bool areAcute(const Vector< T > &a, const Vector< T > &b)
Test if angle between a and b is acute: a.b > 0.
Definition: vectorTools.H:83
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
bool areParallel(const Vector< T > &a, const Vector< T > &b, const T &tolerance=small)
Test if a and b are parallel: a^b = 0.
Definition: vectorTools.H:54
dimensioned< scalar > mag(const dimensioned< Type > &)
Namespace for OpenFOAM.