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 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531
|
<?xml version="1.0" encoding="UTF-8"?>
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.1 plus MathML 2.0//EN" "http://www.w3.org/Math/DTD/mathml2/xhtml-math11-f.dtd">
<html xmlns="http://www.w3.org/1999/xhtml" lang="en-US">
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8"/>
<meta name="viewport" content="width=device-width, initial-scale=0.4"/>
<meta name="google" content="notranslate" />
<link rel="canonical" href="https://sleef.org/dft.xhtml" />
<link rel="icon" href="favicon.png" />
<link rel="stylesheet" type="text/css" href="texlike.css"/>
<link rel="stylesheet" type="text/css" href="//fonts.googleapis.com/css?family=Ubuntu" />
<link rel="stylesheet" type="text/css" href="sleef.css"/>
<title>SLEEF - DFT library reference</title>
</head>
<body translate="no" class="notranslate">
<h1>SLEEF - DFT library reference</h1>
<h2>Table of contents</h2>
<ul class="none" style="font-family: arial, sansserif; padding-left: 0.5cm;">
<li><a class="underlined" href="index.xhtml">Introduction</a></li>
<li><a class="underlined" href="compile.xhtml">Compiling and installing the library</a></li>
<li><a class="underlined" href="purec.xhtml">Math library reference</a></li>
<li> </li>
<li><a class="underlined" href="dft.xhtml">DFT library reference</a></li>
<ul class="disc">
<li><a href="dft.xhtml#tutorial">Tutorial</a></li>
<!--<li><a href="dft.xhtml#compatibility">Compatibility with other libraries</a></li>-->
<li><a href="dft.xhtml#reference">Function reference</a></li>
</ul>
<li> </li>
<li><a class="underlined" href="misc.xhtml">Other tools included in the package</a></li>
<li><a class="underlined" href="benchmark.xhtml">Benchmark results</a></li>
<li><a class="underlined" href="additional.xhtml">Additional notes</a></li>
</ul>
<h2 id="tutorial">Tutorial</h2>
<p class="noindent">
In this section, how to use the DFT library is explained with an
example source code shown below.
<a class="underlined" href="tutorial.c">This source code</a> is
included in the distribution package under src/dft-tester directory.
</p>
<pre class="code">
<code>// gcc tutorial.c -lsleef -lsleefdft -lm -lstdc++</code>
<code>#include <stdio.h></code>
<code>#include <stdlib.h></code>
<code>#include <stdint.h></code>
<code>#include <math.h></code>
<code>#include <complex.h></code>
<code></code>
<code>#include "sleef.h"</code>
<code>#include "sleefdft.h"</code>
<code></code>
<code>#define THRES 1e-4</code>
<code></code>
<code>typedef double complex cmpl;</code>
<code></code>
<code>cmpl omega(double n, double kn) {</code>
<code> return cexp((-2 * M_PI * _Complex_I / n) * kn);</code>
<code>}</code>
<code></code>
<code>void forward(cmpl *ts, cmpl *fs, int len) {</code>
<code> for(int k=0;k<len;k++) {</code>
<code> fs[k] = 0;</code>
<code> for(int n=0;n<len;n++) fs[k] += ts[n] * omega(len, n*k);</code>
<code> }</code>
<code>}</code>
<code></code>
<code>int main(int argc, char **argv) {</code>
<code> int n = 256;</code>
<code> if (argc == 2) n = 1 << atoi(argv[1]);</code>
<code></code>
<code> SleefDFT_setPlanFilePath("tut.plan", NULL, SLEEF_PLAN_AUTOMATIC);</code>
<code></code>
<code> double *sx = (double *)Sleef_malloc(n*2 * sizeof(double));</code>
<code> double *sy = (double *)Sleef_malloc(n*2 * sizeof(double));</code>
<code></code>
<code> struct SleefDFT *p = SleefDFT_double_init1d(n, sx, sy, SLEEF_MODE_FORWARD);</code>
<code></code>
<code> if (p == NULL) {</code>
<code> printf("SleefDFT initialization failed\n");</code>
<code> exit(-1);</code>
<code> }</code>
<code></code>
<code> cmpl *ts = (cmpl *)malloc(sizeof(cmpl)*n);</code>
<code> cmpl *fs = (cmpl *)malloc(sizeof(cmpl)*n);</code>
<code></code>
<code> for(int i=0;i<n;i++) {</code>
<code> ts[i] =</code>
<code> (2.0 * (rand() / (double)RAND_MAX) - 1) * 1.0 +</code>
<code> (2.0 * (rand() / (double)RAND_MAX) - 1) * _Complex_I;</code>
<code></code>
<code> sx[(i*2+0)] = creal(ts[i]);</code>
<code> sx[(i*2+1)] = cimag(ts[i]);</code>
<code> }</code>
<code></code>
<code> forward(ts, fs, n);</code>
<code></code>
<code> SleefDFT_double_execute(p, NULL, NULL);</code>
<code></code>
<code> int success = 1;</code>
<code></code>
<code> for(int i=0;i<n;i++) {</code>
<code> if ((fabs(sy[(i*2+0)] - creal(fs[i])) > THRES) ||</code>
<code> (fabs(sy[(i*2+1)] - cimag(fs[i])) > THRES)) {</code>
<code> success = 0;</code>
<code> }</code>
<code> }</code>
<code></code>
<code> printf("%s\n", success ? "OK" : "NG");</code>
<code></code>
<code> free(fs); free(ts);</code>
<code> Sleef_free(sy); Sleef_free(sx);</code>
<code></code>
<code> SleefDFT_dispose(p);</code>
<code></code>
<code> exit(success);</code>
<code>}</code>
</pre>
<p style="text-align:center;">
Fig. 4.1: <a href="tutorial.c">Test code for DFT subroutines</a>
</p>
<div style="margin-top: 1.0cm;"></div>
<p>
You can compile the source code with the following command. Note
that you need to link with libstdc++.
</p>
<pre class="command">$ gcc tutorial.c -lsleef -lsleefdft -lm -lstdc++
</pre>
<p>
This program takes one integer argument <i class="var">n</i>, and
executes forward complex transform with size
2<sup><i class="var">n</i></sup>. It then compares the results of a
naive transform and the transform performed by the library. If the
two results match, it prints OK.
</p>
<p>
For the first execution, this program takes a few seconds to
finish. This is because the library constructs an execution plan by
measuring computation time with many different configurations. The
best plan is saved to "tut.plan", which is specified at line
30. Later executions will finish instantly as the library reads the
plan from this file. Instead of specifying the file name in the
program, the file can be specified by SLEEFDFTPLAN environment
variable. Instead of constructing a plan, the library can estimate a
modestly good configuration, if you specify SLEEF_MODE_ESTIMATE flag
at line 30.
</p>
<p>
This library executes transforms using the most suitable SIMD
instructions available on the computer and utilizing
multi-threading. The library requires the input and output arrays to
be aligned to certain boundaries so that the data in the arrays can
be accessed efficiently by SIMD instructions. You can allocate arrays with
<b class="func">Sleef_malloc</b> to make them aligned to the
boundaries, as seen in line 32 and 33. Arrays allocated
with <b class="func">Sleef_malloc</b> have to be freed
with <b class="func">Sleef_free</b>, as seen in line 70. Instead of
allocating arrays with <b class="func">Sleef_malloc</b>, you can
allocate an aligned memory region yourself, and pass the pointer to
the library.
</p>
<p>
The real and imaginary parts of the <i>k</i>-th number
are stored in (2<i>k</i>)-th and
(2<i>k</i>+1)-th elements of the input and output array,
respectively. At line 56, the transform is executed by the
library. You can specify the same array as the input and output.
</p>
<p>
Under src/dft-tester directory, there are other examples showing how
to execute transforms in a way that you get equivalent results to
other libraries. The API of this library is designed to be easy to
migrate from FFTW, and you do not need to change the contents of the
input or output arrays upon migration.
</p>
<h2 id="reference">Function reference</h2>
<p class="funcname"><b class="func">Sleef_malloc</b> - allocate aligned memory</p>
<p class="header">Synopsis</p>
<p class="synopsis">
#include <stdlib.h><br/>
#include <sleef.h><br/>
<br/>
<b class="type">void *</b> <b class="func">Sleef_malloc</b>(<b class="type">size_t</b> <i class="var">z</i>);<br/>
<br/>
Link with -lsleef.
</p>
<p class="header">Description</p>
<p class="noindent">
<b class="func">Sleef_malloc</b> allocates <i class="var">z</i>
bytes of aligned memory region, and return the pointer to the
allocated region. The returned pointer points the address of a
memory region that can be effciently accessed by all SIMD load and
store instructions available on that computer. Memory regions
allocated by <b class="func">Sleef_malloc</b> have to be freed
with <b class="func">Sleef_free</b>.
</p>
<hr/>
<p class="funcname"><b class="func">Sleef_free</b> - free memory allocated by <b class="func">Sleef_malloc</b></p>
<p class="header">Synopsis</p>
<p class="synopsis">
#include <stdlib.h><br/>
#include <sleef.h><br/>
<br/>
<b class="type">void</b> <b class="func">Sleef_free</b>(<b class="type">void *</b><i class="var">ptr</i>);<br/>
<br/>
Link with -lsleef.
</p>
<p class="header">Description</p>
<p class="noindent">
A memory region pointed by <i class="var">ptr</i> that is allocated
by <b class="func">Sleef_malloc</b> can be freed
with <b class="func">Sleef_free</b>.
</p>
<hr/>
<p class="funcname"><b class="func">SleefDFT_setPlanFilePath</b> - set the file path for storing execution plans</p>
<p class="header">Synopsis</p>
<p class="synopsis">
#include <stdint.h><br/>
#include <sleefdft.h><br/>
<br/>
<b class="type">void</b> <b class="func">SleefDFT_setPlanFilePath</b>(<b class="type">const char *</b><i class="var">path</i>, <b class="type">const char *</b><i class="var">arch</i>, <b class="type">uint64_t</b> <i class="var">mode</i>);<br/>
<br/>
Link with -lsleefdft -lsleef.
</p>
<p class="header">Description</p>
<p class="noindent">
File name for storing execution plan can be specified by this
function. If NULL is specified as <i class="var">path</i>, the file
name is read from SLEEFDFTPLAN environment variable. A string for
identifying system micro architecture can be also given. The library
will try to automatically detect the micro architecture if NULL is
given as <i class="var">arch</i>. Management options for the plan file
can be specified by the <i class="var">mode</i> parameter, as shown
below.
</p>
<div style="margin-top: 1.0cm;"></div>
<table style="text-align:center;" align="center">
<tr align="center">
<td class="caption">Table 4.2: Mode flags for SleefFT_setPlanFilePath</td>
</tr>
<tr align="center">
<td>
<table class="lt">
<tr>
<td class="lt-hl"></td>
<td class="lt-hl"></td>
</tr>
<tr>
<td class="lt-br" align="center">Flag</td>
<td class="lt-b" align="center">Meaning</td>
</tr>
<tr>
<td class="lt-hl"></td>
<td class="lt-hl"></td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_PLAN_AUTOMATIC</td>
<td class="lt-" align="left">Execution plans are automatically loaded and saved. Plans are generated if it does not exist.</td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_PLAN_READONLY</td>
<td class="lt-" align="left">Execution plans are automatically loaded, but not saved.</td>
</tr>
<tr>
<td class="lt-br" align="left">SLEEF_PLAN_RESET</td>
<td class="lt-b" align="left">Existing execution plans are reset and constructed from the beginning.</td>
</tr>
</table>
</td>
</tr>
</table>
<hr/>
<p class="funcname"><b class="func">SleefDFT_savePlan</b> - save plan to the specified file</p>
<p class="header">Synopsis</p>
<p class="synopsis">
#include <sleef.h><br/>
<br/>
<b class="type">void</b> <b class="func">SleefDFT_savePlan</b>(<b class="type">const char *</b><i class="var">path</i>);<br/>
<br/>
Link with -lsleefdft -lsleef.
</p>
<p class="header">Description</p>
<p class="noindent">
By calling this function, the generated plan is saved to the file
specified by <i class="var">path</i>. You can use this function to
manually save the plan when SLEEF_PLAN_AUTOMATIC is not specified as
the mode with <b class="func">SleefDFT_setPlanFilePath</b>.
</p>
<hr/>
<p class="funcname"><b class="func">SleefDFT_double_init1d</b>, <b class="func">SleefDFT_float_init1d</b> - initialize the tables for 1D transform</p>
<p class="header">Synopsis</p>
<p class="synopsis">
#include <stdint.h><br/>
#include <sleefdft.h><br/>
<br/>
<b class="type">struct SleefDFT *</b> <b class="func">SleefDFT_double_init1d</b>(<b class="type">uint32_t</b> <i class="var">n</i>, <b class="type">const double *</b><i class="var">in</i>, <b class="type">double *</b><i class="var">out</i>, <b class="type">uint64_t</b> <i class="var">mode</i>);<br/>
<b class="type">struct SleefDFT *</b> <b class="func">SleefDFT_float_init1d</b>(<b class="type">uint32_t</b> <i class="var">n</i>, <b class="type">const float *</b><i class="var">in</i>, <b class="type">float *</b><i class="var">out</i>, <b class="type">uint64_t</b> <i class="var">mode</i>);<br/>
<br/>
Link with -lsleefdft -lsleef.
</p>
<p class="header">Description</p>
<p class="noindent">
These functions generate and initialize the tables that are used for
1D transform, and return the pointer. The size of transform can be
specified by <i class="var">n</i>. Currently, power-of-two sizes can
be only specified. The list of the flags that can be passed
to <i class="var">mode</i> is shown below.
</p>
<div style="margin-top: 1.0cm;"></div>
<table style="text-align:center;" align="center">
<tr align="center">
<td class="caption">Table 4.3: Mode flags for SleefDFT_double_init</td>
</tr>
<tr align="center">
<td>
<table class="lt">
<tr>
<td class="lt-hl"></td>
<td class="lt-hl"></td>
</tr>
<tr>
<td class="lt-br" align="center">Flag</td>
<td class="lt-b" align="center">Meaning</td>
</tr>
<tr>
<td class="lt-hl"></td>
<td class="lt-hl"></td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_MODE_FORWARD</td>
<td class="lt-" align="left">Tables are initialized for forward transforms.</td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_MODE_BACKWARD</td>
<td class="lt-" align="left">Tables are initialized for backward transforms.</td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_MODE_COMPLEX</td>
<td class="lt-" align="left">Tables are initialized for complex transforms.</td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_MODE_REAL</td>
<td class="lt-" align="left">Tables are initialized for real transforms.</td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_MODE_ALT</td>
<td class="lt-" align="left">Tables are initialized for alternative real transforms.</td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_MODE_ESTIMATE</td>
<td class="lt-" align="left">Execution plans are estimated.</td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_MODE_MEASURE</td>
<td class="lt-" align="left">Execution plans are measured when they are needed.</td>
</tr>
<tr>
<td class="lt-r" align="left">SLEEF_MODE_VERBOSE</td>
<td class="lt-" align="left">Messages are displayed.</td>
</tr>
<tr>
<td class="lt-br" align="left">SLEEF_MODE_NO_MT</td>
<td class="lt-b" align="left">Multithreading will be disabled in the computation for transforms.</td>
</tr>
</table>
</td>
</tr>
</table>
<p class="header">Return value</p>
<p class="noindent">
These functions return a pointer to the data that is used for 1D DFT
computation, or NULL if an error occurred.
</p>
<hr/>
<p class="funcname"><b class="func">SleefDFT_double_init2d</b>, <b class="func">SleefDFT_float_init2d</b> - initialize the tables for 2D transform</p>
<p class="header">Synopsis</p>
<p class="synopsis">
#include <stdint.h><br/>
#include <sleefdft.h><br/>
<br/>
<b class="type">struct SleefDFT *</b> <b class="func">SleefDFT_double_init2d</b>(<b class="type">uint32_t</b> <i class="var">n</i>, <b class="type">uint32_t</b> <i class="var">m</i>, <b class="type">const double *</b><i class="var">in</i>, <b class="type">double *</b><i class="var">out</i>, <b class="type">uint64_t</b> <i class="var">mode</i>);<br/>
<b class="type">struct SleefDFT *</b> <b class="func">SleefDFT_float_init2d</b>(<b class="type">uint32_t</b> <i class="var">n</i>, <b class="type">uint32_t</b> <i class="var">m</i>, <b class="type">const float *</b><i class="var">in</i>, <b class="type">float *</b><i class="var">out</i>, <b class="type">uint64_t</b> <i class="var">mode</i>);<br/>
<br/>
Link with -lsleefdft -lsleef.
</p>
<p class="header">Description</p>
<p class="noindent">
These functions generate and initilize the tables that are used for
2D transform, and return the pointer. The size of transform can be
specified by <i class="var">n</i>
and <i class="var">m</i>. Currently, power-of-two sizes can be only
specified. The flags that can be passed to <i class="var">mode</i>
are the same as the flags listed in the section
for <b class="func">SleefDFT_double_init1d</b>.
</p>
<p class="header">Return value</p>
<p class="noindent">
These functions return a pointer to the data that is used for 2D DFT
computation, or NULL if an error occurred.
</p>
<hr/>
<p class="funcname"><b class="func">SleefDFT_double_execute</b>, <b class="func">SleefDFT_float_execute</b>, <b class="func">SleefDFT_execute</b> - execute a transform</p>
<p class="header">Synopsis</p>
<p class="synopsis">
#include <stdint.h><br/>
#include <sleefdft.h><br/>
<br/>
<b class="type">void</b> <b class="func">SleefDFT_double_execute</b>(<b class="type">struct SleefDFT *</b><i class="var">ptr</i>, <b class="type">const double *</b><i class="var">in</i>, <b class="type">double *</b><i class="var">out</i>);<br/>
<b class="type">void</b> <b class="func">SleefDFT_float_execute</b>(<b class="type">struct SleefDFT *</b><i class="var">ptr</i>, <b class="type">const float *</b><i class="var">in</i>, <b class="type">float *</b><i class="var">out</i>);<br/>
<b class="type">void</b> <b class="func">SleefDFT_execute</b>(<b class="type">struct SleefDFT *</b><i class="var">ptr</i>, <b class="type">const void *</b><i class="var">in</i>, <b class="type">void *</b><i class="var">out</i>);<br/>
<br/>
Link with -lsleefdft -lsleef.
</p>
<p class="header">Description</p>
<p class="noindent">
By calling these functions, the transform with the plan specified by
<i class="var">ptr</i> is executed. If <i class="var">in</i>
or <i class="var">out</i> is NULL, the pointer specified upon
initialization is used. Otherwise, the pointers must be pointers
returned from Sleef_malloc function. You can specify the same pointer
to <i class="var">in</i> and <i class="var">out</i>.
</p>
<hr/>
<p class="funcname"><b class="func">SleefDFT_dispose</b> - dispose the tables for transforms</p>
<p class="header">Synopsis</p>
<p class="synopsis">
#include <stdint.h><br/>
#include <sleefdft.h><br/>
<br/>
<b class="type">void</b> <b class="func">SleefDFT_dispose</b>(<b class="type">struct SleefDFT *</b><i class="var">ptr</i>);<br/>
<br/>
Link with -lsleefdft -lsleef.
</p>
<p class="header">Description</p>
<p class="noindent">
This function frees a plan returned
by <b class="func">SleefDFT_double_init1d</b>, <b class="func">SleefDFT_float_init1d</b>, <b class="func">SleefDFT_double_init2d</b>, or <b class="func">SleefDFT_float_init2d</b> functions.
</p>
<p class="footer">
Copyright © 2010-2025 SLEEF Project, Naoki Shibata and contributors.<br/>
SLEEF is open-source software and is distributed under the Boost Software License, Version 1.0.
</p>
<script type="text/javascript">
var sc_project=13098265;
var sc_invisible=1;
var sc_security="518de45e";
</script>
<script type="text/javascript"
src="https://www.statcounter.com/counter/counter.js"
async="async"></script>
<noscript><div class="statcounter"><a title="Web Analytics"
href="https://statcounter.com/" target="_blank"><img
class="statcounter"
src="https://c.statcounter.com/13098265/0/518de45e/1/"
alt="Web Analytics"
referrerPolicy="no-referrer-when-downgrade"></img></a></div></noscript>
</body>
</html>
|