quadrature.py 10.3 KB
Newer Older
Boris Bonev's avatar
Boris Bonev committed
1
2
3
4
# coding=utf-8

# SPDX-FileCopyrightText: Copyright (c) 2022 The torch-harmonics Authors. All rights reserved.
# SPDX-License-Identifier: BSD-3-Clause
5
#
Boris Bonev's avatar
Boris Bonev committed
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
# Redistribution and use in source and binary forms, with or without
# modification, are permitted provided that the following conditions are met:
#
# 1. Redistributions of source code must retain the above copyright notice, this
# list of conditions and the following disclaimer.
#
# 2. Redistributions in binary form must reproduce the above copyright notice,
# this list of conditions and the following disclaimer in the documentation
# and/or other materials provided with the distribution.
#
# 3. Neither the name of the copyright holder nor the names of its
# contributors may be used to endorse or promote products derived from
# this software without specific prior written permission.
#
# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
# AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
# IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
# DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE
# FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
# DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR
# SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
# CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
# OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
# OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#

Thorsten Kurth's avatar
Thorsten Kurth committed
32
33
34
from typing import Tuple, Optional
from torch_harmonics.cache import lru_cache
import math
Boris Bonev's avatar
Boris Bonev committed
35
import numpy as np
Thorsten Kurth's avatar
Thorsten Kurth committed
36
import torch
Boris Bonev's avatar
Boris Bonev committed
37

Thorsten Kurth's avatar
Thorsten Kurth committed
38
39
def _precompute_grid(n: int, grid: Optional[str]="equidistant", a: Optional[float]=0.0, b: Optional[float]=1.0,
                     periodic: Optional[bool]=False) -> Tuple[torch.Tensor, torch.Tensor]:
apaaris's avatar
apaaris committed
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
    r"""
    Precompute grid points and weights for various quadrature rules.
    
    Parameters
    -----------
    n : int
        Number of grid points
    grid : str, optional
        Grid type ("equidistant", "legendre-gauss", "lobatto", "equiangular"), by default "equidistant"
    a : float, optional
        Lower bound of interval, by default 0.0
    b : float, optional
        Upper bound of interval, by default 1.0
    periodic : bool, optional
        Whether the grid is periodic (only for equidistant), by default False
        
    Returns
    -------
    Tuple[torch.Tensor, torch.Tensor]
        Grid points and weights
        
    Raises
    ------
    ValueError
        If periodic is True for non-equidistant grids or unknown grid type
    """
Boris Bonev's avatar
Boris Bonev committed
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83

    if (grid != "equidistant") and periodic:
        raise ValueError(f"Periodic grid is only supported on equidistant grids.")

    # compute coordinates
    if grid == "equidistant":
        xlg, wlg = trapezoidal_weights(n, a=a, b=b, periodic=periodic)
    elif grid == "legendre-gauss":
        xlg, wlg = legendre_gauss_weights(n, a=a, b=b)
    elif grid == "lobatto":
        xlg, wlg = lobatto_weights(n, a=a, b=b)
    elif grid == "equiangular":
        xlg, wlg = clenshaw_curtiss_weights(n, a=a, b=b)
    else:
        raise ValueError(f"Unknown grid type {grid}")

    return xlg, wlg

Thorsten Kurth's avatar
Thorsten Kurth committed
84
85
@lru_cache(typed=True, copy=True)
def _precompute_longitudes(nlon: int):
86
   
Thorsten Kurth's avatar
Thorsten Kurth committed
87
88
89
    lons = torch.linspace(0, 2 * math.pi, nlon+1, dtype=torch.float64, requires_grad=False)[:-1]
    return lons

90

Thorsten Kurth's avatar
Thorsten Kurth committed
91
92
93
@lru_cache(typed=True, copy=True)
def _precompute_latitudes(nlat: int, grid: Optional[str]="equiangular") -> Tuple[torch.Tensor, torch.Tensor]:
        
94
    # compute coordinates in the cosine theta domain
Boris Bonev's avatar
Boris Bonev committed
95
    xlg, wlg = _precompute_grid(nlat, grid=grid, a=-1.0, b=1.0, periodic=False)
Thorsten Kurth's avatar
Thorsten Kurth committed
96
    
97
98
    # to perform the quadrature and account for the jacobian of the sphere, the quadrature rule
    # is formulated in the cosine theta domain, which is designed to integrate functions of cos theta
Thorsten Kurth's avatar
Thorsten Kurth committed
99
100
101
    lats = torch.flip(torch.arccos(xlg), dims=(0,)).clone()
    wlg = torch.flip(wlg, dims=(0,)).clone()
    
102
103
    return lats, wlg

104

Thorsten Kurth's avatar
Thorsten Kurth committed
105
def trapezoidal_weights(n: int, a: Optional[float]=-1.0, b: Optional[float]=1.0, periodic: Optional[bool]=False) -> Tuple[torch.Tensor, torch.Tensor]:
Boris Bonev's avatar
Boris Bonev committed
106
107
108
    r"""
    Helper routine which returns equidistant nodes with trapezoidal weights
    on the interval [a, b]
apaaris's avatar
apaaris committed
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126

    Parameters
    -----------
    n: int
        Number of quadrature nodes
    a: Optional[float]
        Lower bound of the interval
    b: Optional[float]
        Upper bound of the interval
    periodic: Optional[bool]
        Whether the grid is periodic

    Returns
    ------- 
    xlg: torch.Tensor
        Tensor of quadrature nodes  
    wlg: torch.Tensor
        Tensor of quadrature weights
Boris Bonev's avatar
Boris Bonev committed
127
128
    """

Thorsten Kurth's avatar
Thorsten Kurth committed
129
    xlg = torch.as_tensor(np.linspace(a, b, n, endpoint=periodic))
Thorsten Kurth's avatar
Thorsten Kurth committed
130
    wlg = (b - a) / (n - periodic * 1) * torch.ones(n, requires_grad=False)
Boris Bonev's avatar
Boris Bonev committed
131
132
133
134
135
136
137

    if not periodic:
        wlg[0] *= 0.5
        wlg[-1] *= 0.5

    return xlg, wlg

138

Thorsten Kurth's avatar
Thorsten Kurth committed
139
def legendre_gauss_weights(n: int, a: Optional[float]=-1.0, b: Optional[float]=1.0) -> Tuple[torch.Tensor, torch.Tensor]:
140
    r"""
Boris Bonev's avatar
Boris Bonev committed
141
142
    Helper routine which returns the Legendre-Gauss nodes and weights
    on the interval [a, b]
apaaris's avatar
apaaris committed
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158

    Parameters
    -----------
    n: int
        Number of quadrature nodes
    a: Optional[float]
        Lower bound of the interval
    b: Optional[float]
        Upper bound of the interval

    Returns
    -------
    xlg: torch.Tensor
        Tensor of quadrature nodes  
    wlg: torch.Tensor
        Tensor of quadrature weights
Boris Bonev's avatar
Boris Bonev committed
159
160
161
    """

    xlg, wlg = np.polynomial.legendre.leggauss(n)
Thorsten Kurth's avatar
Thorsten Kurth committed
162
163
    xlg = torch.as_tensor(xlg).clone()
    wlg = torch.as_tensor(wlg).clone()
Boris Bonev's avatar
Boris Bonev committed
164
165
166
167
168
    xlg = (b - a) * 0.5 * xlg + (b + a) * 0.5
    wlg = wlg * (b - a) * 0.5

    return xlg, wlg

169

Thorsten Kurth's avatar
Thorsten Kurth committed
170
171
def lobatto_weights(n: int, a: Optional[float]=-1.0, b: Optional[float]=1.0,
                    tol: Optional[float]=1e-16, maxiter: Optional[int]=100) -> Tuple[torch.Tensor, torch.Tensor]:
172
    r"""
Boris Bonev's avatar
Boris Bonev committed
173
174
    Helper routine which returns the Legendre-Gauss-Lobatto nodes and weights
    on the interval [a, b]
apaaris's avatar
apaaris committed
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195

    Parameters
    -----------
    n: int
        Number of quadrature nodes
    a: Optional[float]
        Lower bound of the interval
    b: Optional[float]
        Upper bound of the interval
    tol: Optional[float]
        Tolerance for the iteration
    maxiter: Optional[int]
        Maximum number of iterations

    Returns
    -------
    tlg: torch.Tensor
        Tensor of quadrature nodes
    wlg: torch.Tensor
        Tensor of quadrature weights

Boris Bonev's avatar
Boris Bonev committed
196
197
    """

Thorsten Kurth's avatar
Thorsten Kurth committed
198
199
200
    wlg = torch.zeros((n,), dtype=torch.float64, requires_grad=False)
    tlg = torch.zeros((n,), dtype=torch.float64, requires_grad=False)
    tmp = torch.zeros((n,), dtype=torch.float64, requires_grad=False)
Boris Bonev's avatar
Boris Bonev committed
201
202

    # Vandermonde Matrix
Thorsten Kurth's avatar
Thorsten Kurth committed
203
    vdm = torch.zeros((n, n), dtype=torch.float64, requires_grad=False)
204

Boris Bonev's avatar
Boris Bonev committed
205
    # initialize Chebyshev nodes as first guess
206
    for i in range(n):
207
        tlg[i] = -math.cos(math.pi * i / (n - 1))
208

Boris Bonev's avatar
Boris Bonev committed
209
    tmp = 2.0
210

Boris Bonev's avatar
Boris Bonev committed
211
212
    for i in range(maxiter):
        tmp = tlg
213
214
215
216

        vdm[:, 0] = 1.0
        vdm[:, 1] = tlg

Boris Bonev's avatar
Boris Bonev committed
217
        for k in range(2, n):
218
219
220
221
222
223
224
225
            vdm[:, k] = ((2 * k - 1) * tlg * vdm[:, k - 1] - (k - 1) * vdm[:, k - 2]) / k

        tlg = tmp - (tlg * vdm[:, n - 1] - vdm[:, n - 2]) / (n * vdm[:, n - 1])

        if max(abs(tlg - tmp).flatten()) < tol:
            break

    wlg = 2.0 / ((n * (n - 1)) * (vdm[:, n - 1] ** 2))
Boris Bonev's avatar
Boris Bonev committed
226
227
228
229

    # rescale
    tlg = (b - a) * 0.5 * tlg + (b + a) * 0.5
    wlg = wlg * (b - a) * 0.5
230

Boris Bonev's avatar
Boris Bonev committed
231
232
233
    return tlg, wlg


Thorsten Kurth's avatar
Thorsten Kurth committed
234
def clenshaw_curtiss_weights(n: int, a: Optional[float]=-1.0, b: Optional[float]=1.0) -> Tuple[torch.Tensor, torch.Tensor]:
235
    r"""
Boris Bonev's avatar
Boris Bonev committed
236
237
238
    Computation of the Clenshaw-Curtis quadrature nodes and weights.
    This implementation follows

apaaris's avatar
apaaris committed
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
    Parameters
    -----------
    n: int
        Number of quadrature nodes
    a: Optional[float]
        Lower bound of the interval
    b: Optional[float]
        Upper bound of the interval

    Returns
    -------
    tcc: torch.Tensor
        Tensor of quadrature nodes
    wcc: torch.Tensor
        Tensor of quadrature weights

    References
    ----------
Boris Bonev's avatar
Boris Bonev committed
257
258
259
    [1] Joerg Waldvogel, Fast Construction of the Fejer and Clenshaw-Curtis Quadrature Rules; BIT Numerical Mathematics, Vol. 43, No. 1, pp. 001–018.
    """

260
    assert n > 1
Boris Bonev's avatar
Boris Bonev committed
261

Thorsten Kurth's avatar
Thorsten Kurth committed
262
    tcc = torch.cos(torch.linspace(math.pi, 0, n, dtype=torch.float64, requires_grad=False))
Boris Bonev's avatar
Boris Bonev committed
263
264

    if n == 2:
265
        wcc = torch.as_tensor([1.0, 1.0], dtype=torch.float64)
Boris Bonev's avatar
Boris Bonev committed
266
267
268
    else:

        n1 = n - 1
Thorsten Kurth's avatar
Thorsten Kurth committed
269
        N = torch.arange(1, n1, 2, dtype=torch.float64)
Boris Bonev's avatar
Boris Bonev committed
270
271
272
        l = len(N)
        m = n1 - l

Thorsten Kurth's avatar
Thorsten Kurth committed
273
274
275
        v = torch.cat([2 / N / (N - 2), 1 / N[-1:], torch.zeros(m, dtype=torch.float64, requires_grad=False)])
        #v = 0 - v[:-1] - v[-1:0:-1]
        v = 0 - v[:-1] - torch.flip(v[1:], dims=(0,))
Boris Bonev's avatar
Boris Bonev committed
276

Thorsten Kurth's avatar
Thorsten Kurth committed
277
        g0 = -torch.ones(n1, dtype=torch.float64, requires_grad=False)
Boris Bonev's avatar
Boris Bonev committed
278
279
        g0[l] = g0[l] + n1
        g0[m] = g0[m] + n1
280
        g = g0 / (n1**2 - 1 + (n1 % 2))
Thorsten Kurth's avatar
Thorsten Kurth committed
281
282
        wcc = torch.fft.ifft(v + g).real
        wcc = torch.cat((wcc, wcc[:1]))
Boris Bonev's avatar
Boris Bonev committed
283
284
285
286
287
288
289

    # rescale
    tcc = (b - a) * 0.5 * tcc + (b + a) * 0.5
    wcc = wcc * (b - a) * 0.5

    return tcc, wcc

290

Thorsten Kurth's avatar
Thorsten Kurth committed
291
def fejer2_weights(n: int, a: Optional[float]=-1.0, b: Optional[float]=1.0) -> Tuple[torch.Tensor, torch.Tensor]:
292
    r"""
Boris Bonev's avatar
Boris Bonev committed
293
294
    Computation of the Fejer quadrature nodes and weights.

apaaris's avatar
apaaris committed
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
    Parameters
    -----------
    n: int
        Number of quadrature nodes
    a: Optional[float]
        Lower bound of the interval
    b: Optional[float]
        Upper bound of the interval

    Returns
    -------
    tcc: torch.Tensor
        Tensor of quadrature nodes
    wcc: torch.Tensor
        Tensor of quadrature weights

    References
    ----------
Boris Bonev's avatar
Boris Bonev committed
313
314
315
    [1] Joerg Waldvogel, Fast Construction of the Fejer and Clenshaw-Curtis Quadrature Rules; BIT Numerical Mathematics, Vol. 43, No. 1, pp. 001–018.
    """

316
    assert n > 2
Boris Bonev's avatar
Boris Bonev committed
317

Thorsten Kurth's avatar
Thorsten Kurth committed
318
    tcc = torch.cos(torch.linspace(math.pi, 0, n, dtype=torch.float64, requires_grad=False))
Boris Bonev's avatar
Boris Bonev committed
319
320

    n1 = n - 1
Thorsten Kurth's avatar
Thorsten Kurth committed
321
    N = torch.arange(1, n1, 2, dtype=torch.float64)
Boris Bonev's avatar
Boris Bonev committed
322
323
324
    l = len(N)
    m = n1 - l

Thorsten Kurth's avatar
Thorsten Kurth committed
325
326
327
    v = torch.cat([2 / N / (N - 2), 1 / N[-1:], torch.zeros(m, dtype=torch.float64, requires_grad=False)])
    #v = 0 - v[:-1] - v[-1:0:-1]
    v = 0 - v[:-1] - torch.flip(v[1:], dims=(0,))
Boris Bonev's avatar
Boris Bonev committed
328

Thorsten Kurth's avatar
Thorsten Kurth committed
329
330
    wcc = torch.fft.ifft(v).real
    wcc = torch.cat((wcc, wcc[:1]))
Boris Bonev's avatar
Boris Bonev committed
331
332
333
334
335

    # rescale
    tcc = (b - a) * 0.5 * tcc + (b + a) * 0.5
    wcc = wcc * (b - a) * 0.5

336
    return tcc, wcc