Skip to content

Polynomial

beignet.add_polynomial

add_polynomial(input, other)

Returns the sum of two polynomials.

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required
other Tensor

Polynomial coefficients.

required

Returns:

Name Type Description
output Tensor

Polynomial coefficients.

Source code in src/beignet/_add_polynomial.py
 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
def add_polynomial(input: Tensor, other: Tensor) -> Tensor:
    r"""
    Returns the sum of two polynomials.

    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    other : Tensor
        Polynomial coefficients.

    Returns
    -------
    output : Tensor
        Polynomial coefficients.
    """
    input = torch.atleast_1d(input)
    other = torch.atleast_1d(other)

    dtype = torch.promote_types(input.dtype, other.dtype)

    input = input.to(dtype)
    other = other.to(dtype)

    if input.shape[0] > other.shape[0]:
        output = torch.concatenate(
            [
                other,
                torch.zeros(
                    input.shape[0] - other.shape[0],
                    dtype=other.dtype,
                ),
            ],
        )

        output = input + output
    else:
        output = torch.concatenate(
            [
                input,
                torch.zeros(
                    other.shape[0] - input.shape[0],
                    dtype=input.dtype,
                ),
            ]
        )

        output = other + output

    return output

beignet.differentiate_polynomial

differentiate_polynomial(input, order=None, scale=None, dim=0)

Returns the derivative of a polynomial.

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required
order Tensor
None
scale Tensor
None
dim int
0

Returns:

Name Type Description
output Tensor

Polynomial coefficients of the derivative.

Source code in src/beignet/_differentiate_polynomial.py
 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
def differentiate_polynomial(
    input: Tensor,
    order: Tensor | None = None,
    scale: Tensor | None = None,
    dim: int = 0,
) -> Tensor:
    r"""
    Returns the derivative of a polynomial.

    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    order : Tensor, optional

    scale : Tensor, optional

    dim : int, default=0

    Returns
    -------
    output : Tensor
        Polynomial coefficients of the derivative.
    """
    input = torch.atleast_1d(input)

    if order == 0:
        return input

    input = torch.moveaxis(input, dim, 0)

    if order >= input.shape[0]:
        output = torch.zeros_like(input[:1])
    else:
        d = torch.arange(input.shape[0])

        output = input

        for _ in range(0, order):
            output = (d * output.T).T

            output = torch.roll(output, -1, dims=[0]) * scale

            output[-1] = 0.0

        output = output[:-order]

    output = torch.moveaxis(output, 0, dim)

    return output

beignet.divide_polynomial

divide_polynomial(input, other)

Returns the quotient and remainder of two polynomials.

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required
other Tensor

Polynomial coefficients.

required

Returns:

Name Type Description
output Tuple[Tensor, Tensor]

Polynomial coefficients of the quotient and remainder.

Source code in src/beignet/_divide_polynomial.py
  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
def divide_polynomial(
    input: Tensor,
    other: Tensor,
) -> Tuple[Tensor, Tensor]:
    r"""
    Returns the quotient and remainder of two polynomials.

    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    other : Tensor
        Polynomial coefficients.

    Returns
    -------
    output : Tuple[Tensor, Tensor]
        Polynomial coefficients of the quotient and remainder.
    """
    input = torch.atleast_1d(input)
    other = torch.atleast_1d(other)

    dtype = torch.promote_types(input.dtype, other.dtype)

    input = input.to(dtype)
    other = other.to(dtype)

    m = input.shape[0]
    n = other.shape[0]

    if m < n:
        return torch.zeros_like(input[:1]), input

    if n == 1:
        return input / other[-1], torch.zeros_like(input[:1])

    def f(x: Tensor) -> Tensor:
        indicies = torch.flip(x, [0])

        indicies = torch.nonzero(indicies, as_tuple=False)

        if indicies.shape[0] > 1:
            indicies = indicies[:1]

        if indicies.shape[0] < 1:
            indicies = torch.concatenate(
                [
                    indicies,
                    torch.full(
                        [
                            1 - indicies.shape[0],
                            indicies.shape[1],
                        ],
                        0,
                    ),
                ],
                0,
            )

        return x.shape[0] - 1 - indicies[0][0]

    quotient = torch.zeros(m - n + 1, dtype=input.dtype)

    ridx = input.shape[0] - 1

    size = m - f(other) - 1

    y = torch.zeros(m + n + 1, dtype=input.dtype)

    y[size] = 1.0

    x = quotient, input, y, ridx

    for index in range(0, size):
        quotient, remainder, y2, ridx1 = x

        j = size - index

        p = multiply_polynomial(y2, other)

        pidx = f(p)

        t = remainder[ridx1] / p[pidx]

        remainder_modified = remainder.clone()
        remainder_modified[ridx1] = 0.0

        a = remainder_modified

        p_modified = p.clone()
        p_modified[pidx] = 0.0

        b = t * p_modified

        a = torch.atleast_1d(a)
        b = torch.atleast_1d(b)

        dtype = torch.promote_types(a.dtype, b.dtype)

        a = a.to(dtype)
        b = b.to(dtype)

        if a.shape[0] > b.shape[0]:
            output = -b

            output = torch.concatenate(
                [
                    output,
                    torch.zeros(
                        a.shape[0] - b.shape[0],
                        dtype=b.dtype,
                    ),
                ],
            )
            output = a + output
        else:
            output = -b

            output = torch.concatenate(
                [
                    output[: a.shape[0]] + a,
                    output[a.shape[0] :],
                ],
            )

        remainder = output

        remainder = remainder[: remainder.shape[0]]

        quotient[j] = t

        ridx1 = ridx1 - 1

        y2 = torch.roll(y2, -1)

        x = quotient, remainder, y2, ridx1

    quotient, remainder, _, _ = x

    return quotient, remainder

beignet.evaluate_polynomial

evaluate_polynomial(input, coefficients, tensor=True)

Parameters:

Name Type Description Default
input Tensor
required
coefficients Tensor
required
tensor bool
True

Returns:

Name Type Description
output Tensor
Source code in src/beignet/_evaluate_polynomial.py
 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
def evaluate_polynomial(
    input: Tensor, coefficients: Tensor, tensor: bool = True
) -> Tensor:
    r"""
    Parameters
    ----------
    input : Tensor

    coefficients : Tensor

    tensor : bool

    Returns
    -------
    output : Tensor
    """
    coefficients = torch.atleast_1d(coefficients)

    if tensor:
        coefficients = torch.reshape(
            coefficients,
            coefficients.shape + (1,) * input.ndim,
        )

    output = coefficients[-1] + torch.zeros_like(input)

    for i in range(2, coefficients.shape[0] + 1):
        output = coefficients[-i] + output * input

    return output

beignet.evaluate_polynomial_2d

evaluate_polynomial_2d(x, y, coefficients)

Parameters:

Name Type Description Default
x Tensor
required
y Tensor
required
coefficients Tensor
required

Returns:

Name Type Description
output Tensor
Source code in src/beignet/_evaluate_polynomial_2d.py
 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
def evaluate_polynomial_2d(x: Tensor, y: Tensor, coefficients: Tensor) -> Tensor:
    r"""
    Parameters
    ----------
    x : Tensor

    y : Tensor

    coefficients : Tensor

    Returns
    -------
    output : Tensor
    """
    points = [x, y]

    if not all(a.shape == points[0].shape for a in points[1:]):
        match len(points):
            case 2:
                raise ValueError
            case 3:
                raise ValueError
            case _:
                raise ValueError

    points = iter(points)

    output = evaluate_polynomial(
        next(points),
        coefficients,
    )

    for x in points:
        output = evaluate_polynomial(
            x,
            output,
            tensor=False,
        )

    return output

beignet.evaluate_polynomial_3d

evaluate_polynomial_3d(x, y, z, coefficients)

Parameters:

Name Type Description Default
x Tensor
required
y Tensor
required
z Tensor
required
coefficients Tensor
required

Returns:

Name Type Description
output Tensor
Source code in src/beignet/_evaluate_polynomial_3d.py
 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
def evaluate_polynomial_3d(
    x: Tensor,
    y: Tensor,
    z: Tensor,
    coefficients: Tensor,
) -> Tensor:
    r"""
    Parameters
    ----------
    x : Tensor

    y : Tensor

    z : Tensor

    coefficients : Tensor

    Returns
    -------
    output : Tensor
    """
    points = [x, y, z]

    if not all(a.shape == points[0].shape for a in points[1:]):
        match len(points):
            case 2:
                raise ValueError
            case 3:
                raise ValueError
            case _:
                raise ValueError

    points = iter(points)

    output = evaluate_polynomial(
        next(points),
        coefficients,
    )

    for x in points:
        output = evaluate_polynomial(
            x,
            output,
            tensor=False,
        )

    return output

beignet.evaluate_polynomial_cartesian_2d

evaluate_polynomial_cartesian_2d(x, y, coefficients)

Parameters:

Name Type Description Default
x Tensor
required
y Tensor
required
coefficients Tensor
required

Returns:

Name Type Description
output Tensor
Source code in src/beignet/_evaluate_polynomial_cartesian_2d.py
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
def evaluate_polynomial_cartesian_2d(
    x: Tensor, y: Tensor, coefficients: Tensor
) -> Tensor:
    r"""
    Parameters
    ----------
    x : Tensor

    y : Tensor

    coefficients : Tensor

    Returns
    -------
    output : Tensor
    """
    for input in [x, y]:
        coefficients = evaluate_polynomial(input, coefficients)

    return coefficients

beignet.evaluate_polynomial_cartesian_3d

evaluate_polynomial_cartesian_3d(x, y, z, coefficients)

Parameters:

Name Type Description Default
x Tensor
required
y Tensor
required
z Tensor
required
coefficients Tensor
required

Returns:

Name Type Description
out Tensor
Source code in src/beignet/_evaluate_polynomial_cartesian_3d.py
 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
def evaluate_polynomial_cartesian_3d(
    x: Tensor, y: Tensor, z: Tensor, coefficients: Tensor
) -> Tensor:
    r"""
    Parameters
    ----------
    x : Tensor

    y : Tensor

    z : Tensor

    coefficients : Tensor

    Returns
    -------
    out : Tensor
    """
    for input in [x, y, z]:
        coefficients = evaluate_polynomial(
            input,
            coefficients,
        )

    return coefficients

beignet.evaluate_polynomial_from_roots

evaluate_polynomial_from_roots(input, other, tensor=True)
Source code in src/beignet/_evaluate_polynomial_from_roots.py
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
def evaluate_polynomial_from_roots(
    input: Tensor,
    other: Tensor,
    tensor: bool = True,
) -> Tensor:
    if other.ndim == 0:
        other = torch.ravel(other)

    if tensor:
        other = torch.reshape(other, other.shape + (1,) * input.ndim)

    if input.ndim >= other.ndim:
        raise ValueError

    output = torch.prod(input - other, dim=0)

    return output

beignet.fit_polynomial

fit_polynomial(input, other, degree, relative_condition=None, full=False, weight=None)

Parameters:

Name Type Description Default
input Tensor

Independent variable.

required
other Tensor

Dependent variable.

required
degree Tensor or int

Degree of the fitting polynomial.

required
relative_condition float

Relative condition number.

None
full bool

Return additional information.

False
weight Tensor

Weights.

None

Returns:

Name Type Description
output Tensor

Polynomial coefficients of the fit.

Source code in src/beignet/_fit_polynomial.py
  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
def fit_polynomial(
    input: Tensor,
    other: Tensor,
    degree: Tensor | int,
    relative_condition: float | None = None,
    full: bool = False,
    weight: Tensor | None = None,
):
    r"""
    Parameters
    ----------
    input : Tensor
        Independent variable.

    other : Tensor
        Dependent variable.

    degree : Tensor or int
        Degree of the fitting polynomial.

    relative_condition : float, optional
        Relative condition number.

    full : bool, default=False
        Return additional information.

    weight : Tensor, optional
        Weights.

    Returns
    -------
    output : Tensor
        Polynomial coefficients of the fit.
    """
    input = torch.tensor(input)
    other = torch.tensor(other)
    degree = torch.tensor(degree)
    if degree.ndim > 1:
        raise TypeError
    # if deg.dtype.kind not in "iu":
    #     raise TypeError
    if math.prod(degree.shape) == 0:
        raise TypeError
    if degree.min() < 0:
        raise ValueError
    if input.ndim != 1:
        raise TypeError
    if input.size == 0:
        raise TypeError
    if other.ndim < 1 or other.ndim > 2:
        raise TypeError
    if len(input) != len(other):
        raise TypeError
    if degree.ndim == 0:
        lmax = int(degree)
        van = polynomial_vandermonde(input, lmax)
    else:
        degree, _ = torch.sort(degree)
        lmax = int(degree[-1])
        van = polynomial_vandermonde(input, lmax)[:, degree]
    # set up the least squares matrices in transposed form
    lhs = van.T
    rhs = other.T
    if weight is not None:
        if weight.ndim != 1:
            raise TypeError("expected 1D vector for w")

        if len(input) != len(weight):
            raise TypeError("expected x and w to have same length")

        # apply weights. Don't use inplace operations as they
        # can cause problems with NA.
        lhs = lhs * weight
        rhs = rhs * weight
    # set rcond
    if relative_condition is None:
        relative_condition = len(input) * torch.finfo(input.dtype).eps
    # Determine the norms of the design matrix columns.
    if torch.is_complex(lhs):
        scl = torch.sqrt((torch.square(lhs.real) + torch.square(lhs.imag)).sum(1))
    else:
        scl = torch.sqrt(torch.square(lhs).sum(1))
    scl = torch.where(scl == 0, 1, scl)
    # Solve the least squares problem.
    c, resids, rank, s = torch.linalg.lstsq(lhs.T / scl, rhs.T, relative_condition)
    c = (c.T / scl).T
    # Expand c to include non-fitted coefficients which are set to zero
    if degree.ndim > 0:
        if c.ndim == 2:
            cc = torch.zeros((lmax + 1, c.shape[1]), dtype=c.dtype)
        else:
            cc = torch.zeros(lmax + 1, dtype=c.dtype)

        cc[degree] = c

        c = cc
    if full:
        result = c, [resids, rank, s, relative_condition]
    else:
        result = c
    return result

beignet.integrate_polynomial

integrate_polynomial(input, order=1, k=None, lower_bound=0, scale=1, dim=0)

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required
order int
1
k int
None
lower_bound float
0
scale float
1
dim int
0

Returns:

Name Type Description
output Tensor

Polynomial coefficients of the integral.

Source code in src/beignet/_integrate_polynomial.py
  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
def integrate_polynomial(
    input: Tensor,
    order: int = 1,
    k=None,
    lower_bound: float = 0,
    scale: float = 1,
    dim: int = 0,
) -> Tensor:
    r"""
    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    order : int, default=1

    k : int, optional

    lower_bound : float, default=0

    scale : float, default=1

    dim : int, default=0

    Returns
    -------
    output : Tensor
        Polynomial coefficients of the integral.
    """
    if k is None:
        k = []

    input = torch.atleast_1d(input)

    lower_bound, scale = map(torch.tensor, (lower_bound, scale))

    if not numpy.iterable(k):
        k = [k]

    if order < 0:
        raise ValueError

    if len(k) > order:
        raise ValueError

    if lower_bound.ndim != 0:
        raise ValueError

    if scale.ndim != 0:
        raise ValueError

    if order == 0:
        return input

    k = torch.tensor(list(k) + [0] * (order - len(k)))
    k = torch.atleast_1d(k)

    n = input.shape[dim]

    padding = torch.tensor([0, order])
    input = torch.moveaxis(input, dim, 0)
    if padding[0] < 0:
        input = input[torch.abs(padding[0]) :]

        padding = (0, padding[1])
    if padding[1] < 0:
        input = input[: -torch.abs(padding[1])]

        padding = (padding[0], 0)
    npad = torch.tensor([(0, 0)] * input.ndim)
    npad[0] = padding
    output = torch._numpy._funcs_impl.pad(
        input,
        pad_width=npad,
        mode="constant",
        constant_values=0,
    )
    input = torch.moveaxis(output, 0, dim)

    input = torch.moveaxis(input, dim, 0)

    d = torch.arange(n + order) + 1

    for i in range(0, order):
        input = input * scale

        input = (input.T / d).T

        input = torch.roll(input, 1, dims=[0])

        input[0] = 0.0

        input[0] += k[i] - evaluate_polynomial(lower_bound, input)

    return torch.moveaxis(input, 0, dim)

beignet.linear_polynomial

linear_polynomial(input, other)
Source code in src/beignet/_linear_polynomial.py
5
6
def linear_polynomial(input: float, other: float) -> Tensor:
    return torch.tensor([input, other])

beignet.multiply_polynomial

multiply_polynomial(input, other, mode='full')

Returns the product of two polynomials.

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required
other Tensor

Polynomial coefficients.

required

Returns:

Name Type Description
output Tensor

Polynomial coefficients of the product.

Source code in src/beignet/_multiply_polynomial.py
 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
def multiply_polynomial(
    input: Tensor,
    other: Tensor,
    mode: Literal["full", "same", "valid"] = "full",
) -> Tensor:
    r"""
    Returns the product of two polynomials.

    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    other : Tensor
        Polynomial coefficients.

    Returns
    -------
    output : Tensor
        Polynomial coefficients of the product.
    """
    input = torch.atleast_1d(input)
    other = torch.atleast_1d(other)

    dtype = torch.promote_types(input.dtype, other.dtype)

    input = input.to(dtype)
    other = other.to(dtype)

    output = convolve(input, other)

    if mode == "same":
        output = output[: max(input.shape[0], other.shape[0])]

    return output

beignet.multiply_polynomial_by_x

multiply_polynomial_by_x(input, mode='full')

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required
mode Literal['full', 'same', 'valid']
'full'

Returns:

Name Type Description
output Tensor

Polynomial coefficients of the product of the polynomial and the independent variable.

Source code in src/beignet/_multiply_polynomial_by_x.py
 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
def multiply_polynomial_by_x(
    input: Tensor,
    mode: Literal["full", "same", "valid"] = "full",
) -> Tensor:
    r"""
    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    mode : Literal["full", "same", "valid"]

    Returns
    -------
    output : Tensor
        Polynomial coefficients of the product of the polynomial and the
        independent variable.
    """
    input = torch.atleast_1d(input)

    output = torch.zeros(input.shape[0] + 1, dtype=input.dtype)

    output[1:] = input

    if mode == "same":
        output = output[: input.shape[0]]

    return output

beignet.polynomial_companion

polynomial_companion(input)

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required

Returns:

Name Type Description
output Tensor, shape=(degree, degree)

Companion matrix.

Source code in src/beignet/_polynomial_companion.py
 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
def polynomial_companion(input: Tensor) -> Tensor:
    r"""
    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    Returns
    -------
    output : Tensor, shape=(degree, degree)
        Companion matrix.
    """
    input = torch.atleast_1d(input)

    if input.shape[0] < 2:
        raise ValueError

    if input.shape[0] == 2:
        output = torch.tensor([[-input[0] / input[1]]])
    else:
        n = input.shape[0] - 1

        output = torch.reshape(torch.zeros([n, n], dtype=input.dtype), [-1])

        output[n :: n + 1] = 1.0

        output = torch.reshape(output, [n, n])

        output[:, -1] = output[:, -1] + (-input[:-1] / input[-1])

    return output

beignet.polynomial_domain module-attribute

polynomial_domain = tensor([-1.0, 1.0])

beignet.polynomial_from_roots

polynomial_from_roots(input)

Parameters:

Name Type Description Default
input Tensor

Roots.

required

Returns:

Name Type Description
output Tensor

Polynomial coefficients.

Source code in src/beignet/_polynomial_from_roots.py
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
def polynomial_from_roots(input: Tensor) -> Tensor:
    r"""
    Parameters
    ----------
    input : Tensor
        Roots.

    Returns
    -------
    output : Tensor
        Polynomial coefficients.
    """
    f = linear_polynomial
    g = multiply_polynomial
    if math.prod(input.shape) == 0:
        return torch.ones([1])

    input, _ = torch.sort(input)

    ys = []

    for x in input:
        a = torch.zeros(input.shape[0] + 1, dtype=x.dtype)
        b = f(-x, 1)

        a = torch.atleast_1d(a)
        b = torch.atleast_1d(b)

        dtype = torch.promote_types(a.dtype, b.dtype)

        a = a.to(dtype)
        b = b.to(dtype)

        if a.shape[0] > b.shape[0]:
            y = torch.concatenate(
                [
                    b,
                    torch.zeros(
                        a.shape[0] - b.shape[0],
                        dtype=b.dtype,
                    ),
                ],
            )

            y = a + y
        else:
            y = torch.concatenate(
                [
                    a,
                    torch.zeros(
                        b.shape[0] - a.shape[0],
                        dtype=a.dtype,
                    ),
                ]
            )

            y = b + y

        ys = [*ys, y]

    p = torch.stack(ys)

    m = p.shape[0]

    x = m, p

    while x[0] > 1:
        m, r = divmod(x[0], 2)

        z = x[1]

        previous = torch.zeros([len(p), input.shape[0] + 1])

        y = previous

        for i in range(0, m):
            y[i] = g(z[i], z[i + m])[: input.shape[0] + 1]

        previous = y

        if r:
            previous[0] = g(previous[0], z[2 * m])[: input.shape[0] + 1]

        x = m, previous

    _, output = x

    return output[0]

beignet.polynomial_one module-attribute

polynomial_one = tensor([1.0])

beignet.polynomial_power

polynomial_power(input, exponent, maximum_exponent=16.0)

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required
exponent float or Tensor
required
maximum_exponent float or Tensor
16.0

Returns:

Name Type Description
output Tensor

Polynomial coefficients of the power.

Source code in src/beignet/_polynomial_power.py
 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
def polynomial_power(
    input: Tensor,
    exponent: float | Tensor,
    maximum_exponent: float | Tensor = 16.0,
) -> Tensor:
    r"""
    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    exponent : float or Tensor

    maximum_exponent : float or Tensor, default=16.0

    Returns
    -------
    output : Tensor
        Polynomial coefficients of the power.
    """
    input = torch.atleast_1d(input)
    _exponent = int(exponent)
    if _exponent != exponent or _exponent < 0:
        raise ValueError
    if maximum_exponent is not None and _exponent > maximum_exponent:
        raise ValueError
    match _exponent:
        case 0:
            output = torch.tensor([1], dtype=input.dtype)
        case 1:
            output = input
        case _:
            output = torch.zeros(input.shape[0] * exponent, dtype=input.dtype)

            input = torch.atleast_1d(input)
            output = torch.atleast_1d(output)

            dtype = torch.promote_types(input.dtype, output.dtype)

            input = input.to(dtype)
            output = output.to(dtype)

            if output.shape[0] > input.shape[0]:
                input = torch.concatenate(
                    [
                        input,
                        torch.zeros(
                            output.shape[0] - input.shape[0],
                            dtype=input.dtype,
                        ),
                    ],
                )

                output = output + input
            else:
                output = torch.concatenate(
                    [
                        output,
                        torch.zeros(
                            input.shape[0] - output.shape[0],
                            dtype=output.dtype,
                        ),
                    ]
                )

                output = input + output

            for _ in range(2, _exponent + 1):
                output = multiply_polynomial(output, input, mode="same")
    return output

beignet.polynomial_roots

polynomial_roots(input)

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required

Returns:

Name Type Description
output Tensor

Roots.

Source code in src/beignet/_polynomial_roots.py
 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
def polynomial_roots(input: Tensor) -> Tensor:
    r"""
    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    Returns
    -------
    output : Tensor
        Roots.
    """
    input = torch.atleast_1d(input)

    if input.shape[0] < 2:
        return torch.tensor([], dtype=input.dtype)

    if input.shape[0] == 2:
        return torch.tensor([-input[0] / input[1]])

    if input.shape[0] < 2:
        raise ValueError

    if input.shape[0] == 2:
        output = torch.tensor([[-input[0] / input[1]]])
    else:
        n = input.shape[0] - 1

        output = torch.reshape(torch.zeros([n, n], dtype=input.dtype), [-1])

        output[n :: n + 1] = 1.0

        output = torch.reshape(output, [n, n])

        output[:, -1] = output[:, -1] + (-input[:-1] / input[-1])

    output = torch.flip(output, dims=[0])
    output = torch.flip(output, dims=[1])

    output = torch.linalg.eigvals(output)

    output, _ = torch.sort(output.real)

    return output

beignet.polynomial_to_chebyshev_polynomial

polynomial_to_chebyshev_polynomial(input)
Source code in src/beignet/_polynomial_to_chebyshev_polynomial.py
 8
 9
10
11
12
13
14
15
16
17
18
def polynomial_to_chebyshev_polynomial(input: Tensor) -> Tensor:
    input = torch.atleast_1d(input)

    output = torch.zeros_like(input)

    for i in range(0, input.shape[0] - 1 + 1):
        output = multiply_chebyshev_polynomial_by_x(output, mode="same")

        output = add_chebyshev_polynomial(output, input[input.shape[0] - 1 - i])

    return output

beignet.polynomial_to_laguerre_polynomial

polynomial_to_laguerre_polynomial(input)
Source code in src/beignet/_polynomial_to_laguerre_polynomial.py
 8
 9
10
11
12
13
14
15
16
17
18
def polynomial_to_laguerre_polynomial(input: Tensor) -> Tensor:
    input = torch.atleast_1d(input)

    output = torch.zeros_like(input)

    for i in range(0, input.shape[0]):
        output = multiply_laguerre_polynomial_by_x(output, mode="same")

        output = add_laguerre_polynomial(output, torch.flip(input, dims=[0])[i])

    return output

beignet.polynomial_to_legendre_polynomial

polynomial_to_legendre_polynomial(input)
Source code in src/beignet/_polynomial_to_legendre_polynomial.py
 8
 9
10
11
12
13
14
15
16
17
18
def polynomial_to_legendre_polynomial(input: Tensor) -> Tensor:
    input = torch.atleast_1d(input)

    output = torch.zeros_like(input)

    for i in range(0, input.shape[0] - 1 + 1):
        output = multiply_legendre_polynomial_by_x(output, mode="same")

        output = add_legendre_polynomial(output, input[input.shape[0] - 1 - i])

    return output

beignet.polynomial_to_physicists_hermite_polynomial

polynomial_to_physicists_hermite_polynomial(input)
Source code in src/beignet/_polynomial_to_physicists_hermite_polynomial.py
10
11
12
13
14
15
16
17
18
19
20
21
22
def polynomial_to_physicists_hermite_polynomial(input: Tensor) -> Tensor:
    input = torch.atleast_1d(input)

    output = torch.zeros_like(input)

    for index in range(0, input.shape[0] - 1 + 1):
        output = multiply_physicists_hermite_polynomial_by_x(output, mode="same")

        output = add_physicists_hermite_polynomial(
            output, input[input.shape[0] - 1 - index]
        )

    return output

beignet.polynomial_to_probabilists_hermite_polynomial

polynomial_to_probabilists_hermite_polynomial(input)
Source code in src/beignet/_polynomial_to_probabilists_hermite_polynomial.py
10
11
12
13
14
15
16
17
18
19
20
21
22
def polynomial_to_probabilists_hermite_polynomial(input: Tensor) -> Tensor:
    input = torch.atleast_1d(input)

    output = torch.zeros_like(input)

    for i in range(0, input.shape[0] - 1 + 1):
        output = multiply_probabilists_hermite_polynomial_by_x(output, mode="same")

        output = add_probabilists_hermite_polynomial(
            output, input[input.shape[0] - 1 - i]
        )

    return output

beignet.polynomial_vandermonde

polynomial_vandermonde(input, degree)

Parameters:

Name Type Description Default
input Tensor
required
degree Tensor
required

Returns:

Name Type Description
output Tensor
Source code in src/beignet/_polynomial_vandermonde.py
 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
def polynomial_vandermonde(input: Tensor, degree: Tensor) -> Tensor:
    r"""
    Parameters
    ----------
    input : Tensor

    degree : Tensor

    Returns
    -------
    output : Tensor
    """
    if degree < 0:
        raise ValueError

    degree = int(degree)

    input = torch.atleast_1d(input)

    output = torch.empty([degree + 1, *input.shape], dtype=input.dtype)

    output[0] = torch.ones_like(input)

    for i in range(1, degree + 1):
        output[i] = output[i - 1] * input

    output = torch.moveaxis(output, 0, -1)

    return output

beignet.polynomial_vandermonde_2d

polynomial_vandermonde_2d(x, y, degree)

Parameters:

Name Type Description Default
x Tensor
required
y Tensor
required
degree Tensor
required

Returns:

Name Type Description
output Tensor
Source code in src/beignet/_polynomial_vandermonde_2d.py
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
def polynomial_vandermonde_2d(x: Tensor, y: Tensor, degree: Tensor) -> Tensor:
    r"""
    Parameters
    ----------
    x : Tensor

    y : Tensor

    degree : Tensor

    Returns
    -------
    output : Tensor
    """
    functions = (
        polynomial_vandermonde,
        polynomial_vandermonde,
    )

    n = len(functions)

    if n != len([x, y]):
        raise ValueError

    if n != len(degree):
        raise ValueError

    if n == 0:
        raise ValueError

    matrices = []

    for i in range(n):
        matrix = functions[i]((x, y)[i], degree[i])

        matrices = [
            *matrices,
            matrix[(..., *tuple(slice(None) if j == i else None for j in range(n)))],
        ]

    vandermonde = functools.reduce(
        operator.mul,
        matrices,
    )

    return torch.reshape(
        vandermonde,
        [*vandermonde.shape[: -len(degree)], -1],
    )

beignet.polynomial_vandermonde_3d

polynomial_vandermonde_3d(x, y, z, degree)

Parameters:

Name Type Description Default
x Tensor
required
y Tensor
required
z Tensor
required
degree Tensor
required

Returns:

Name Type Description
output Tensor
Source code in src/beignet/_polynomial_vandermonde_3d.py
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
def polynomial_vandermonde_3d(
    x: Tensor, y: Tensor, z: Tensor, degree: Tensor
) -> Tensor:
    r"""
    Parameters
    ----------
    x : Tensor

    y : Tensor

    z : Tensor

    degree : Tensor

    Returns
    -------
    output : Tensor
    """
    vandermonde_functions = (
        polynomial_vandermonde,
        polynomial_vandermonde,
        polynomial_vandermonde,
    )

    n = len(vandermonde_functions)

    if n != len([x, y, z]):
        raise ValueError

    if n != len(degree):
        raise ValueError

    if n == 0:
        raise ValueError

    matrices = []

    for i in range(n):
        matrix = vandermonde_functions[i]((x, y, z)[i], degree[i])

        matrices = [
            *matrices,
            matrix[(..., *tuple(slice(None) if j == i else None for j in range(n)))],
        ]

    vandermonde = functools.reduce(
        operator.mul,
        matrices,
    )

    return torch.reshape(
        vandermonde,
        [*vandermonde.shape[: -len(degree)], -1],
    )

beignet.polynomial_x module-attribute

polynomial_x = tensor([0.0, 1.0])

beignet.polynomial_zero module-attribute

polynomial_zero = tensor([0.0])

beignet.subtract_polynomial

subtract_polynomial(input, other)

Returns the difference of two polynomials.

Parameters:

Name Type Description Default
input Tensor

Polynomial coefficients.

required
other Tensor

Polynomial coefficients.

required

Returns:

Name Type Description
output Tensor

Polynomial coefficients of the difference.

Source code in src/beignet/_subtract_polynomial.py
 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
def subtract_polynomial(input: Tensor, other: Tensor) -> Tensor:
    r"""
    Returns the difference of two polynomials.

    Parameters
    ----------
    input : Tensor
        Polynomial coefficients.

    other : Tensor
        Polynomial coefficients.

    Returns
    -------
    output : Tensor
        Polynomial coefficients of the difference.
    """
    input = torch.atleast_1d(input)
    other = torch.atleast_1d(other)

    dtype = torch.promote_types(input.dtype, other.dtype)

    input = input.to(dtype)
    other = other.to(dtype)

    if input.shape[0] > other.shape[0]:
        output = -other

        output = torch.concatenate(
            [
                output,
                torch.zeros(
                    input.shape[0] - other.shape[0],
                    dtype=other.dtype,
                ),
            ],
        )
        output = input + output
    else:
        output = -other

        output = torch.concatenate(
            [
                output[: input.shape[0]] + input,
                output[input.shape[0] :],
            ],
        )

    return output

beignet.trim_polynomial_coefficients

trim_polynomial_coefficients(input, tol=0.0)
Source code in src/beignet/_trim_polynomial_coefficients.py
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
def trim_polynomial_coefficients(
    input: Tensor,
    tol: float = 0.0,
) -> Tensor:
    if tol < 0:
        raise ValueError

    input = torch.atleast_1d(input)

    indices = torch.nonzero(torch.abs(input) > tol)

    if indices.shape[0] == 0:
        output = input[:1] * 0
    else:
        output = input[: indices[-1] + 1]

    return output