README.md 6.57 KB
Newer Older
Boris Bonev's avatar
Boris Bonev committed
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
<!-- 
SPDX-FileCopyrightText: Copyright (c) 2022 The torch-harmonics Authors. All rights reserved.

SPDX-License-Identifier: BSD-3-Clause
 
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.
-->

Boris Bonev's avatar
Boris Bonev committed
31
32
33
34
35
36
37
38
<p align="center">
    <img src="./images/logo/logo.png"  width="568">
</p>

<!-- # torch-harmonics: differentiable harmonic transforms -->

<!-- ## What is torch-harmonics? -->

Boris Bonev's avatar
Boris Bonev committed
39
`torch_harmonics` is a differentiable implementation of the Spherical Harmonic transform in PyTorch. It uses quadrature rules to compute the projection onto the associated Legendre polynomials and FFTs for the projection onto the harmonic basis. This algorithm tends to outperform others with better asymptotic scaling for most practical purposes.
Boris Bonev's avatar
Boris Bonev committed
40
41
42

`torch_harmonics` uses PyTorch primitives to implement these operations, making it fully differentiable. Moreover, the quadrature can be distributed onto multiple ranks making it spatially distributed.

Boris Bonev's avatar
Boris Bonev committed
43
`torch_harmonics` has been used to implement a variety of differentiable PDE solvers which generated the animations below. Moreover, it has enabled the development of spherical Fourier Neural Operators (SFNOs).
Boris Bonev's avatar
Boris Bonev committed
44
45
46
47


<table border="0" cellspacing="0" cellpadding="0">
    <tr>
Boris Bonev's avatar
Boris Bonev committed
48
        <td><img src="./images/sfno.gif"  width="240"></td>
Boris Bonev's avatar
Boris Bonev committed
49
50
51
        <td><img src="./images/zonal_jet.gif"  width="240"></td>
        <td><img src="./images/ginzburg-landau.gif"  width="240"></td>
        <td><img src="./images/allen-cahn.gif"  width="240"></td>
Boris Bonev's avatar
Boris Bonev committed
52
    </tr> 
Boris Bonev's avatar
Boris Bonev committed
53
<!--     <tr>
Boris Bonev's avatar
Boris Bonev committed
54
55
56
        <td style="text-align:center; border-style : hidden!important;">Shallow Water Eqns.</td>
        <td style="text-align:center; border-style : hidden!important;">Ginzburg-Landau Eqn.</td>
        <td style="text-align:center; border-style : hidden!important;">Allen-Cahn Eqn.</td>
Boris Bonev's avatar
Boris Bonev committed
57
    </tr>  -->
Boris Bonev's avatar
Boris Bonev committed
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
</table>


## Installation
Build in your environment using the Python package:

```
git clone git@github.com:NVIDIA/torch-harmonics.git
pip install ./torch_harmonics
```

Alternatively, use the Dockerfile to build your custom container after cloning:


```
git clone git@github.com:NVIDIA/torch-harmonics.git
cd torch_harmonics
docker build . -t torch_harmonics
docker run --gpus all -it --rm --ipc=host --ulimit memlock=-1 --ulimit stack=67108864 torch_harmonics
```

## Contributors

 - Boris Bonev (bbonev@nvidia.com)
 - Christian Hundt (chundt@nvidia.com)
 - Thorsten Kurth (tkurth@nvidia.com)
Boris Bonev's avatar
Boris Bonev committed
84
 - Nikola Kovachki (nkovachki@nvidia.com)
Boris Bonev's avatar
Boris Bonev committed
85
 - Jean Kossaifi (jkossaifi@nvidia.com)
Boris Bonev's avatar
Boris Bonev committed
86
87

## Implementation
Boris Bonev's avatar
Boris Bonev committed
88
The implementation follows the algorithm as presented in [2].
Boris Bonev's avatar
Boris Bonev committed
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

### Spherical harmonic transform

The truncated series expansion of a function $f$ defined on the surface of a sphere can be written as

$$
f(\theta, \lambda) = \sum_{m=-M}^{M} \exp(im\lambda) \sum_{n=|m|}^{M} F_n^m \bar{P}_n^m (\cos \theta),
$$

where $\theta$ is the colatitude, $\lambda$ the longitude, $\bar{P}_n^m$ the normalized, associated Legendre polynomials and $F_n^m$, the expansion coefficient associated to the mode $(m,n)$.

A direct spherical harmonic transform can be accomplished by a Fourier transform

$$
F^m(\theta) = \frac{1}{2 \pi} \int_{0}^{2\pi} f(\theta, \lambda) \exp(-im\lambda)  \mathrm{d}\lambda
$$

in longitude and a Legendre transform

$$
F_n^m = \frac{1}{2} \int_{-1}^1 F^m(\theta) \bar{P}_n^m(\cos \theta)  \mathrm{d} \cos \theta
$$

in latitude.

### Discrete Legendre transform

in order to apply the Legendre transfor, we shall use Gauss-Legendre points in the latitudinal direction. The integral

$$
F_n^m = \int_{0}^\pi F^m(\theta) \bar{P}_n^m(\cos \theta) \sin \theta \mathrm{d} \theta
$$

is approximated by the sum

$$
F_n^m = \sum_{j=1}^{N_\theta} F^m(\theta_j) \bar{P}_n^m(\cos \theta_j) w_j
$$

## Usage

### Getting started

The main functionality of `torch_harmonics` is provided in the form of `torch.nn.Modules` for composability. A minimum example is given by:

```python
import torch
Boris Bonev's avatar
Boris Bonev committed
136
import torch_harmonics as th
Boris Bonev's avatar
Boris Bonev committed
137
138
139
140
141
142
143
144
145

device = torch.device('cuda' if torch.cuda.is_available() else 'cpu')

nlat = 512
nlon = 2*nlat
batch_size = 32
signal = torch.randn(batch_size, nlat, nlon)

# transform data on an equiangular grid
Boris Bonev's avatar
Boris Bonev committed
146
sht = th.RealSHT(nlat, nlon, grid="equiangular").to(device).float()
Boris Bonev's avatar
Boris Bonev committed
147
148
149

coeffs = sht(signal)
```
Boris Bonev's avatar
Boris Bonev committed
150

Boris Bonev's avatar
Boris Bonev committed
151
152
`torch_harmonics` also implements a distributed variant of the SHT located in `torch-harmonics.distributed`.

Boris Bonev's avatar
Boris Bonev committed
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
### Cite us

If you use `torch-harmonics` in an academic paper, please cite [1]

```
@misc{bonev2023spherical,
      title={Spherical Fourier Neural Operators: Learning Stable Dynamics on the Sphere}, 
      author={Boris Bonev and Thorsten Kurth and Christian Hundt and Jaideep Pathak and Maximilian Baust and Karthik Kashinath and Anima Anandkumar},
      year={2023},
      eprint={2306.03838},
      archivePrefix={arXiv},
      primaryClass={cs.LG}
}
```

Boris Bonev's avatar
Boris Bonev committed
168
169
170
## References

<a id="1">[1]</a> 
Boris Bonev's avatar
Boris Bonev committed
171
172
173
174
175
176
Bonev B., Kurth T., Hundt C., Pathak, J., Baust M., Kashinath K., Anandkumar A.;
Spherical Fourier Neural Operators: Learning Stable Dynamics on the Sphere;
arXiv 2306.0383, 2023.

<a id="1">[2]</a> 
Schaeffer N.;
Boris Bonev's avatar
Boris Bonev committed
177
Efficient spherical harmonic transforms aimed at pseudospectral numerical simulations;
Boris Bonev's avatar
Boris Bonev committed
178
G3: Geochemistry, Geophysics, Geosystems, 2013.
Boris Bonev's avatar
Boris Bonev committed
179

Boris Bonev's avatar
Boris Bonev committed
180
181
<a id="1">[3]</a> 
Wang B., Wang L., Xie Z.;
Boris Bonev's avatar
Boris Bonev committed
182
Accurate calculation of spherical and vector spherical harmonic expansions via spectral element grids;
Boris Bonev's avatar
Boris Bonev committed
183
Adv Comput Math, 2018.