complexI.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 | Copyright (C) 2011 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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
27 
28 namespace Foam
29 {
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 {}
35 
36 
37 inline complex::complex(const scalar Re, const scalar Im)
38 :
39  re(Re),
40  im(Im)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
45 
46 inline scalar complex::Re() const
47 {
48  return re;
49 }
50 
51 
52 inline scalar complex::Im() const
53 {
54  return im;
55 }
56 
57 
58 inline scalar& complex::Re()
59 {
60  return re;
61 }
62 
63 
64 inline scalar& complex::Im()
65 {
66  return im;
67 }
68 
69 
71 {
72  return complex(re, -im);
73 }
74 
75 
76 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
77 
78 inline const complex& complex::operator=(const complex& c)
79 {
80  re = c.re;
81  im = c.im;
82  return *this;
83 }
84 
85 
86 inline void complex::operator+=(const complex& c)
87 {
88  re += c.re;
89  im += c.im;
90 }
91 
92 
93 inline void complex::operator-=(const complex& c)
94 {
95  re -= c.re;
96  im -= c.im;
97 }
98 
99 
100 inline void complex::operator*=(const complex& c)
101 {
102  *this = (*this)*c;
103 }
104 
105 
106 inline void complex::operator/=(const complex& c)
107 {
108  *this = *this/c;
109 }
110 
111 
112 inline const complex& complex::operator=(const scalar s)
113 {
114  re = s;
115  im = 0.0;
116  return *this;
117 }
118 
119 
120 inline void complex::operator+=(const scalar s)
121 {
122  re += s;
123 }
124 
125 
126 inline void complex::operator-=(const scalar s)
127 {
128  re -= s;
129 }
130 
131 
132 inline void complex::operator*=(const scalar s)
133 {
134  re *= s;
135  im *= s;
136 }
137 
138 
139 inline void complex::operator/=(const scalar s)
140 {
141  re /= s;
142  im /= s;
143 }
144 
145 
147 {
148  return conjugate();
149 }
150 
151 
152 inline bool complex::operator==(const complex& c) const
153 {
154  return (equal(re, c.re) && equal(im, c.im));
155 }
156 
157 
158 inline bool complex::operator!=(const complex& c) const
159 {
160  return !operator==(c);
161 }
162 
163 
164 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
165 
166 
167 inline scalar magSqr(const complex& c)
168 {
169  return (c.re*c.re + c.im*c.im);
170 }
171 
172 
173 inline complex sqr(const complex& c)
174 {
175  return c * c;
176 }
177 
178 
179 inline scalar mag(const complex& c)
180 {
181  return sqrt(magSqr(c));
182 }
183 
184 
185 inline const complex& max(const complex& c1, const complex& c2)
186 {
187  if (mag(c1) > mag(c2))
188  {
189  return c1;
190  }
191  else
192  {
193  return c2;
194  }
195 }
196 
197 
198 inline const complex& min(const complex& c1, const complex& c2)
199 {
200  if (mag(c1) < mag(c2))
201  {
202  return c1;
203  }
204  else
205  {
206  return c2;
207  }
208 }
209 
210 
211 inline complex limit(const complex& c1, const complex& c2)
212 {
213  return complex(limit(c1.re, c2.re), limit(c1.im, c2.im));
214 }
215 
216 
217 inline const complex& sum(const complex& c)
218 {
219  return c;
220 }
221 
222 
223 template<class Cmpt>
224 class Tensor;
225 
226 inline complex transform(const Tensor<scalar>&, const complex c)
227 {
228  return c;
229 }
230 
231 
232 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
233 
234 inline complex operator+(const complex& c1, const complex& c2)
235 {
236  return complex
237  (
238  c1.re + c2.re,
239  c1.im + c2.im
240  );
241 }
242 
243 
244 inline complex operator-(const complex& c)
245 {
246  return complex
247  (
248  -c.re,
249  -c.im
250  );
251 }
252 
253 
254 inline complex operator-(const complex& c1, const complex& c2)
255 {
256  return complex
257  (
258  c1.re - c2.re,
259  c1.im - c2.im
260  );
261 }
262 
263 
264 inline complex operator*(const complex& c1, const complex& c2)
265 {
266  return complex
267  (
268  c1.re*c2.re - c1.im*c2.im,
269  c1.im*c2.re + c1.re*c2.im
270  );
271 }
272 
273 
274 inline complex operator/(const complex& c1, const complex& c2)
275 {
276  scalar sqrC2 = magSqr(c2);
277 
278  return complex
279  (
280  (c1.re*c2.re + c1.im*c2.im)/sqrC2,
281  (c1.im*c2.re - c1.re*c2.im)/sqrC2
282  );
283 }
284 
285 
286 inline complex operator*(const scalar s, const complex& c)
287 {
288  return complex(s*c.re, s*c.im);
289 }
290 
291 
292 inline complex operator*(const complex& c, const scalar s)
293 {
294  return complex(s*c.re, s*c.im);
295 }
296 
297 
298 inline complex operator/(const complex& c, const scalar s)
299 {
300  return complex(c.re/s, c.im/s);
301 }
302 
303 
304 inline complex operator/(const scalar s, const complex& c)
305 {
306  return complex(s/c.re, s/c.im);
307 }
308 
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 } // End namespace Foam
313 
314 // ************************************************************************* //
Foam::Im
scalarField Im(const UList< complex > &cf)
Definition: complexFields.C:110
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::complex::conjugate
complex conjugate() const
Definition: complexI.H:70
Foam::complex::operator!
complex operator!() const
Definition: complexI.H:146
Foam::complex::re
scalar re
Real and imaginary parts of the complex number.
Definition: complex.H:81
Foam::Re
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Foam::complex::operator=
const complex & operator=(const complex &)
Definition: complexI.H:78
Foam::complex::operator!=
bool operator!=(const complex &) const
Definition: complexI.H:158
Foam::complex::im
scalar im
Definition: complex.H:81
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::complex::operator*=
void operator*=(const complex &)
Definition: complexI.H:100
Foam::complex::operator/=
void operator/=(const complex &)
Definition: complexI.H:106
Foam::transform
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
Foam::complex::operator+=
void operator+=(const complex &)
Definition: complexI.H:86
Foam::operator-
tmp< fvMatrix< Type > > operator-(const fvMatrix< Type > &)
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::complex::operator==
bool operator==(const complex &) const
Definition: complexI.H:152
Foam::constant::atomic::re
const dimensionedScalar re
Classical electron radius: default SI units: [m].
Foam::complex::Re
scalar Re() const
Definition: complexI.H:46
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::complex
Extension to the c++ complex library type.
Definition: complex.H:76
Foam::constant::physicoChemical::c2
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::operator*
tmp< fvMatrix< Type > > operator*(const DimensionedField< scalar, volMesh > &, const fvMatrix< Type > &)
Foam::operator/
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
Definition: dimensionedScalar.C:65
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::complex::operator-=
void operator-=(const complex &)
Definition: complexI.H:93
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::limit
complex limit(const complex &, const complex &)
Definition: complexI.H:211
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::complex::Im
scalar Im() const
Definition: complexI.H:52
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::equal
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:78
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::operator+
tmp< fvMatrix< Type > > operator+(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::complex::complex
complex()
Construct null.
Definition: complexI.H:33