Loading...
Searching...
No Matches
complexVectorField.C
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 OpenFOAM Foundation
9 Copyright (C) 2019-2023 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
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
38}
39
40
41// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
42
43void Foam::zip
44(
45 complexVectorField& result,
46 const UList<vector>& realValues,
47 const UList<vector>& imagValues
48)
49{
50 const label len = result.size();
51
52 #ifdef FULLDEBUG
53 if (len != realValues.size() || len != imagValues.size())
54 {
56 << "Components sizes do not match: " << len << " ("
57 << realValues.size() << ' ' << imagValues.size() << ')' << nl
58 << abort(FatalError);
59 }
60 #endif
61
62 for (label i=0; i < len; ++i)
63 {
64 result[i] = Foam::zip(realValues[i], imagValues[i]);
65 }
66}
67
68
69void Foam::zip
70(
71 complexVectorField& result,
72 const UList<vector>& realValues,
73 const vector& imagValue
74)
75{
76 const label len = result.size();
77
78 #ifdef FULLDEBUG
79 if (len != realValues.size())
80 {
82 << "Components sizes do not match: " << len << " != "
83 << realValues.size() << nl
84 << abort(FatalError);
85 }
86 #endif
87
88 for (label i=0; i < len; ++i)
89 {
90 result[i] = Foam::zip(realValues[i], imagValue);
91 }
92}
93
94
95void Foam::zip
96(
97 complexVectorField& result,
98 const vector& realValue,
99 const UList<vector>& imagValues
100)
101{
102 const label len = result.size();
103
104 #ifdef FULLDEBUG
105 if (len != imagValues.size())
106 {
108 << "Components sizes do not match: " << len << " != "
109 << imagValues.size() << nl
110 << abort(FatalError);
111 }
112 #endif
113
114 for (label i=0; i < len; ++i)
115 {
116 result[i] = Foam::zip(realValue, imagValues[i]);
117 }
118}
119
120
122(
123 const UList<vector>& realValues,
124 const UList<vector>& imagValues
125)
126{
127 complexVectorField result(realValues.size());
129 Foam::zip(result, realValues, imagValues);
130
131 return result;
132}
133
134
136(
137 const UList<vector>& realValues,
138 const vector& imagValue
139)
140{
141 complexVectorField result(realValues.size());
143 Foam::zip(result, realValues, imagValue);
144
145 return result;
146}
147
148
150(
151 const vector& realValue,
152 const UList<vector>& imagValues
153)
154{
155 complexVectorField result(imagValues.size());
157 Foam::zip(result, realValue, imagValues);
158
159 return result;
160}
161
162
163Foam::vectorField Foam::ReImSum(const UList<complexVector>& cmplx)
164{
165 vectorField result(cmplx.size());
166
167 std::transform
168 (
169 cmplx.cbegin(),
170 cmplx.cend(),
171 result.begin(),
172 [](const complexVector& c)
173 {
174 return vector(c.x().cmptSum(), c.y().cmptSum(), c.z().cmptSum());
176 );
177
178 return result;
179}
180
181
182Foam::vectorField Foam::Re(const UList<complexVector>& cmplx)
183{
184 vectorField result(cmplx.size());
185
186 std::transform
187 (
188 cmplx.cbegin(),
189 cmplx.cend(),
190 result.begin(),
191 [](const complexVector& c)
192 {
193 return vector(c.x().real(), c.y().real(), c.z().real());
195 );
196
197 return result;
198}
199
200
201Foam::vectorField Foam::Im(const UList<complexVector>& cmplx)
202{
203 vectorField result(cmplx.size());
204
205 std::transform
206 (
207 cmplx.cbegin(),
208 cmplx.cend(),
209 result.begin(),
210 [](const complexVector& c)
211 {
212 return vector(c.x().imag(), c.y().imag(), c.z().imag());
214 );
215
216 return result;
217}
218
219
220Foam::complexVectorField Foam::operator^
221(
222 const UList<vector>& vec,
223 const UList<complexVector>& cmplx
224)
225{
226 const label len = cmplx.size();
227
228 #ifdef FULLDEBUG
229 if (len != vec.size())
230 {
232 << "Parameter sizes do not match: " << vec.size()
233 << " != " << len << nl
234 << abort(FatalError);
235 }
236 #endif
237
238 complexVectorField result(len);
239
240 for (label i=0; i < len; ++i)
241 {
242 result[i] = (vec[i] ^ cmplx[i]);
243 }
244
245 return result;
246}
247
248
249// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
const_iterator cend() const noexcept
Return const_iterator to end traversing the constant UList.
Definition UListI.H:468
const_iterator cbegin() const noexcept
Return const_iterator to begin traversing the constant UList.
Definition UListI.H:424
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for OpenFOAM.
scalarField ReImSum(const UList< complex > &cmplx)
Sum real and imag components.
scalarField Im(const UList< complex > &cmplx)
Extract imag component.
Field< complexVector > complexVectorField
Specialisation of Field<T> for complexVector.
Vector< complex > complexVector
A Vector of complex values with 'scalar' precision.
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
complexField ComplexField(const UList< scalar > &realValues, const UList< scalar > &imagValues)
Create complex field by zipping two lists of real/imag values.
void zip(FieldField< Field, SphericalTensor< Cmpt > > &result, const FieldField< Field, Cmpt > &ii)
Zip together sphericalTensor field field from components.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addCompoundToRunTimeSelectionTable(Type, Tag)
Add compound to selection tables, lookup using typeName.
Definition token.H:1463
#define defineCompoundTypeName(Type, UnusedTag)
Define compound using Type for its name.
Definition token.H:1451