Report a bug
If you spot a problem with this page, click here to create a Bugzilla issue.
Improve this page
Quickly fork, edit online, and submit a pull request for this page. Requires a signed-in GitHub account. This works well for small changes. If you'd like to make larger changes you may want to consider using a local clone.

std.numeric

This module is a port of a growing fragment of the numeric header in Alexander Stepanov's Standard Template Library, with a few additions.
Authors:
Andrei Alexandrescu, Don Clugston, Robert Jacques, Ilya Yaroshenko
enum CustomFloatFlags: int;
Format flags for CustomFloat.
signed
Adds a sign bit to allow for signed numbers.
storeNormalized
Store values in normalized form by default. The actual precision of the significand is extended by 1 bit by assuming an implicit leading bit of 1 instead of 0. i.e. 1.nnnn instead of 0.nnnn. True for all IEE754 types
allowDenorm
Stores the significand in IEEE754 denormalized form when the exponent is 0. Required to express the value 0.
infinity
Allows the storage of IEEE754 infinity values.
nan
Allows the storage of IEEE754 Not a Number values.
probability
If set, select an exponent bias such that max_exp = 1. i.e. so that the maximum value is >= 1.0 and < 2.0. Ignored if the exponent bias is manually specified.
negativeUnsigned
If set, unsigned custom floats are assumed to be negative.
allowDenormZeroOnly
If set, 0 is the only allowed IEEE754 denormalized number. Requires allowDenorm and storeNormalized.
ieee
Include all of the IEEE754 options.
none
Include none of the above options.
template CustomFloat(uint bits) if (bits == 8 || bits == 16 || bits == 32 || bits == 64 || bits == 80)

template CustomFloat(uint precision, uint exponentWidth, CustomFloatFlags flags = CustomFloatFlags.ieee) if (((flags & flags.signed) + precision + exponentWidth) % 8 == 0 && precision + exponentWidth > 0)

struct CustomFloat(uint precision, uint exponentWidth, CustomFloatFlags flags, uint bias) if (((flags & flags.signed) + precision + exponentWidth) % 8 == 0 && precision + exponentWidth > 0);
Allows user code to define custom floating-point formats. These formats are for storage only; all operations on them are performed by first implicitly extracting them to real first. After the operation is completed the result can be stored in a custom floating-point value via assignment.
Examples:
import std.math : sin, cos;

// Define a 16-bit floating point values
CustomFloat!16                                x;     // Using the number of bits
CustomFloat!(10, 5)                           y;     // Using the precision and exponent width
CustomFloat!(10, 5,CustomFloatFlags.ieee)     z;     // Using the precision, exponent width and format flags
CustomFloat!(10, 5,CustomFloatFlags.ieee, 15) w;     // Using the precision, exponent width, format flags and exponent offset bias

// Use the 16-bit floats mostly like normal numbers
w = x*y - 1;

// Functions calls require conversion
z = sin(+x)           + cos(+y);                     // Use unary plus to concisely convert to a real
z = sin(x.get!float)  + cos(y.get!float);            // Or use get!T
z = sin(cast(float) x) + cos(cast(float) y);           // Or use cast(T) to explicitly convert

// Define a 8-bit custom float for storing probabilities
alias Probability = CustomFloat!(4, 4, CustomFloatFlags.ieee^CustomFloatFlags.probability^CustomFloatFlags.signed );
auto p = Probability(0.5);
template FPTemporary(F) if (isFloatingPoint!F)
Defines the fastest type to use when storing temporaries of a calculation intended to ultimately yield a result of type F (where F must be one of float, double, or real). When doing a multi-step computation, you may want to store intermediate results as FPTemporary!F.
The necessity of FPTemporary stems from the optimized floating-point operations and registers present in virtually all processors. When adding numbers in the example above, the addition may in fact be done in real precision internally. In that case, storing the intermediate result in double format is not only less precise, it is also (surprisingly) slower, because a conversion from real to double is performed every pass through the loop. This being a lose-lose situation, FPTemporary!F has been defined as the fastest type to use for calculations at precision F. There is no need to define a type for the most accurate calculations, as that is always real.
Finally, there is no guarantee that using FPTemporary!F will always be fastest, as the speed of floating-point calculations depends on very many factors.
Examples:
import std.math : approxEqual;

// Average numbers in an array
double avg(in double[] a)
{
    if (a.length == 0) return 0;
    FPTemporary!double result = 0;
    foreach (e; a) result += e;
    return result / a.length;
}

auto a = [1.0, 2.0, 3.0];
assert(approxEqual(avg(a), 2));
template secantMethod(alias fun)
Implements the secant method for finding a root of the function fun starting from points [xn_1, x_n] (ideally close to the root). Num may be float, double, or real.
Examples:
import std.math : approxEqual, cos;

float f(float x)
{
    return cos(x) - x*x*x;
}
auto x = secantMethod!(f)(0f, 1f);
assert(approxEqual(x, 0.865474));
T findRoot(T, DF, DT)(scope DF f, in T a, in T b, scope DT tolerance)
if (isFloatingPoint!T && is(typeof(tolerance(T.init, T.init)) : bool) && is(typeof(f(T.init)) == R, R) && isFloatingPoint!R);

T findRoot(T, DF)(scope DF f, in T a, in T b);
Find a real root of a real function f(x) via bracketing.
Given a function f and a range [a .. b] such that f(a) and f(b) have opposite signs or at least one of them equals ±0, returns the value of x in the range which is closest to a root of f(x). If f(x) has more than one root in the range, one will be chosen arbitrarily. If f(x) returns NaN, NaN will be returned; otherwise, this algorithm is guaranteed to succeed.
Uses an algorithm based on TOMS748, which uses inverse cubic interpolation whenever possible, otherwise reverting to parabolic or secant interpolation. Compared to TOMS748, this implementation improves worst-case performance by a factor of more than 100, and typical performance by a factor of 2. For 80-bit reals, most problems require 8 to 15 calls to f(x) to achieve full machine precision. The worst-case performance (pathological cases) is approximately twice the number of bits.

References "On Enclosing Simple Roots of Nonlinear Equations", G. Alefeld, F.A. Potra, Yixun Shi, Mathematics of Computation 61, pp733-744 (1993). Fortran code available from www.netlib.org as algorithm TOMS478.

Tuple!(T, T, R, R) findRoot(T, R, DF, DT)(scope DF f, in T ax, in T bx, in R fax, in R fbx, scope DT tolerance)
if (isFloatingPoint!T && is(typeof(tolerance(T.init, T.init)) : bool) && is(typeof(f(T.init)) == R) && isFloatingPoint!R);

Tuple!(T, T, R, R) findRoot(T, R, DF)(scope DF f, in T ax, in T bx, in R fax, in R fbx);

T findRoot(T, R)(scope R delegate(T) f, in T a, in T b, scope bool delegate(T lo, T hi) tolerance = (T a, T b) => false);
Find root of a real function f(x) by bracketing, allowing the termination condition to be specified.
Parameters:
DF f Function to be analyzed
T ax Left bound of initial range of f known to contain the root.
T bx Right bound of initial range of f known to contain the root.
R fax Value of f(ax).
R fbx Value of f(bx). fax and fbx should have opposite signs. (f(ax) and f(bx) are commonly known in advance.)
DT tolerance Defines an early termination condition. Receives the current upper and lower bounds on the root. The delegate must return true when these bounds are acceptable. If this function always returns false, full machine precision will be achieved.
Returns:
A tuple consisting of two ranges. The first two elements are the range (in x) of the root, while the second pair of elements are the corresponding function values at those points. If an exact root was found, both of the first two elements will contain the root, and the second pair of elements will be 0.
Tuple!(T, "x", Unqual!(ReturnType!DF), "y", T, "error") findLocalMin(T, DF)(scope DF f, in T ax, in T bx, in T relTolerance = sqrt(T.epsilon), in T absTolerance = sqrt(T.epsilon))
if (isFloatingPoint!T && __traits(compiles, () { T _ = DF.init(T.init); } ));
Find a real minimum of a real function f(x) via bracketing. Given a function f and a range (ax .. bx), returns the value of x in the range which is closest to a minimum of f(x). f is never evaluted at the endpoints of ax and bx. If f(x) has more than one minimum in the range, one will be chosen arbitrarily. If f(x) returns NaN or -Infinity, (x, f(x), NaN) will be returned; otherwise, this algorithm is guaranteed to succeed.
Parameters:
DF f Function to be analyzed
T ax Left bound of initial range of f known to contain the minimum.
T bx Right bound of initial range of f known to contain the minimum.
T relTolerance Relative tolerance.
T absTolerance Absolute tolerance.

Preconditions ax and bx shall be finite reals.
relTolerance shall be normal positive real.
absTolerance shall be normal positive real no less then T.epsilon*2.

Returns:
A tuple consisting of x, y = f(x) and error = 3 * (absTolerance * fabs(x) + relTolerance).
The method used is a combination of golden section search and successive parabolic interpolation. Convergence is never much slower than that for a Fibonacci search.

References "Algorithms for Minimization without Derivatives", Richard Brent, Prentice-Hall, Inc. (1973)

Examples:
import std.math : approxEqual;

auto ret = findLocalMin((double x) => (x-4)^^2, -1e7, 1e7);
assert(ret.x.approxEqual(4.0));
assert(ret.y.approxEqual(0.0));
CommonType!(ElementType!Range1, ElementType!Range2) euclideanDistance(Range1, Range2)(Range1 a, Range2 b)
if (isInputRange!Range1 && isInputRange!Range2);

CommonType!(ElementType!Range1, ElementType!Range2) euclideanDistance(Range1, Range2, F)(Range1 a, Range2 b, F limit)
if (isInputRange!Range1 && isInputRange!Range2);
Computes Euclidean distance between input ranges a and b. The two ranges must have the same length. The three-parameter version stops computation as soon as the distance is greater than or equal to limit (this is useful to save computation if a small distance is sought).
CommonType!(ElementType!Range1, ElementType!Range2) dotProduct(Range1, Range2)(Range1 a, Range2 b)
if (isInputRange!Range1 && isInputRange!Range2 && !(isArray!Range1 && isArray!Range2));

CommonType!(F1, F2) dotProduct(F1, F2)(in F1[] avector, in F2[] bvector);
Computes the dot product of input ranges a and b. The two ranges must have the same length. If both ranges define length, the check is done once; otherwise, it is done at each iteration.
CommonType!(ElementType!Range1, ElementType!Range2) cosineSimilarity(Range1, Range2)(Range1 a, Range2 b)
if (isInputRange!Range1 && isInputRange!Range2);
Computes the cosine similarity of input ranges a and b. The two ranges must have the same length. If both ranges define length, the check is done once; otherwise, it is done at each iteration. If either range has all-zero elements, return 0.
bool normalize(R)(R range, ElementType!R sum = 1)
if (isForwardRange!R);
Normalizes values in range by multiplying each element with a number chosen such that values sum up to sum. If elements in range sum to zero, assigns sum / range.length to all. Normalization makes sense only if all elements in range are positive. normalize assumes that is the case without checking it.
Returns:
true if normalization completed normally, false if all elements in range were zero or if range is empty.
Examples:
double[] a = [];
assert(!normalize(a));
a = [ 1.0, 3.0 ];
assert(normalize(a));
writeln(a); // [0.25, 0.75]
a = [ 0.0, 0.0 ];
assert(!normalize(a));
writeln(a); // [0.5, 0.5]
ElementType!Range sumOfLog2s(Range)(Range r)
if (isInputRange!Range && isFloatingPoint!(ElementType!Range));
Compute the sum of binary logarithms of the input range r. The error of this method is much smaller than with a naive sum of log2.
Examples:
import std.math : isNaN;

writeln(sumOfLog2s(new double[0])); // 0
writeln(sumOfLog2s([0.0L])); // -real.infinity
writeln(sumOfLog2s([-0.0L])); // -real.infinity
writeln(sumOfLog2s([2.0L])); // 1
assert(sumOfLog2s([-2.0L]).isNaN());
assert(sumOfLog2s([real.nan]).isNaN());
assert(sumOfLog2s([-real.nan]).isNaN());
writeln(sumOfLog2s([real.infinity])); // real.infinity
assert(sumOfLog2s([-real.infinity]).isNaN());
writeln(sumOfLog2s([0.25, 0.25, 0.25, 0.125])); // -9
ElementType!Range entropy(Range)(Range r)
if (isInputRange!Range);

ElementType!Range entropy(Range, F)(Range r, F max)
if (isInputRange!Range && !is(CommonType!(ElementType!Range, F) == void));
Computes entropy of input range r in bits. This function assumes (without checking) that the values in r are all in [0, 1]. For the entropy to be meaningful, often r should be normalized too (i.e., its values should sum to 1). The two-parameter version stops evaluating as soon as the intermediate result is greater than or equal to max.
CommonType!(ElementType!Range1, ElementType!Range2) kullbackLeiblerDivergence(Range1, Range2)(Range1 a, Range2 b)
if (isInputRange!Range1 && isInputRange!Range2);
Computes the Kullback-Leibler divergence between input ranges a and b, which is the sum ai * log(ai / bi). The base of logarithm is 2. The ranges are assumed to contain elements in [0, 1]. Usually the ranges are normalized probability distributions, but this is not required or checked by kullbackLeiblerDivergence. If any element bi is zero and the corresponding element ai nonzero, returns infinity. (Otherwise, if ai == 0 && bi == 0, the term ai * log(ai / bi) is considered zero.) If the inputs are normalized, the result is positive.
Examples:
import std.math : approxEqual;

double[] p = [ 0.0, 0, 0, 1 ];
writeln(kullbackLeiblerDivergence(p, p)); // 0
double[] p1 = [ 0.25, 0.25, 0.25, 0.25 ];
writeln(kullbackLeiblerDivergence(p1, p1)); // 0
writeln(kullbackLeiblerDivergence(p, p1)); // 2
writeln(kullbackLeiblerDivergence(p1, p)); // double.infinity
double[] p2 = [ 0.2, 0.2, 0.2, 0.4 ];
assert(approxEqual(kullbackLeiblerDivergence(p1, p2), 0.0719281));
assert(approxEqual(kullbackLeiblerDivergence(p2, p1), 0.0780719));
CommonType!(ElementType!Range1, ElementType!Range2) jensenShannonDivergence(Range1, Range2)(Range1 a, Range2 b)
if (isInputRange!Range1 && isInputRange!Range2 && is(CommonType!(ElementType!Range1, ElementType!Range2)));

CommonType!(ElementType!Range1, ElementType!Range2) jensenShannonDivergence(Range1, Range2, F)(Range1 a, Range2 b, F limit)
if (isInputRange!Range1 && isInputRange!Range2 && is(typeof(CommonType!(ElementType!Range1, ElementType!Range2).init >= F.init) : bool));
Computes the Jensen-Shannon divergence between a and b, which is the sum (ai * log(2 * ai / (ai + bi)) + bi * log(2 * bi / (ai + bi))) / 2. The base of logarithm is 2. The ranges are assumed to contain elements in [0, 1]. Usually the ranges are normalized probability distributions, but this is not required or checked by jensenShannonDivergence. If the inputs are normalized, the result is bounded within [0, 1]. The three-parameter version stops evaluations as soon as the intermediate result is greater than or equal to limit.
Examples:
import std.math : approxEqual;

double[] p = [ 0.0, 0, 0, 1 ];
writeln(jensenShannonDivergence(p, p)); // 0
double[] p1 = [ 0.25, 0.25, 0.25, 0.25 ];
writeln(jensenShannonDivergence(p1, p1)); // 0
assert(approxEqual(jensenShannonDivergence(p1, p), 0.548795));
double[] p2 = [ 0.2, 0.2, 0.2, 0.4 ];
assert(approxEqual(jensenShannonDivergence(p1, p2), 0.0186218));
assert(approxEqual(jensenShannonDivergence(p2, p1), 0.0186218));
assert(approxEqual(jensenShannonDivergence(p2, p1, 0.005), 0.00602366));
F gapWeightedSimilarity(alias comp = "a == b", R1, R2, F)(R1 s, R2 t, F lambda)
if (isRandomAccessRange!R1 && hasLength!R1 && isRandomAccessRange!R2 && hasLength!R2);
The so-called "all-lengths gap-weighted string kernel" computes a similarity measure between s and t based on all of their common subsequences of all lengths. Gapped subsequences are also included.
To understand what gapWeightedSimilarity(s, t, lambda) computes, consider first the case lambda = 1 and the strings s = ["Hello", "brave", "new", "world"] and t = ["Hello", "new", "world"]. In that case, gapWeightedSimilarity counts the following matches:
  1. three matches of length 1, namely "Hello", "new", and "world";
  2. three matches of length 2, namely ("Hello", "new"), ("Hello", "world"), and ("new", "world");
  3. one match of length 3, namely ("Hello", "new", "world").
The call gapWeightedSimilarity(s, t, 1) simply counts all of these matches and adds them up, returning 7.
string[] s = ["Hello", "brave", "new", "world"];
string[] t = ["Hello", "new", "world"];
assert(gapWeightedSimilarity(s, t, 1) == 7);
Note how the gaps in matching are simply ignored, for example ("Hello", "new") is deemed as good a match as ("new", "world"). This may be too permissive for some applications. To eliminate gapped matches entirely, use lambda = 0:
string[] s = ["Hello", "brave", "new", "world"];
string[] t = ["Hello", "new", "world"];
assert(gapWeightedSimilarity(s, t, 0) == 4);
The call above eliminated the gapped matches ("Hello", "new"), ("Hello", "world"), and ("Hello", "new", "world") from the tally. That leaves only 4 matches.
The most interesting case is when gapped matches still participate in the result, but not as strongly as ungapped matches. The result will be a smooth, fine-grained similarity measure between the input strings. This is where values of lambda between 0 and 1 enter into play: gapped matches are exponentially penalized with the number of gaps with base lambda. This means that an ungapped match adds 1 to the return value; a match with one gap in either string adds lambda to the return value; ...; a match with a total of n gaps in both strings adds pow(lambda, n) to the return value. In the example above, we have 4 matches without gaps, 2 matches with one gap, and 1 match with three gaps. The latter match is ("Hello", "world"), which has two gaps in the first string and one gap in the second string, totaling to three gaps. Summing these up we get 4 + 2 * lambda + pow(lambda, 3).
string[] s = ["Hello", "brave", "new", "world"];
string[] t = ["Hello", "new", "world"];
assert(gapWeightedSimilarity(s, t, 0.5) == 4 + 0.5 * 2 + 0.125);
gapWeightedSimilarity is useful wherever a smooth similarity measure between sequences allowing for approximate matches is needed. The examples above are given with words, but any sequences with elements comparable for equality are allowed, e.g. characters or numbers. gapWeightedSimilarity uses a highly optimized dynamic programming implementation that needs 16 * min(s.length, t.length) extra bytes of memory and Ο(s.length * t.length) time to complete.
Select!(isFloatingPoint!F, F, double) gapWeightedSimilarityNormalized(alias comp = "a == b", R1, R2, F)(R1 s, R2 t, F lambda, F sSelfSim = F.init, F tSelfSim = F.init)
if (isRandomAccessRange!R1 && hasLength!R1 && isRandomAccessRange!R2 && hasLength!R2);
The similarity per gapWeightedSimilarity has an issue in that it grows with the lengths of the two strings, even though the strings are not actually very similar. For example, the range ["Hello", "world"] is increasingly similar with the range ["Hello", "world", "world", "world",...] as more instances of "world" are appended. To prevent that, gapWeightedSimilarityNormalized computes a normalized version of the similarity that is computed as gapWeightedSimilarity(s, t, lambda) / sqrt(gapWeightedSimilarity(s, t, lambda) * gapWeightedSimilarity(s, t, lambda)). The function gapWeightedSimilarityNormalized (a so-called normalized kernel) is bounded in [0, 1], reaches 0 only for ranges that don't match in any position, and 1 only for identical ranges.
The optional parameters sSelfSim and tSelfSim are meant for avoiding duplicate computation. Many applications may have already computed gapWeightedSimilarity(s, s, lambda) and/or gapWeightedSimilarity(t, t, lambda). In that case, they can be passed as sSelfSim and tSelfSim, respectively.
Examples:
import std.math : approxEqual, sqrt;

string[] s = ["Hello", "brave", "new", "world"];
string[] t = ["Hello", "new", "world"];
writeln(gapWeightedSimilarity(s, s, 1)); // 15
writeln(gapWeightedSimilarity(t, t, 1)); // 7
writeln(gapWeightedSimilarity(s, t, 1)); // 7
assert(approxEqual(gapWeightedSimilarityNormalized(s, t, 1),
                7.0 / sqrt(15.0 * 7), 0.01));
struct GapWeightedSimilarityIncremental(Range, F = double) if (isRandomAccessRange!Range && hasLength!Range);

GapWeightedSimilarityIncremental!(R, F) gapWeightedSimilarityIncremental(R, F)(R r1, R r2, F penalty);
Similar to gapWeightedSimilarity, just works in an incremental manner by first revealing the matches of length 1, then gapped matches of length 2, and so on. The memory requirement is Ο(s.length * t.length). The time complexity is Ο(s.length * t.length) time for computing each step. Continuing on the previous example:
The implementation is based on the pseudocode in Fig. 4 of the paper "Efficient Computation of Gapped Substring Kernels on Large Alphabets" by Rousu et al., with additional algorithmic and systems-level optimizations.
Examples:
string[] s = ["Hello", "brave", "new", "world"];
string[] t = ["Hello", "new", "world"];
auto simIter = gapWeightedSimilarityIncremental(s, t, 1.0);
assert(simIter.front == 3); // three 1-length matches
simIter.popFront();
assert(simIter.front == 3); // three 2-length matches
simIter.popFront();
assert(simIter.front == 1); // one 3-length match
simIter.popFront();
assert(simIter.empty);     // no more match
this(Range s, Range t, F lambda);
Constructs an object given two ranges s and t and a penalty lambda. Constructor completes in Ο(s.length * t.length) time and computes all matches of length 1.
ref GapWeightedSimilarityIncremental opSlice();
Returns:
this.
void popFront();
Computes the match of the popFront length. Completes in Ο(s.length * t.length) time.
@property F front();
Returns:
The gapped similarity at the current match length (initially 1, grows with each call to popFront).
@property bool empty();
Returns:
Whether there are more matches.
T gcd(T)(T a, T b)
if (isIntegral!T);

T gcd(T)(T a, T b)
if (!isIntegral!T && is(typeof(T.init % T.init)) && is(typeof(T.init == 0 || T.init > 0)));
Computes the greatest common divisor of a and b by using an efficient algorithm such as Euclid's or Stein's algorithm.
Parameters:
T Any numerical type that supports the modulo operator %. If bit-shifting << and >> are also supported, Stein's algorithm will be used; otherwise, Euclid's algorithm is used as a fallback.
Returns:
The greatest common divisor of the given arguments.
Examples:
writeln(gcd(2 * 5 * 7 * 7, 5 * 7 * 11)); // 5 * 7
const int a = 5 * 13 * 23 * 23, b = 13 * 59;
writeln(gcd(a, b)); // 13
class Fft;
A class for performing fast Fourier transforms of power of two sizes. This class encapsulates a large amount of state that is reusable when performing multiple FFTs of sizes smaller than or equal to that specified in the constructor. This results in substantial speedups when performing multiple FFTs with a known maximum size. However, a free function API is provided for convenience if you need to perform a one-off FFT.
this(size_t size);
Create an Fft object for computing fast Fourier transforms of power of two sizes of size or smaller. size must be a power of two.
const Complex!F[] fft(F = double, R)(R range)
if (isFloatingPoint!F && isRandomAccessRange!R);
Compute the Fourier transform of range using the Ο(N log N) Cooley-Tukey Algorithm. range must be a random-access range with slicing and a length equal to size as provided at the construction of this object. The contents of range can be either numeric types, which will be interpreted as pure real values, or complex types with properties or members .re and .im that can be read.

Note Pure real FFTs are automatically detected and the relevant optimizations are performed.

Returns:
An array of complex numbers representing the transformed data in the frequency domain.

Conventions The exponent is negative and the factor is one, i.e., output[j] := sum[ exp(-2 PI i j k / N) input[k] ].

const void fft(Ret, R)(R range, Ret buf)
if (isRandomAccessRange!Ret && isComplexLike!(ElementType!Ret) && hasSlicing!Ret);
Same as the overload, but allows for the results to be stored in a user- provided buffer. The buffer must be of the same length as range, must be a random-access range, must have slicing, and must contain elements that are complex-like. This means that they must have a .re and a .im member or property that can be both read and written and are floating point numbers.
const Complex!F[] inverseFft(F = double, R)(R range)
if (isRandomAccessRange!R && isComplexLike!(ElementType!R) && isFloatingPoint!F);
Computes the inverse Fourier transform of a range. The range must be a random access range with slicing, have a length equal to the size provided at construction of this object, and contain elements that are either of type std.complex.Complex or have essentially the same compile-time interface.
Returns:
The time-domain signal.

Conventions The exponent is positive and the factor is 1/N, i.e., output[j] := (1 / N) sum[ exp(+2 PI i j k / N) input[k] ].

const void inverseFft(Ret, R)(R range, Ret buf)
if (isRandomAccessRange!Ret && isComplexLike!(ElementType!Ret) && hasSlicing!Ret);
Inverse FFT that allows a user-supplied buffer to be provided. The buffer must be a random access range with slicing, and its elements must be some complex-like type.
Complex!F[] fft(F = double, R)(R range);

void fft(Ret, R)(R range, Ret buf);

Complex!F[] inverseFft(F = double, R)(R range);

void inverseFft(Ret, R)(R range, Ret buf);
Convenience functions that create an Fft object, run the FFT or inverse FFT and return the result. Useful for one-off FFTs.

Note In addition to convenience, these functions are slightly more efficient than manually creating an Fft object for a single use, as the Fft object is deterministically destroyed before these functions return.