forked from DmitryLyakh/TAL_SH
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtalsh_complex.h
460 lines (437 loc) · 11.9 KB
/
talsh_complex.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
/** ExaTensor::TAL-SH: Complex arithmetic header.
REVISION: 2019/01/04
Copyright (C) 2014-2019 Dmitry I. Lyakh (Liakh)
Copyright (C) 2014-2019 Oak Ridge National Laboratory (UT-Battelle)
This file is part of ExaTensor.
ExaTensor is free software: you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as published
by the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
ExaTensor is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License
along with ExaTensor. If not, see <http://www.gnu.org/licenses/>.
------------------------------------------------------------------------
**/
#ifndef TALSH_COMPLEX_H_
#define TALSH_COMPLEX_H_
#include <math.h>
#ifdef __cplusplus
#include <complex>
#endif
#ifndef NO_GPU
#include <cuComplex.h>
#endif
//DECLARATIONS:
// Complex number:
#ifndef NO_GPU
typedef cuFloatComplex talshComplex4;
typedef cuDoubleComplex talshComplex8;
#else
#ifdef __cplusplus
typedef std::complex<float> talshComplex4;
typedef std::complex<double> talshComplex8;
#else
typedef struct{float real; float imag;} talshComplex4;
typedef struct{double real; double imag;} talshComplex8;
#endif
#endif /*NO_GPU*/
/* TAL-SH complex arithmetic headers:
inline talshComplex4 talshComplex4Set(float real, float imag);
inline talshComplex8 talshComplex8Set(double real, double imag);
inline float talshComplex4Real(talshComplex4 cmplx);
inline double talshComplex8Real(talshComplex8 cmplx);
inline float talshComplex4Imag(talshComplex4 cmplx);
inline double talshComplex8Imag(talshComplex8 cmplx);
inline talshComplex4 talshComplex4Conjg(talshComplex4 cmplx);
inline talshComplex8 talshComplex8Conjg(talshComplex8 cmplx);
inline float talshComplex4Abs(talshComplex4 cmplx);
inline double talshComplex8Abs(talshComplex8 cmplx);
inline float talshComplex4Asq(talshComplex4 cmplx);
inline double talshComplex8Asq(talshComplex8 cmplx);
inline talshComplex4 talshComplex4Add(talshComplex4 x, talshComplex4 y);
inline talshComplex8 talshComplex8Add(talshComplex8 x, talshComplex8 y);
inline void talshComplex4AddEq(talshComplex4 * x, talshComplex4 y);
inline void talshComplex8AddEq(talshComplex8 * x, talshComplex8 y);
inline talshComplex4 talshComplex4Sub(talshComplex4 x, talshComplex4 y);
inline talshComplex8 talshComplex8Sub(talshComplex8 x, talshComplex8 y);
inline void talshComplex4SubEq(talshComplex4 * x, talshComplex4 y);
inline void talshComplex8SubEq(talshComplex8 * x, talshComplex8 y);
inline talshComplex4 talshComplex4Mul(talshComplex4 x, talshComplex4 y);
inline talshComplex8 talshComplex8Mul(talshComplex8 x, talshComplex8 y);
inline void talshComplex4MulEq(talshComplex4 * x, talshComplex4 y);
inline void talshComplex8MulEq(talshComplex8 * x, talshComplex8 y);
inline talshComplex4 talshComplex4Div(talshComplex4 x, talshComplex4 y);
inline talshComplex8 talshComplex8Div(talshComplex8 x, talshComplex8 y);
inline void talshComplex4DivEq(talshComplex4 * x, talshComplex4 y);
inline void talshComplex8DivEq(talshComplex8 * x, talshComplex8 y);
*/
//DEFINITIONS:
//Construct a complex number:
#ifndef NO_GPU
__host__ __device__ __forceinline__ talshComplex4 talshComplex4Set(float real, float imag)
{
return make_cuFloatComplex(real,imag);
}
__host__ __device__ __forceinline__ talshComplex8 talshComplex8Set(double real, double imag)
{
return make_cuDoubleComplex(real,imag);
}
#else
#ifdef __cplusplus
inline talshComplex4 talshComplex4Set(float real, float imag)
{
return talshComplex4(real,imag);
}
inline talshComplex8 talshComplex8Set(double real, double imag)
{
return talshComplex8(real,imag);
}
#else
talshComplex4 talshComplex4Set(float real, float imag)
{
talshComplex4 result = {real,imag};
return result;
}
talshComplex8 talshComplex8Set(double real, double imag)
{
talshComplex8 result = {real,imag};
return result;
}
#endif
#endif
//Get the real component of a complex number:
#ifndef NO_GPU
__host__ __device__ __forceinline__ float talshComplex4Real(talshComplex4 cmplx)
{
return cuCrealf(cmplx);
}
__host__ __device__ __forceinline__ double talshComplex8Real(talshComplex8 cmplx)
{
return cuCreal(cmplx);
}
#else
#ifdef __cplusplus
inline float talshComplex4Real(talshComplex4 cmplx)
{
return cmplx.real();
}
inline double talshComplex8Real(talshComplex8 cmplx)
{
return cmplx.real();
}
#else
float talshComplex4Real(talshComplex4 cmplx)
{
return cmplx.real;
}
double talshComplex8Real(talshComplex8 cmplx)
{
return cmplx.real;
}
#endif
#endif
//Get the imaginary component of a complex number:
#ifndef NO_GPU
__host__ __device__ __forceinline__ float talshComplex4Imag(talshComplex4 cmplx)
{
return cuCimagf(cmplx);
}
__host__ __device__ __forceinline__ double talshComplex8Imag(talshComplex8 cmplx)
{
return cuCimag(cmplx);
}
#else
#ifdef __cplusplus
inline float talshComplex4Imag(talshComplex4 cmplx)
{
return cmplx.imag();
}
inline double talshComplex8Imag(talshComplex8 cmplx)
{
return cmplx.imag();
}
#else
float talshComplex4Imag(talshComplex4 cmplx)
{
return cmplx.imag;
}
double talshComplex8Imag(talshComplex8 cmplx)
{
return cmplx.imag;
}
#endif
#endif
//Get the complex conjugate:
#ifndef NO_GPU
__host__ __device__ __forceinline__ talshComplex4 talshComplex4Conjg(talshComplex4 cmplx)
{
return cuConjf(cmplx);
}
__host__ __device__ __forceinline__ talshComplex8 talshComplex8Conjg(talshComplex8 cmplx)
{
return cuConj(cmplx);
}
#else
#ifdef __cplusplus
inline talshComplex4 talshComplex4Conjg(talshComplex4 cmplx)
{
return std::conj(cmplx);
}
inline talshComplex8 talshComplex8Conjg(talshComplex8 cmplx)
{
return std::conj(cmplx);
}
#else
talshComplex4 talshComplex4Conjg(talshComplex4 cmplx)
{
talshComplex4 result = {cmplx.real,-cmplx.imag};
return result;
}
talshComplex8 talshComplex8Conjg(talshComplex8 cmplx)
{
talshComplex8 result = {cmplx.real,-cmplx.imag};
return result;
}
#endif
#endif
//Get the absolute magnitude of a complex number:
#ifndef NO_GPU
__host__ __device__ __forceinline__ float talshComplex4Abs(talshComplex4 cmplx)
{
return cuCabsf(cmplx);
}
__host__ __device__ __forceinline__ double talshComplex8Abs(talshComplex8 cmplx)
{
return cuCabs(cmplx);
}
#else
#ifdef __cplusplus
inline float talshComplex4Abs(talshComplex4 cmplx)
{
return std::abs(cmplx);
}
inline double talshComplex8Abs(talshComplex8 cmplx)
{
return std::abs(cmplx);
}
#else
float talshComplex4Abs(talshComplex4 cmplx)
{
return (float)sqrt((double)((cmplx.real)*(cmplx.real)) + (double)((cmplx.imag)*(cmplx.imag)));
}
double talshComplex8Abs(talshComplex8 cmplx)
{
return sqrt(((cmplx.real)*(cmplx.real)) + ((cmplx.imag)*(cmplx.imag)));
}
#endif
#endif
//Get the squared magnitude of a complex number:
#ifndef NO_GPU
__host__ __device__ __forceinline__ float talshComplex4Asq(talshComplex4 cmplx)
{
float rl = cuCrealf(cmplx); float im = cuCimagf(cmplx);
return (rl*rl + im*im);
}
__host__ __device__ __forceinline__ double talshComplex8Asq(talshComplex8 cmplx)
{
double rl = cuCreal(cmplx); double im = cuCimag(cmplx);
return (rl*rl + im*im);
}
#else
#ifdef __cplusplus
inline float talshComplex4Asq(talshComplex4 cmplx)
{
float rl = cmplx.real(); float im = cmplx.imag();
return (rl*rl + im*im);
}
inline double talshComplex8Asq(talshComplex8 cmplx)
{
double rl = cmplx.real(); double im = cmplx.imag();
return (rl*rl + im*im);
}
#else
float talshComplex4Asq(talshComplex4 cmplx)
{
return ((cmplx.real)*(cmplx.real) + (cmplx.imag)*(cmplx.imag));
}
double talshComplex8Asq(talshComplex8 cmplx)
{
return ((cmplx.real)*(cmplx.real) + (cmplx.imag)*(cmplx.imag));
}
#endif
#endif
//Add two complex numbers:
#ifndef NO_GPU
__host__ __device__ __forceinline__ talshComplex4 talshComplex4Add(talshComplex4 x, talshComplex4 y)
{
return cuCaddf(x,y);
}
__host__ __device__ __forceinline__ talshComplex8 talshComplex8Add(talshComplex8 x, talshComplex8 y)
{
return cuCadd(x,y);
}
#else
#ifdef __cplusplus
inline talshComplex4 talshComplex4Add(talshComplex4 x, talshComplex4 y)
{
return x+y;
}
inline talshComplex8 talshComplex8Add(talshComplex8 x, talshComplex8 y)
{
return x+y;
}
#else
talshComplex4 talshComplex4Add(talshComplex4 x, talshComplex4 y)
{
return talshComplex4Set(x.real+y.real,x.imag+y.imag);
}
talshComplex8 talshComplex8Add(talshComplex8 x, talshComplex8 y)
{
return talshComplex8Set(x.real+y.real,x.imag+y.imag);
}
#endif
#endif
//Add two complex numbers in-place:
#ifndef NO_GPU
__host__ __device__ __forceinline__ void talshComplex4AddEq(talshComplex4 * x, talshComplex4 y)
{
*x = cuCaddf(*x,y);
return;
}
__host__ __device__ __forceinline__ void talshComplex8AddEq(talshComplex8 * x, talshComplex8 y)
{
*x = cuCadd(*x,y);
return;
}
#else
#ifdef __cplusplus
inline void talshComplex4AddEq(talshComplex4 * x, talshComplex4 y)
{
*x = *x + y;
return;
}
inline void talshComplex8AddEq(talshComplex8 * x, talshComplex8 y)
{
*x = *x + y;
return;
}
#else
void talshComplex4AddEq(talshComplex4 * x, talshComplex4 y)
{
*x = talshComplex4Set(x->real+y.real,x->imag+y.imag);
return;
}
void talshComplex8AddEq(talshComplex8 * x, talshComplex8 y)
{
*x = talshComplex8Set(x->real+y.real,x->imag+y.imag);
return;
}
#endif
#endif
//Subtract two complex numbers:
#ifndef NO_GPU
__host__ __device__ __forceinline__ talshComplex4 talshComplex4Sub(talshComplex4 x, talshComplex4 y)
{
return cuCsubf(x,y);
}
__host__ __device__ __forceinline__ talshComplex8 talshComplex8Sub(talshComplex8 x, talshComplex8 y)
{
return cuCsub(x,y);
}
#else
#ifdef __cplusplus
inline talshComplex4 talshComplex4Sub(talshComplex4 x, talshComplex4 y)
{
return x-y;
}
inline talshComplex8 talshComplex8Sub(talshComplex8 x, talshComplex8 y)
{
return x-y;
}
#else
talshComplex4 talshComplex4Sub(talshComplex4 x, talshComplex4 y)
{
return talshComplex4Set(x.real-y.real,x.imag-y.imag);
}
talshComplex8 talshComplex8Sub(talshComplex8 x, talshComplex8 y)
{
return talshComplex8Set(x.real-y.real,x.imag-y.imag);
}
#endif
#endif
//Multiply two complex numbers:
#ifndef NO_GPU
__host__ __device__ __forceinline__ talshComplex4 talshComplex4Mul(talshComplex4 x, talshComplex4 y)
{
return cuCmulf(x,y);
}
__host__ __device__ __forceinline__ talshComplex8 talshComplex8Mul(talshComplex8 x, talshComplex8 y)
{
return cuCmul(x,y);
}
#else
#ifdef __cplusplus
inline talshComplex4 talshComplex4Mul(talshComplex4 x, talshComplex4 y)
{
return x*y;
}
inline talshComplex8 talshComplex8Mul(talshComplex8 x, talshComplex8 y)
{
return x*y;
}
#else
talshComplex4 talshComplex4Mul(talshComplex4 x, talshComplex4 y)
{
float rlx = x.real; float imx = x.imag;
float rly = y.real; float imy = y.imag;
return talshComplex4Set(rlx*rly-imx*imy,rlx*imy+imx*rly);
}
talshComplex8 talshComplex8Mul(talshComplex8 x, talshComplex8 y)
{
double rlx = x.real; double imx = x.imag;
double rly = y.real; double imy = y.imag;
return talshComplex8Set(rlx*rly-imx*imy,rlx*imy+imx*rly);
}
#endif
#endif
//Divide two complex numbers:
#ifndef NO_GPU
__host__ __device__ __forceinline__ talshComplex4 talshComplex4Div(talshComplex4 x, talshComplex4 y)
{
return cuCdivf(x,y);
}
__host__ __device__ __forceinline__ talshComplex8 talshComplex8Div(talshComplex8 x, talshComplex8 y)
{
return cuCdiv(x,y);
}
#else
#ifdef __cplusplus
inline talshComplex4 talshComplex4Div(talshComplex4 x, talshComplex4 y)
{
return x/y;
}
inline talshComplex8 talshComplex8Div(talshComplex8 x, talshComplex8 y)
{
return x/y;
}
#else
talshComplex4 talshComplex4Div(talshComplex4 x, talshComplex4 y)
{
float rlx = x.real; float imx = x.imag;
float rly = y.real; float imy = y.imag;
float dny = 1.0f/(rly*rly + imy*imy);
return talshComplex4Set((rlx*rly+imx*imy)*dny,(imx*rly-rlx*imy)*dny);
}
talshComplex8 talshComplex8Div(talshComplex8 x, talshComplex8 y)
{
double rlx = x.real; double imx = x.imag;
double rly = y.real; double imy = y.imag;
double dny = 1.0/(rly*rly + imy*imy);
return talshComplex8Set((rlx*rly+imx*imy)*dny,(imx*rly-rlx*imy)*dny);
}
#endif
#endif
#endif /*TALSH_COMPLEX_H_*/