Loading...
Searching...
No Matches
BarycentricTensorI.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) 2017 OpenFOAM Foundation
9 Copyright (C) 2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
31template<class Cmpt>
33:
35{}
36
37
38template<class Cmpt>
40(
41 const Barycentric<Cmpt>& x,
42 const Barycentric<Cmpt>& y,
43 const Barycentric<Cmpt>& z
44)
45{
46 this->v_[XA] = x.a();
47 this->v_[XB] = x.b();
48 this->v_[XC] = x.c();
49 this->v_[XD] = x.d();
50
51 this->v_[YA] = y.a();
52 this->v_[YB] = y.b();
53 this->v_[YC] = y.c();
54 this->v_[YD] = y.d();
55
56 this->v_[ZA] = z.a();
57 this->v_[ZB] = z.b();
58 this->v_[ZC] = z.c();
59 this->v_[ZD] = z.d();
60}
61
62
63template<class Cmpt>
65(
66 const Vector<Cmpt>& a,
67 const Vector<Cmpt>& b,
68 const Vector<Cmpt>& c,
69 const Vector<Cmpt>& d
70)
71{
72 this->v_[XA] = a.x();
73 this->v_[XB] = b.x();
74 this->v_[XC] = c.x();
75 this->v_[XD] = d.x();
76
77 this->v_[YA] = a.y();
78 this->v_[YB] = b.y();
79 this->v_[YC] = c.y();
80 this->v_[YD] = d.y();
81
82 this->v_[ZA] = a.z();
83 this->v_[ZB] = b.z();
84 this->v_[ZC] = c.z();
85 this->v_[ZD] = d.z();
86}
87
88
89// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90
91template<class Cmpt>
93{
94 return
96 (
97 this->v_[XA],
98 this->v_[XB],
99 this->v_[XC],
100 this->v_[XD]
101 );
102}
103
104
105template<class Cmpt>
107{
108 return
110 (
111 this->v_[YA],
112 this->v_[YB],
113 this->v_[YC],
114 this->v_[YD]
115 );
116}
117
118
119template<class Cmpt>
121{
122 return
124 (
125 this->v_[ZA],
126 this->v_[ZB],
127 this->v_[ZC],
128 this->v_[ZD]
129 );
130}
131
132
133template<class Cmpt>
135{
136 return Vector<Cmpt>(this->v_[XA], this->v_[YA], this->v_[ZA]);
137}
138
139
140template<class Cmpt>
142{
143 return Vector<Cmpt>(this->v_[XB], this->v_[YB], this->v_[ZB]);
144}
145
146
147template<class Cmpt>
149{
150 return Vector<Cmpt>(this->v_[XC], this->v_[YC], this->v_[ZC]);
151}
152
153
154template<class Cmpt>
156{
157 return Vector<Cmpt>(this->v_[XD], this->v_[YD], this->v_[ZD]);
158}
159
160
161// NB: same workaround for gcc (11+) failure on (tensor dot vector)
162// - not sure it will indeed be required here as well.
163template<class Cmpt>
164#if defined(__GNUC__) && !defined(__clang__)
165__attribute__((optimize("no-tree-vectorize")))
166#endif
169{
170 const BarycentricTensor<Cmpt>& t = *this;
171
172 return Vector<Cmpt>
173 (
174 (bry.a()*t.xa() + bry.b()*t.xb() + bry.c()*t.xc() + bry.d()*t.xd()),
175 (bry.a()*t.ya() + bry.b()*t.yb() + bry.c()*t.yc() + bry.d()*t.yd()),
176 (bry.a()*t.za() + bry.b()*t.zb() + bry.c()*t.zc() + bry.d()*t.zd())
177 );
178}
179
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183namespace Foam
184{
185
186// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
187
188// Transform Barycentric coordinates to Vector
189template<class Cmpt>
190inline Vector<Cmpt> operator&
191(
193 const Barycentric<Cmpt>& b
194)
195{
196 return t.inner(b);
197}
198
199
200// Transform Vector to Barycentric coordinates.
201// Caution: the tensor must be inverse one (see particle.C)
202template<class Cmpt>
203inline Barycentric<Cmpt> operator&
204(
205 const Vector<Cmpt>& v,
206 const BarycentricTensor<Cmpt>& T
207)
208{
209 return Barycentric<Cmpt>(v & T.a(), v & T.b(), v & T.c(), v & T.d());
210}
211
212
213// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214
215} // End namespace Foam
216
217// ************************************************************************* //
scalar y
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
const Cmpt & ya() const noexcept
const Cmpt & xc() const noexcept
const Cmpt & zd() const noexcept
Barycentric< Cmpt > y() const
Vector< Cmpt > d() const
const Cmpt & xb() const noexcept
Vector< Cmpt > b() const
Vector< Cmpt > c() const
const Cmpt & xa() const noexcept
Barycentric< Cmpt > z() const
BarycentricTensor()=default
Default construct.
const Cmpt & yc() const noexcept
Vector< Cmpt > inner(const Barycentric< Cmpt > &bry) const
Tensor/barycentric inner product.
const Cmpt & zc() const noexcept
const Cmpt & zb() const noexcept
const Cmpt & yd() const noexcept
const Cmpt & za() const noexcept
const Cmpt & xd() const noexcept
const Cmpt & yb() const noexcept
Barycentric< Cmpt > x() const
Vector< Cmpt > a() const
Templated 3D Barycentric derived from VectorSpace. Has 4 components, one of which is redundant.
Definition Barycentric.H:53
const Cmpt & b() const noexcept
const Cmpt & d() const noexcept
const Cmpt & a() const noexcept
const Cmpt & c() const noexcept
MatrixSpace< BarycentricTensor< Cmpt >, Cmpt, Mrows, Ncols > msType
Definition MatrixSpace.H:65
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
Cmpt inner(const Vector< Cmpt > &v2) const
Scalar-product of this with another Vector.
Definition VectorI.H:146
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
const volScalarField & T
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & b