PSCF v1.2
prdc/cpu/complex.h
1#ifndef PRDC_CPU_COMPLEX_H
2#define PRDC_CPU_COMPLEX_H
3
4/*
5* PSCF Package - Polymer Self-Consistent Field Theory
6*
7* Copyright 2016 - 2022, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include <pscf/math/complex.h>
12
13#include <fftw3.h>
14#include <complex>
15#include <iostream>
16
17namespace Pscf {
18
19 /*
20 * Types Cpu::Complex and Cpu::Real are defined in prdc/cpu/types.h
21 * as aliases for fftw_complex and double, respectively.
22 */
23
24 // Real and imaginary components
25
33 template <> inline
34 double real<fftw_complex, double>(fftw_complex const& a)
35 { return a[0]; }
36
44 template <> inline
45 double imag<fftw_complex, double>(fftw_complex const& a)
46 { return a[1]; }
47
48 // Absolute magnitude
49
57 template <> inline
58 double abs<fftw_complex, double>(fftw_complex const& a)
59 { return sqrt(a[0] * a[0] + a[1] * a[1]); }
60
68 template <> inline
69 double absSq<fftw_complex, double>(fftw_complex const& a)
70 { return (a[0] * a[0] + a[1] * a[1]); }
71
72 // fftw_complex Conjugation
73
82 template <> inline
83 void conj(fftw_complex& z, fftw_complex const& a)
84 {
85 z[0] = a[0];
86 z[1] = -a[1];
87 }
88
96 template <> inline
97 void conj(fftw_complex& a)
98 {
99 a[0] = a[0];
100 a[1] = -a[1];
101 }
102
103 // Assignment
104
114 template <> inline
115 void assign(fftw_complex& z, double const& a, double const& b)
116 {
117 z[0] = a;
118 z[1] = b;
119 }
120
129 template <> inline
130 void assign(fftw_complex& z, double const& a)
131 {
132 z[0] = a;
133 z[1] = 0;
134 }
135
144 template <> inline
145 void assign(fftw_complex& z, fftw_complex const& a)
146 {
147 z[0] = a[0];
148 z[1] = a[1];
149 }
150
159 template <> inline
160 void assign(fftw_complex & z, std::complex<double> const& a)
161 {
162 z[0] = a.real();
163 z[1] = a.imag();
164 }
165
174 template <> inline
175 void assign(std::complex<double> & z, fftw_complex const& a)
176 { z = std::complex<double>(a[0], a[1]); }
177
178 // Addition
179
189 template <> inline
190 void add(fftw_complex& z, fftw_complex const& a, fftw_complex const& b)
191 {
192 z[0] = a[0] + b[0];
193 z[1] = a[1] + b[1];
194 }
195
205 template <> inline
206 void add(fftw_complex& z, fftw_complex const& a, double const& b)
207 {
208 z[0] = a[0] + b;
209 z[1] = a[1];
210 }
211
220 template <> inline
221 void addEq(fftw_complex& a, fftw_complex const& b)
222 {
223 a[0] += b[0];
224 a[1] += b[1];
225 }
226
235 template <> inline
236 void addEq(fftw_complex& a, double const& b)
237 {
238 a[0] += b;
239 }
240
241 // Subtraction
242
252 template <> inline
253 void sub(fftw_complex& z, fftw_complex const& a, fftw_complex const& b)
254 {
255 z[0] = a[0] - b[0];
256 z[1] = a[1] - b[1];
257 }
258
268 template <> inline
269 void sub(fftw_complex& z, fftw_complex const& a, double const& b)
270 {
271 z[0] = a[0] - b;
272 z[1] = a[1];
273 }
274
283 template <> inline
284 void subEq(fftw_complex & a, fftw_complex const& b)
285 {
286 a[0] -= b[0];
287 a[1] -= b[1];
288 }
289
298 template <> inline
299 void subEq(fftw_complex & a, double const& b)
300 {
301 a[0] -= b;
302 }
303
314 template <> inline
315 double absSqDiff(fftw_complex const& a, fftw_complex const& b)
316 {
317 fftw_complex z;
318 sub(z, a, b);
320 }
321
322 // Multiplication
323
333 template <> inline
334 void mul(fftw_complex& z, fftw_complex const& a, fftw_complex const& b)
335 {
336 z[0] = a[0] * b[0] - a[1] * b[1];
337 z[1] = a[1] * b[0] + a[0] * b[1];
338 }
339
349 template <> inline
350 void mul(fftw_complex& z, fftw_complex const& a, double const& b)
351 {
352 z[0] = a[0]*b;
353 z[1] = a[1]*b;
354 }
355
364 template <> inline
365 void mulEq(fftw_complex & a, fftw_complex const& b)
366 {
367 double a0;
368 a0 = a[0] * b[0] - a[1] * b[1];
369 a[1] = a[1] * b[0] + a[0] * b[1];
370 a[0] = a0;
371 }
372
381 template <> inline
382 void mulEq(fftw_complex & a, double const& b)
383 {
384 a[0] *= b;
385 a[1] *= b;
386 }
387
396 template <> inline
397 void square(fftw_complex& z, fftw_complex const& a)
398 {
399 z[0] = a[0] * a[0] - a[1] * a[1];
400 z[1] = 2.0 * a[1] * a[0];
401 }
402
403 // Division
404
414 template <> inline
415 void div(fftw_complex& z, fftw_complex const& a, fftw_complex const& b)
416 {
417 double bSq = b[0] * b[0] + b[1] * b[1];
418 z[0] = (a[0] * b[0] + a[1] * b[1])/bSq;
419 z[1] = (a[1] * b[0] - a[0] * b[1])/bSq;
420 }
421
431 template <> inline
432 void div(fftw_complex& z, fftw_complex const& a, double const& b)
433 {
434 z[0] = a[0]/b;
435 z[1] = a[1]/b;
436 }
437
446 template <> inline
447 void divEq(fftw_complex & a, fftw_complex const & b)
448 {
449 double bSq = b[0] * b[0] + b[1] * b[1];
450 double a0 = (a[0] * b[0] + a[1] * b[1])/bSq;
451 a[1] = (a[1] * b[0] - a[0] * b[1])/bSq;
452 a[0] = a0;
453 }
454
463 template <> inline
464 void divEq(fftw_complex & a, double const& b)
465 {
466 a[0] /= b;
467 a[1] /= b;
468 }
469
476 std::istream& operator >> (std::istream& is, fftw_complex & z);
477
485 std::ostream& operator << (std::ostream& os, fftw_complex const & z);
486
487 #if 0
497 template <typename Archive>
498 void serialize(Archive& ar, double z[2], const unsigned int version = 0)
499 {
500 serialize(ar, z[0], version);
501 serialize(ar, z[1], version);
502 }
503 #endif
504
505} // namespace Pscf
506#endif
double imag< fftw_complex, double >(fftw_complex const &a)
Return the imaginary part of a complex number.
double absSq< fftw_complex, double >(fftw_complex const &a)
Return square of absolute magnitude of a complex number.
double real< fftw_complex, double >(fftw_complex const &a)
Return the real part of a complex number.
double abs< fftw_complex, double >(fftw_complex const &a)
Return absolute magnitude of a complex number.
void div(CT &z, CT const &a, CT const &b)
Division of two complex numbers, z = a / b .
void square(CT &z, CT const &a)
Compute complex square of a complex number, z = a * a.
void conj(CT &z, CT const &a)
Compute complex conjugate, z = a^*.
void divEq(CT &a, CT const &b)
In place division of two complex number, a /= b.
void subEq(CT &a, CT const &b)
In place subtraction of two complex numbers, a -= b.
void mulEq(CT &a, CT const &b)
In place multiplication of two complex number, a *= b.
void assign(CT &z, RT const &a, RT const &b)
Create a complex number from real and imaginary parts, z = a + ib.
void addEq(CT &a, CT const &b)
In place addition of complex numbers, a += b.
RT absSqDiff(CT const &a, CT const &b)
Return square of the absolute magnitude of a complex difference.
void add(CT &z, CT const &a, CT const &b)
Addition of two complex numbers, z = a + b.
void sub(CT &z, CT const &a, CT const &b)
Subtraction of two complex numbers, z = a - b.
void mul(CT &z, CT const &a, CT const &b)
Multiplication of two complex numbers, z = a * b.
PSCF package top-level namespace.
Definition param_pc.dox:1
std::istream & operator>>(std::istream &in, Pair< Data > &pair)
Input a Pair from an istream.
Definition Pair.h:44
std::ostream & operator<<(std::ostream &out, const Pair< Data > &pair)
Output a Pair to an ostream, without line breaks.
Definition Pair.h:57