linear_interpolation.py 5.2 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# This program is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with this program.  If not, see <http://www.gnu.org/licenses/>.
#
# Copyright(C) 2013-2018 Max-Planck-Society
#
# NIFTy is being developed at the Max-Planck-Institut fuer Astrophysik
# and financially supported by the Studienstiftung des deutschen Volkes.

from __future__ import absolute_import, division, print_function

Philipp Arras's avatar
Cleanup    
Philipp Arras committed
21
22
23
24
from numpy import (abs, arange, array, int64, mgrid, prod, ravel,
                   ravel_multi_index, zeros)
from scipy.sparse import coo_matrix
from scipy.sparse.linalg import aslinearoperator
25

Philipp Arras's avatar
Cleanup    
Philipp Arras committed
26
27
28
from ..compat import *
from ..domains.unstructured_domain import UnstructuredDomain
from ..field import Field
29
from ..sugar import makeDomain
30
31
32
33
34
35
36
37
38
39
from .linear_operator import LinearOperator


class LinearInterpolator(LinearOperator):
    def __init__(self, domain, positions):
        """

        :param domain:
            RGSpace
        :param positions:
Philipp Arras's avatar
Cleanup    
Philipp Arras committed
40
            positions at which to interpolate, shape (dim, ndata)
41
42
43
44
45
46
47
48
49
        """
        self._domain = makeDomain(domain)
        N_points = positions.shape[1]
        self._target = makeDomain(UnstructuredDomain(N_points))
        self._capability = self.TIMES | self.ADJOINT_TIMES
        self._build_mat(positions, N_points)

    def _build_mat(self, positions, N_points):
        ndim = positions.shape[0]
50
        mg = mgrid[(slice(0, 2),)*ndim]
51
        mg = array(list(map(ravel, mg)))
52
        dist = array(self.domain[0].distances).reshape((-1, 1))
53
54
55
56
57
58
59
        pos = positions/dist
        excess = pos-pos.astype(int64)
        pos = pos.astype(int64)
        data = zeros((len(mg[0]), N_points))
        ii = zeros((len(mg[0]), N_points), dtype=int64)
        jj = zeros((len(mg[0]), N_points), dtype=int64)
        for i in range(len(mg[0])):
60
61
62
            factor = prod(abs(1-mg[:, i].reshape((-1, 1))-excess), axis=0)
            data[i, :] = factor
            fromi = pos+mg[:, i].reshape((-1, 1))
63
64
            ii[i, :] = arange(N_points)
            jj[i, :] = ravel_multi_index(fromi, self.domain.shape)
65
66
67
        self._mat = coo_matrix((data.reshape(-1),
                               (ii.reshape(-1), jj.reshape(-1))),
                               (N_points, prod(self.domain.shape)))
68
        self._mat = aslinearoperator(self._mat)
69

70
71
72
73
74
75
76
77
78
    def apply(self, x, mode):
        self._check_input(x, mode)
        x_val = x.to_global_data()
        if mode == self.TIMES:
            res = self._mat.matvec(x_val.reshape((-1,)))
            return Field.from_global_data(self.target, res)
        res = self._mat.rmatvec(x_val).reshape(self.domain.shape)
        return Field.from_global_data(self.domain, res)

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

# import numpy as np
# from ..domains.rg_space import RGSpace
# import itertools
#
#
# class LinearInterpolator(LinearOperator):
#     def __init__(self, domain, positions):
#         """
#         :param domain:
#             RGSpace
#         :param target:
#             UnstructuredDomain, shape (ndata,)
#         :param positions:
#             positions at which to interpolate
#             Field with UnstructuredDomain, shape (dim, ndata)
#         """
#         if not isinstance(domain, RGSpace):
#             raise TypeError("RGSpace needed")
#         if np.any(domain.shape < 2):
#             raise ValueError("RGSpace shape too small")
#         if positions.ndim != 2:
#             raise ValueError("positions must be a 2D array")
#         ndim = len(domain.shape)
#         if positions.shape[0] != ndim:
#             raise ValueError("shape mismatch")
#         self._domain = makeDomain(domain)
#         N_points = positions.shape[1]
#         dist = np.array(domain.distances).reshape((ndim, -1))
#         self._pos = positions/dist
#         shp = np.array(domain.shape, dtype=int).reshape((ndim, -1))
#         self._idx = np.maximum(0, np.minimum(shp-2, self._pos.astype(int)))
#         self._pos -= self._idx
#         tmp = tuple([0, 1] for i in range(ndim))
#         self._corners = np.array(list(itertools.product(*tmp)))
#         self._target = makeDomain(UnstructuredDomain(N_points))
#         self._capability = self.TIMES | self.ADJOINT_TIMES
#
#     def apply(self, x, mode):
#         self._check_input(x, mode)
#         x = x.to_global_data()
#         ndim = len(self._domain.shape)
#
#         res = np.zeros(self._tgt(mode).shape, dtype=x.dtype)
#         for corner in self._corners:
#             corner = corner.reshape(ndim, -1)
#             idx = self._idx+corner
#             idx2 = tuple(idx[t, :] for t in range(idx.shape[0]))
#             wgt = np.prod(self._pos*corner+(1-self._pos)*(1-corner), axis=0)
#             if mode == self.TIMES:
#                 res += wgt*x[idx2]
#             else:
#                 np.add.at(res, idx2, wgt*x)
#         return Field.from_global_data(self._tgt(mode), res)