pde_sphere.py 6.29 KB
Newer Older
Boris Bonev's avatar
Boris Bonev committed
1
2
# coding=utf-8

Boris Bonev's avatar
Boris Bonev committed
3
# SPDX-FileCopyrightText: Copyright (c) 2022 The torch-harmonics Authors. All rights reserved.
Boris Bonev's avatar
Boris Bonev committed
4
# 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
32
33
34
# 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.
#


import torch
import torch.nn as nn
Boris Bonev's avatar
Boris Bonev committed
35
import torch_harmonics as th
Thorsten Kurth's avatar
Thorsten Kurth committed
36
from torch_harmonics.quadrature import _precompute_longitudes
Boris Bonev's avatar
Boris Bonev committed
37

Thorsten Kurth's avatar
Thorsten Kurth committed
38
import math
Boris Bonev's avatar
Boris Bonev committed
39
40
41
42
43
44
45
46
47
import numpy as np


class SphereSolver(nn.Module):
    """
    Solver class on the sphere. Can solve the following PDEs:
    - Allen-Cahn eq
    """

48
    def __init__(self, nlat, nlon, dt, lmax=None, mmax=None, grid="equiangular", radius=1.0, coeff=0.001):
Boris Bonev's avatar
Boris Bonev committed
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
        super().__init__()

        # time stepping param
        self.dt = dt

        # grid parameters
        self.nlat = nlat
        self.nlon = nlon
        self.grid = grid

        # physical sonstants
        self.register_buffer('radius', torch.as_tensor(radius, dtype=torch.float64))
        self.register_buffer('coeff', torch.as_tensor(coeff, dtype=torch.float64))

        # SHT
Boris Bonev's avatar
Boris Bonev committed
64
65
        self.sht = th.RealSHT(nlat, nlon, lmax=lmax, mmax=mmax, grid=grid, csphase=False)
        self.isht = th.InverseRealSHT(nlat, nlon, lmax=lmax, mmax=mmax, grid=grid, csphase=False)
Boris Bonev's avatar
Boris Bonev committed
66
67
68
69
70
71

        self.lmax = lmax or self.sht.lmax
        self.mmax = lmax or self.sht.mmax

        # compute gridpoints
        if self.grid == "legendre-gauss":
Boris Bonev's avatar
Boris Bonev committed
72
            cost, _ = th.quadrature.legendre_gauss_weights(self.nlat, -1, 1)
Boris Bonev's avatar
Boris Bonev committed
73
        elif self.grid == "lobatto":
Boris Bonev's avatar
Boris Bonev committed
74
            cost, _ = th.quadrature.lobatto_weights(self.nlat, -1, 1)
Boris Bonev's avatar
Boris Bonev committed
75
        elif self.grid == "equiangular":
Boris Bonev's avatar
Boris Bonev committed
76
            cost, _ = th.quadrature.clenshaw_curtiss_weights(self.nlat, -1, 1)
Boris Bonev's avatar
Boris Bonev committed
77
78

        # apply cosine transform and flip them
Thorsten Kurth's avatar
Thorsten Kurth committed
79
80
        lats = -torch.arcsin(cost)
        lons = _precompute_longitudes(self.nlon)
Boris Bonev's avatar
Boris Bonev committed
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100

        self.lmax = self.sht.lmax
        self.mmax = self.sht.mmax

        l = torch.arange(0, self.lmax).reshape(self.lmax, 1).cdouble()
        l = l.expand(self.lmax, self.mmax)
        # the laplace operator acting on the coefficients is given by l (l + 1)
        lap = - l * (l + 1) / self.radius**2
        invlap = - self.radius**2 / l / (l + 1)
        invlap[0] = 0.

        # register all
        self.register_buffer('lats', lats)
        self.register_buffer('lons', lons)
        self.register_buffer('l', l)
        self.register_buffer('lap', lap)
        self.register_buffer('invlap', invlap)

    def grid2spec(self, u):
        """spectral coefficients from spatial data"""
101

Boris Bonev's avatar
Boris Bonev committed
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
        return self.sht(u)

    def spec2grid(self, uspec):
        """spatial data from spectral coefficients"""

        return self.isht(uspec)

    def dudtspec(self, uspec, pde='allen-cahn'):

        if pde == 'allen-cahn':
            ugrid = self.spec2grid(uspec)
            u3spec  = self.grid2spec(ugrid**3)
            dudtspec = self.coeff*self.lap*uspec + uspec - u3spec
        elif pde == 'ginzburg-landau':
            ugrid = self.spec2grid(uspec)
            u3spec  = self.grid2spec(ugrid**3)
            dudtspec = uspec + (1. + 2.j)*self.coeff*self.lap*uspec - (1. + 2.j)*u3spec
        else:
            NotImplementedError
121

Boris Bonev's avatar
Boris Bonev committed
122
123
124
125
126
127
128
129
130
131
132
133
134
        return dudtspec

    def randspec(self):
        """random data on the sphere"""

        rspec = torch.randn_like(self.lap) / 4 / torch.pi
        return rspec


    def plot_griddata(self, data, fig, cmap='twilight_shifted', vmax=None, vmin=None, projection='3d', title=None, antialiased=False):
        """
        plotting routine for data on the grid. Requires cartopy for 3d plots.
        """
Boris Bonev's avatar
Boris Bonev committed
135
        import matplotlib.pyplot as plt
Boris Bonev's avatar
Boris Bonev committed
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

        lons = self.lons.squeeze() - torch.pi
        lats = self.lats.squeeze()

        if data.is_cuda:
            data = data.cpu()
            lons = lons.cpu()
            lats = lats.cpu()

        Lons, Lats = np.meshgrid(lons, lats)

        if projection == 'mollweide':

            #ax = plt.gca(projection=projection)
            ax = fig.add_subplot(projection=projection)
            im = ax.pcolormesh(Lons, Lats, data, cmap=cmap, vmax=vmax, vmin=vmin)
            # ax.set_title("Elevation map of mars")
            ax.grid(True)
            ax.set_xticklabels([])
            ax.set_yticklabels([])
            plt.colorbar(im, orientation='horizontal')
            plt.title(title)

        elif projection == '3d':

Boris Bonev's avatar
Boris Bonev committed
161
            import cartopy.crs as ccrs
Boris Bonev's avatar
Boris Bonev committed
162

163
            proj = ccrs.Orthographic(central_longitude=0.0, central_latitude=25.0)
Boris Bonev's avatar
Boris Bonev committed
164
165
166

            #ax = plt.gca(projection=proj, frameon=True)
            ax = fig.add_subplot(projection=proj)
Thorsten Kurth's avatar
Thorsten Kurth committed
167
168
            Lons = Lons*180/math.pi
            Lats = Lats*180/math.pi
Boris Bonev's avatar
Boris Bonev committed
169
170
171
172
173
174
175
176
177
178
179

            # contour data over the map.
            im = ax.pcolormesh(Lons, Lats, data, cmap=cmap, transform=ccrs.PlateCarree(), antialiased=antialiased, vmax=vmax, vmin=vmin)
            plt.title(title, y=1.05)

        else:
            raise NotImplementedError

        return im

    def plot_specdata(self, data, fig, **kwargs):
Thorsten Kurth's avatar
Thorsten Kurth committed
180
        return self.plot_griddata(self.isht(data), fig, **kwargs)