BarycentricTensor2DI.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) 2021 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 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
28 template<class Cmpt>
30 {}
31 
32 
33 template<class Cmpt>
35 :
37 {}
38 
39 
40 template<class Cmpt>
42 (
43  const Barycentric2D<Cmpt>& x,
44  const Barycentric2D<Cmpt>& y,
45  const Barycentric2D<Cmpt>& z
46 )
47 {
48  this->v_[XA] = x.a();
49  this->v_[XB] = x.b();
50  this->v_[XC] = x.c();
51 
52  this->v_[YA] = y.a();
53  this->v_[YB] = y.b();
54  this->v_[YC] = y.c();
55 
56  this->v_[ZA] = z.a();
57  this->v_[ZB] = z.b();
58  this->v_[ZC] = z.c();
59 }
60 
61 
62 template<class Cmpt>
64 (
65  const Vector<Cmpt>& a,
66  const Vector<Cmpt>& b,
67  const Vector<Cmpt>& c
68 )
69 {
70  this->v_[XA] = a.x();
71  this->v_[XB] = b.x();
72  this->v_[XC] = c.x();
73 
74  this->v_[YA] = a.y();
75  this->v_[YB] = b.y();
76  this->v_[YC] = c.y();
77 
78  this->v_[ZA] = a.z();
79  this->v_[ZB] = b.z();
80  this->v_[ZC] = c.z();
81 }
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
86 template<class Cmpt>
88 {
89  return
91  (
92  this->v_[XA],
93  this->v_[XB],
94  this->v_[XC]
95  );
96 }
97 
98 
99 template<class Cmpt>
101 {
102  return
104  (
105  this->v_[YA],
106  this->v_[YB],
107  this->v_[YC]
108  );
109 }
110 
111 
112 template<class Cmpt>
114 {
115  return
117  (
118  this->v_[ZA],
119  this->v_[ZB],
120  this->v_[ZC]
121  );
122 }
123 
124 
125 template<class Cmpt>
127 {
128  return Vector<Cmpt>(this->v_[XA], this->v_[YA], this->v_[ZA]);
129 }
130 
131 
132 template<class Cmpt>
134 {
135  return Vector<Cmpt>(this->v_[XB], this->v_[YB], this->v_[ZB]);
136 }
137 
138 
139 template<class Cmpt>
141 {
142  return Vector<Cmpt>(this->v_[XC], this->v_[YC], this->v_[ZC]);
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 
148 namespace Foam
149 {
150 
151 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
152 
153 template<class Cmpt>
154 inline Vector<Cmpt> operator&
155 (
157  const Barycentric2D<Cmpt>& b
158 )
159 {
160  return Vector<Cmpt>(T.x() & b, T.y() & b, T.z() & b);
161 }
162 
163 
164 template<class Cmpt>
165 inline Barycentric2D<Cmpt> operator&
166 (
167  const Vector<Cmpt>& v,
169 )
170 {
171  return Barycentric2D<Cmpt>(v & T.a(), v & T.b(), v & T.c());
172 }
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 } // End namespace Foam
178 
179 // ************************************************************************* //
scalar y
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant.
Definition: Barycentric2D.H:53
const Cmpt & a() const
const Cmpt & b() const
const Cmpt & c() const
Templated 3x3 tensor derived from VectorSpace. Has 9 components. Can represent a barycentric transfor...
Barycentric2D< Cmpt > z() const
Vector< Cmpt > b() const
Vector< Cmpt > c() const
Barycentric2D< Cmpt > y() const
Barycentric2D< Cmpt > x() const
BarycentricTensor2D()
Construct null.
Vector< Cmpt > a() const
Templated matrix space.
Definition: MatrixSpace.H:58
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:60
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
volScalarField & b
Definition: createFields.H:25
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)