los_response.py 9.17 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
# 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/>.
#
14
# Copyright(C) 2013-2019 Max-Planck-Society
15
#
16
# NIFTy is being developed at the Max-Planck-Institut fuer Astrophysik.
Philipp Arras's avatar
Philipp Arras committed
17

Martin Reinecke's avatar
Martin Reinecke committed
18
19
20
import numpy as np
from scipy.sparse import coo_matrix
from scipy.sparse.linalg import aslinearoperator
Philipp Arras's avatar
Philipp Arras committed
21
22
from scipy.special import erfc

Martin Reinecke's avatar
Martin Reinecke committed
23
24
25
26
from ..domain_tuple import DomainTuple
from ..domains.rg_space import RGSpace
from ..domains.unstructured_domain import UnstructuredDomain
from ..field import Field
Philipp Arras's avatar
Philipp Arras committed
27
from ..operators.linear_operator import LinearOperator
Martin Reinecke's avatar
Martin Reinecke committed
28
29


30
def _gaussian_sf(x):
31
    return 0.5*erfc(x/np.sqrt(2.))
Martin Reinecke's avatar
Martin Reinecke committed
32
33


34
def _comp_traverse(start, end, shp, dist, lo, mid, hi, sig, erf):
Martin Reinecke's avatar
Martin Reinecke committed
35
36
    ndim = start.shape[0]
    nlos = start.shape[1]
37
    inc = np.full(len(shp), 1, dtype=np.int64)
Martin Reinecke's avatar
Martin Reinecke committed
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
    for i in range(-2, -len(shp)-1, -1):
        inc[i] = inc[i+1]*shp[i+1]

    pmax = np.array(shp)

    out = [None]*nlos
    for i in range(nlos):
        dir = end[:, i]-start[:, i]
        dirx = np.where(dir == 0., 1e-12, dir)
        d0 = np.where(dir == 0., ((start[:, i] > 0)-0.5)*1e12,
                      -start[:, i]/dirx)
        d1 = np.where(dir == 0., ((start[:, i] < pmax)-0.5)*-1e12,
                      (pmax-start[:, i])/dirx)
        (dmin, dmax) = (np.minimum(d0, d1), np.maximum(d0, d1))
        dmin = dmin.max()
        dmax = dmax.min()
        dmin = np.maximum(0., dmin)
        dmax = np.minimum(1., dmax)
        dmax = np.maximum(dmin, dmax)
        # hack: move away from potential grid crossings
        dmin += 1e-7
        dmax -= 1e-7
Martin Reinecke's avatar
Martin Reinecke committed
60
        if dmin >= dmax:  # no intersection
61
            out[i] = (np.full(0, 0, dtype=np.int64), np.full(0, 0.))
Martin Reinecke's avatar
Martin Reinecke committed
62
63
64
65
66
67
68
69
70
71
72
73
74
            continue
        # determine coordinates of first cell crossing
        c_first = np.ceil(start[:, i]+dir*dmin)
        c_first = np.where(dir > 0., c_first, c_first-1.)
        c_first = (c_first-start[:, i])/dirx
        pos1 = np.asarray((start[:, i]+dmin*dir), dtype=np.int)
        pos1 = np.sum(pos1*inc)
        cdist = np.empty(0, dtype=np.float64)
        add = np.empty(0, dtype=np.int)
        for j in range(ndim):
            if dir[j] != 0:
                step = inc[j] if dir[j] > 0 else -inc[j]
                tmp = np.arange(start=c_first[j], stop=dmax,
Martin Reinecke's avatar
Martin Reinecke committed
75
                                step=abs(1./dir[j]))
Martin Reinecke's avatar
Martin Reinecke committed
76
                cdist = np.append(cdist, tmp)
77
                add = np.append(add, np.full(len(tmp), step, dtype=np.int64))
Martin Reinecke's avatar
Martin Reinecke committed
78
79
80
        idx = np.argsort(cdist)
        cdist = cdist[idx]
        add = add[idx]
81
82
        cdist = np.append(np.full(1, dmin), cdist)
        cdist = np.append(cdist, np.full(1, dmax))
Martin Reinecke's avatar
Martin Reinecke committed
83
84
85
86
        corfac = np.linalg.norm(dir*dist)
        cdist *= corfac
        wgt = np.diff(cdist)
        mdist = 0.5*(cdist[:-1]+cdist[1:])
87
        wgt = apply_erf(wgt, mdist, lo[i], mid[i], hi[i], sig[i], erf)
Martin Reinecke's avatar
Martin Reinecke committed
88
89
90
91
92
93
        add = np.append(pos1, add)
        add = np.cumsum(add)
        out[i] = (add, wgt)
    return out


94
def apply_erf(wgt, dist, lo, mid, hi, sig, erf):
Martin Reinecke's avatar
Martin Reinecke committed
95
96
97
    wgt = wgt.copy()
    mask = dist > hi
    wgt[mask] = 0.
98
    mask = (dist > lo) & (dist <= hi)
99
    wgt[mask] *= erf((-1/dist[mask]+1/mid)/sig)
Martin Reinecke's avatar
Martin Reinecke committed
100
101
102
103
    return wgt


class LOSResponse(LinearOperator):
Martin Reinecke's avatar
Martin Reinecke committed
104
105
    """Line-of-sight response operator

106
    This operator transforms from a single RGSpace to an UnstructuredDomain
Martin Reinecke's avatar
Martin Reinecke committed
107
108
109
110
111
112
113
114
115
116
117
118
    with as many entries as there were lines of sight passed to the
    constructor. Adjoint application is also provided.

    Parameters
    ----------
    domain : RGSpace or DomainTuple
        The operator's input domain. This must be a single RGSpace.
    starts, ends : numpy.ndarray(float) with two dimensions
        Arrays containing the start and end points of the individual lines
        of sight. The first dimension must have as many entries as `domain`
        has dimensions. The second dimensions must be identical for both arrays
        and indicated the total number of lines of sight.
119
    sigmas: numpy.ndarray(float) (optional)
Martin Reinecke's avatar
Martin Reinecke committed
120
121
122
        If this is not None, the inverse of the lengths of the LOSs are assumed
        to be Gaussian distributed with these sigmas. The start point will
        remain the same, but the endpoint is assumed to be unknown.
123
124
        This is a typical statistical model for astrophysical parallaxes.
        The LOS response then returns the expected integral
Martin Reinecke's avatar
Martin Reinecke committed
125
126
        over the input given that the length of the LOS is unknown and
        therefore the result is averaged over different endpoints.
Philipp Arras's avatar
Philipp Arras committed
127
        Default: None.
128
129
    truncation: float (optional)
        Use only if the sigmas keyword argument is used!
Martin Reinecke's avatar
Martin Reinecke committed
130
131
132
133
        This truncates the probability of the endpoint lying more sigmas away
        than the truncation. Used to speed up computation and to avoid negative
        distances. It should hold that `1./(1./length-sigma*truncation)>0`
        for all lengths of the LOSs and all corresponding sigma of sigmas.
134
        If unsure, leave blank.
Philipp Arras's avatar
Philipp Arras committed
135
        Default: 3.
Martin Reinecke's avatar
Martin Reinecke committed
136
137
138

    Notes
    -----
139
    `starts, `ends`, `sigmas`, and `truncation` have to be identical on
Martin Reinecke's avatar
Martin Reinecke committed
140
141
    every calling MPI task (i.e. the full LOS information has to be provided on
    every task).
Martin Reinecke's avatar
Martin Reinecke committed
142
    """
Philipp Arras's avatar
Philipp Arras committed
143

144
    def __init__(self, domain, starts, ends, sigmas=None, truncation=3.):
Martin Reinecke's avatar
Martin Reinecke committed
145
        self._domain = DomainTuple.make(domain)
Martin Reinecke's avatar
Martin Reinecke committed
146
        self._capability = self.TIMES | self.ADJOINT_TIMES
Martin Reinecke's avatar
Martin Reinecke committed
147
148
149
150
151
152
153
154
155

        if ((not isinstance(self.domain[0], RGSpace)) or
                (len(self._domain) != 1)):
            raise TypeError("The domain must be exactly one RGSpace instance.")

        ndim = len(self.domain[0].shape)
        starts = np.array(starts)
        nlos = starts.shape[1]
        ends = np.array(ends)
156
157
158
        if sigmas is None:
            sigmas = np.zeros(nlos, dtype=np.float32)
        sigmas = np.array(sigmas)
Martin Reinecke's avatar
Martin Reinecke committed
159
160
        if starts.shape[0] != ndim:
            raise TypeError("dimension mismatch")
161
        if nlos != sigmas.shape[0]:
Martin Reinecke's avatar
Martin Reinecke committed
162
163
164
            raise TypeError("dimension mismatch")
        if starts.shape != ends.shape:
            raise TypeError("dimension mismatch")
Martin Reinecke's avatar
Martin Reinecke committed
165
166
167
168

        diffs = ends-starts
        difflen = np.linalg.norm(diffs, axis=0)
        diffs /= difflen
169
        real_distances = 1./(1./difflen - truncation*sigmas)
Martin Reinecke's avatar
Martin Reinecke committed
170
171
172
        if np.any(real_distances < 0):
            raise ValueError("parallax error truncation to high: "
                             "getting negative distances")
173
        real_ends = starts + diffs*real_distances
Martin Reinecke's avatar
Martin Reinecke committed
174
        dist = np.array(self.domain[0].distances).reshape((-1, 1))
Martin Reinecke's avatar
stage 1  
Martin Reinecke committed
175
176
        localized_pixel_starts = starts/dist + 0.5
        localized_pixel_ends = real_ends/dist + 0.5
Martin Reinecke's avatar
Martin Reinecke committed
177
178
179
180
181
182

        # get the shape of the local data slice
        w_i = _comp_traverse(localized_pixel_starts,
                             localized_pixel_ends,
                             self._local_shape,
                             np.array(self.domain[0].distances),
Martin Reinecke's avatar
Martin Reinecke committed
183
184
185
                             1./(1./difflen+truncation*sigmas),
                             difflen,
                             1./(1./difflen-truncation*sigmas),
186
                             sigmas,
187
                             _gaussian_sf)
Martin Reinecke's avatar
Martin Reinecke committed
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212

        boxsz = 16
        nlos = len(w_i)
        npix = np.prod(self._local_shape)
        ntot = 0
        for i in w_i:
            ntot += len(i[1])
        pri = np.empty(ntot, dtype=np.float64)
        ilos = np.empty(ntot, dtype=np.int32)
        iarr = np.empty(ntot, dtype=np.int32)
        xwgt = np.empty(ntot, dtype=np.float32)
        ofs = 0
        cnt = 0
        for i in w_i:
            nval = len(i[1])
            ilos[ofs:ofs+nval] = cnt
            iarr[ofs:ofs+nval] = i[0]
            xwgt[ofs:ofs+nval] = i[1]
            fullidx = np.unravel_index(i[0], self._local_shape)
            tmp = np.zeros(nval, dtype=np.float64)
            fct = 1.
            for j in range(ndim):
                tmp += (fullidx[j]//boxsz)*fct
                fct *= self._local_shape[j]
            tmp += cnt/float(nlos)
213
            tmp += iarr[ofs:ofs+nval]/(float(nlos)*float(npix))
Martin Reinecke's avatar
Martin Reinecke committed
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
            pri[ofs:ofs+nval] = tmp
            ofs += nval
            cnt += 1
        xtmp = np.argsort(pri)
        ilos = ilos[xtmp]
        iarr = iarr[xtmp]
        xwgt = xwgt[xtmp]
        self._smat = aslinearoperator(
            coo_matrix((xwgt, (ilos, iarr)),
                       shape=(nlos, np.prod(self._local_shape))))

        self._target = DomainTuple.make(UnstructuredDomain(nlos))

    def apply(self, x, mode):
        self._check_input(x, mode)
        if mode == self.TIMES:
            result_arr = self._smat.matvec(x.local_data.reshape(-1))
            return Field.from_global_data(self._target, result_arr,
                                          sum_up=True)
        local_input_data = x.to_global_data().reshape(-1)
        res = self._smat.rmatvec(local_input_data).reshape(self._local_shape)
        return Field.from_local_data(self._domain, res)