complexField.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 -------------------------------------------------------------------------------
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 #include "complexField.H"
31 
32 #define TEMPLATE
33 #include "FieldFunctionsM.C"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineCompoundTypeName(List<complex>, complexList);
41 }
42 
43 
44 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
45 
46 void Foam::zip
47 (
48  complexField& result,
49  const UList<scalar>& realValues,
50  const UList<scalar>& imagValues
51 )
52 {
53  const label len = result.size();
54 
55  #ifdef FULLDEBUG
56  if (len != realValues.size() || len != imagValues.size())
57  {
59  << "Components sizes do not match: " << len << " ("
60  << realValues.size() << ' ' << imagValues.size() << ')' << nl
61  << abort(FatalError);
62  }
63  #endif
64 
65  for (label i=0; i < len; ++i)
66  {
67  result[i].real(realValues[i]);
68  result[i].imag(imagValues[i]);
69  }
70 }
71 
72 
73 void Foam::zip
74 (
75  complexField& result,
76  const UList<scalar>& realValues,
77  const scalar imagValue
78 )
79 {
80  const label len = result.size();
81 
82  #ifdef FULLDEBUG
83  if (len != realValues.size())
84  {
86  << "Components sizes do not match: " << len
87  << " != " << realValues.size() << nl
88  << abort(FatalError);
89  }
90  #endif
91 
92  for (label i=0; i < len; ++i)
93  {
94  result[i].real(realValues[i]);
95  result[i].imag(imagValue);
96  }
97 }
98 
99 
100 void Foam::zip
101 (
102  complexField& result,
103  const scalar realValue,
104  const UList<scalar>& imagValues
105 )
106 {
107  const label len = result.size();
108 
109  #ifdef FULLDEBUG
110  if (len != imagValues.size())
111  {
113  << "Components sizes do not match: " << len
114  << " != " << imagValues.size() << nl
115  << abort(FatalError);
116  }
117  #endif
118 
119  for (label i=0; i < len; ++i)
120  {
121  result[i].real(realValue);
122  result[i].imag(imagValues[i]);
123  }
124 }
125 
126 
127 void Foam::unzip
128 (
129  const UList<complex>& input,
130  scalarField& realValues,
131  scalarField& imagValues
132 )
133 {
134  const label len = input.size();
135 
136  #ifdef FULLDEBUG
137  if (len != realValues.size() || len != imagValues.size())
138  {
140  << "Components sizes do not match: " << len << " ("
141  << realValues.size() << ' ' << imagValues.size() << ')' << nl
142  << abort(FatalError);
143  }
144  #endif
145 
146  for (label i=0; i < len; ++i)
147  {
148  realValues[i] = input[i].real();
149  imagValues[i] = input[i].imag();
150  }
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
157 (
158  const UList<scalar>& realValues,
159  const UList<scalar>& imagValues
160 )
161 {
162  complexField result(realValues.size());
163 
164  Foam::zip(result, realValues, imagValues);
165 
166  return result;
167 }
168 
169 
171 (
172  const UList<scalar>& realValues,
173  const scalar imagValue
174 )
175 {
176  complexField result(realValues.size());
177 
178  Foam::zip(result, realValues, imagValue);
179 
180  return result;
181 }
182 
183 
185 (
186  const scalar realValue,
187  const UList<scalar>& imagValues
188 )
189 {
190  complexField result(imagValues.size());
192  Foam::zip(result, realValue, imagValues);
193 
194  return result;
195 }
196 
197 
198 Foam::scalarField Foam::ReImSum(const UList<complex>& cmplx)
199 {
200  scalarField result(cmplx.size());
201 
203  (
204  cmplx.cbegin(),
205  cmplx.cend(),
206  result.begin(),
207  [](const complex& c) { return c.cmptSum(); }
208  );
209 
210  return result;
211 }
212 
213 
214 Foam::scalarField Foam::Re(const UList<complex>& cmplx)
215 {
216  scalarField result(cmplx.size());
217 
219  (
220  cmplx.cbegin(),
221  cmplx.cend(),
222  result.begin(),
223  [](const complex& c) { return c.real(); }
224  );
225 
226  return result;
227 }
228 
229 
230 Foam::scalarField Foam::Im(const UList<complex>& cmplx)
231 {
232  scalarField result(cmplx.size());
233 
235  (
236  cmplx.cbegin(),
237  cmplx.cend(),
238  result.begin(),
239  [](const complex& c) { return c.imag(); }
240  );
241 
242  return result;
243 }
244 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 namespace Foam
249 {
250 
251 template<>
252 complex sumProd(const UList<complex>& f1, const UList<complex>& f2)
253 {
254  complex result(0);
255  if (f1.size() && (f1.size() == f2.size()))
256  {
257  // std::inner_product
258  TFOR_ALL_S_OP_F_OP_F(complex, result, +=, complex, f1, *, complex, f2)
259  }
260  return result;
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
268 
271 
272 
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 
296 
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 } // End namespace Foam
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 #include "undefFieldFunctionsM.H"
305 
306 // ************************************************************************* //
dimensionedScalar tanh(const dimensionedScalar &ds)
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensionedScalar acos(const dimensionedScalar &ds)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar log(const dimensionedScalar &ds)
addCompoundToRunTimeSelectionTable(List< complex >, complexList)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
defineCompoundTypeName(List< complex >, complexList)
Field< complex > complexField
Specialisation of Field<T> for complex.
Definition: complexField.H:50
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void zip(FieldField< Field, SphericalTensor< Cmpt >> &result, const FieldField< Field, Cmpt > &ii)
Zip together sphericalTensor field field from components.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
scalarField ReImSum(const UList< complex > &cmplx)
Sum real and imag components.
Definition: complexField.C:191
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
scalarField Im(const UList< complex > &cmplx)
Extract imag component.
Definition: complexField.C:223
complex sumProd(const UList< complex > &f1, const UList< complex > &f2)
Sum product.
Definition: complexField.C:245
void unzip(const FieldField< Field, SphericalTensor< Cmpt >> &input, FieldField< Field, Cmpt > &ii)
Unzip sphericalTensor field field into components.
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
Generic templated field type.
Definition: Field.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar asinh(const dimensionedScalar &ds)
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:33
complexField ComplexField(const UList< scalar > &realValues, const UList< scalar > &imagValues)
Create complex field by zipping two lists of real/imag values.
Definition: complexField.C:150
dimensionedScalar atanh(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
dimensionedScalar sin(const dimensionedScalar &ds)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
#define TFOR_ALL_S_OP_F_OP_F(typeS, s, OP1, typeF1, f1, OP2, typeF2, f2)
Definition: FieldM.H:540
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow6(const dimensionedScalar &ds)
A complex number, similar to the C++ complex type.
Definition: complex.H:70
dimensionedScalar cosh(const dimensionedScalar &ds)
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
dimensionedScalar log10(const dimensionedScalar &ds)
Namespace for OpenFOAM.