los_response.py 8.94 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
# 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.

19
from __future__ import absolute_import, division, print_function
Philipp Arras's avatar
Philipp Arras committed
20

Martin Reinecke's avatar
Martin Reinecke committed
21
22
23
import numpy as np
from scipy.sparse import coo_matrix
from scipy.sparse.linalg import aslinearoperator
Philipp Arras's avatar
Philipp Arras committed
24
25
26
27
from scipy.special import erfc

from .. import dobj
from ..compat import *
Martin Reinecke's avatar
Martin Reinecke committed
28
29
30
31
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
32
from ..operators.linear_operator import LinearOperator
Martin Reinecke's avatar
Martin Reinecke committed
33
34
35


def _gaussian_error_function(x):
36
    return 0.5*erfc(x/np.sqrt(2.))
Martin Reinecke's avatar
Martin Reinecke committed
37
38
39
40

def _comp_traverse(start, end, shp, dist, lo, mid, hi, erf):
    ndim = start.shape[0]
    nlos = start.shape[1]
41
    inc = np.full(len(shp), 1, dtype=np.int64)
Martin Reinecke's avatar
Martin Reinecke committed
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
    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
64
        if dmin >= dmax:  # no intersection
65
            out[i] = (np.full(0, 0, dtype=np.int64), np.full(0, 0.))
Martin Reinecke's avatar
Martin Reinecke committed
66
67
68
69
70
71
72
73
74
75
76
77
78
            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
79
                                step=abs(1./dir[j]))
Martin Reinecke's avatar
Martin Reinecke committed
80
                cdist = np.append(cdist, tmp)
81
                add = np.append(add, np.full(len(tmp), step, dtype=np.int64))
Martin Reinecke's avatar
Martin Reinecke committed
82
83
84
        idx = np.argsort(cdist)
        cdist = cdist[idx]
        add = add[idx]
85
86
        cdist = np.append(np.full(1, dmin), cdist)
        cdist = np.append(cdist, np.full(1, dmax))
Martin Reinecke's avatar
Martin Reinecke committed
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
        corfac = np.linalg.norm(dir*dist)
        cdist *= corfac
        wgt = np.diff(cdist)
        mdist = 0.5*(cdist[:-1]+cdist[1:])
        wgt = apply_erf(wgt, mdist, lo[i], mid[i], hi[i], erf)
        add = np.append(pos1, add)
        add = np.cumsum(add)
        out[i] = (add, wgt)
    return out


def apply_erf(wgt, dist, lo, mid, hi, erf):
    wgt = wgt.copy()
    mask = dist > hi
    wgt[mask] = 0.
102
103
    mask = (dist > lo) & (dist <= hi)
    sig = (1/lo - 1/mid)/3
104
    wgt[mask] *= erf((-1/dist[mask]+1/mid)/sig)
Martin Reinecke's avatar
Martin Reinecke committed
105
106
107
108
    return wgt


class LOSResponse(LinearOperator):
Martin Reinecke's avatar
Martin Reinecke committed
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
    """Line-of-sight response operator

    This operator transforms from a single RGSpace to an unstructured domain
    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.
    sigmas_low, sigmas_up : numpy.ndarray(float) (optional)
125
        sigmas_low is 1/parallax-1/(parallax+3*parallax_error), where the parallax
126
127
        error is assumed to be Gaussian distributed.
        sigmas_up is the distance at which the weight is truncated.
128
129
        Should be at least 1/(parallax-3*parallax_error)-1/parallax, 
        but could be higher.
130
        If unsure, leave blank.
Martin Reinecke's avatar
Martin Reinecke committed
131
132
133
134
135
136

    Notes
    -----
    `starts, `ends`, `sigmas_low`, and `sigmas_up` have to be identical on
    every calling MPI task (i.e. the full LOS information has to be provided on
    every task).
Martin Reinecke's avatar
Martin Reinecke committed
137
    """
Philipp Arras's avatar
Philipp Arras committed
138

Martin Reinecke's avatar
Martin Reinecke committed
139
140
    def __init__(self, domain, starts, ends, sigmas_low=None, sigmas_up=None):
        self._domain = DomainTuple.make(domain)
Martin Reinecke's avatar
Martin Reinecke committed
141
        self._capability = self.TIMES | self.ADJOINT_TIMES
Martin Reinecke's avatar
Martin Reinecke committed
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156

        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)
        if sigmas_low is None:
            sigmas_low = np.zeros(nlos, dtype=np.float32)
        if sigmas_up is None:
            sigmas_up = np.zeros(nlos, dtype=np.float32)
        sigmas_low = np.array(sigmas_low)
        sigmas_up = np.array(sigmas_up)
Martin Reinecke's avatar
Martin Reinecke committed
157
158
159
160
161
162
163
164
        if starts.shape[0] != ndim:
            raise TypeError("dimension mismatch")
        if nlos != sigmas_low.shape[0]:
            raise TypeError("dimension mismatch")
        if starts.shape != ends.shape:
            raise TypeError("dimension mismatch")
        if sigmas_low.shape != sigmas_up.shape:
            raise TypeError("dimension mismatch")
Martin Reinecke's avatar
Martin Reinecke committed
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211

        self._local_shape = dobj.local_shape(self.domain[0].shape)
        local_zero_point = (np.array(
            dobj.ibegin_from_shape(self.domain[0].shape)) *
            np.array(self.domain[0].distances))

        diffs = ends-starts
        difflen = np.linalg.norm(diffs, axis=0)
        diffs /= difflen
        real_ends = ends + sigmas_up*diffs
        lzp = local_zero_point.reshape((-1, 1))
        dist = np.array(self.domain[0].distances).reshape((-1, 1))
        localized_pixel_starts = (starts-lzp)/dist + 0.5
        localized_pixel_ends = (real_ends-lzp)/dist + 0.5

        # 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),
                             difflen-sigmas_low, difflen, difflen+sigmas_up,
                             _gaussian_error_function)

        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)
212
            tmp += iarr[ofs:ofs+nval]/(float(nlos)*float(npix))
Martin Reinecke's avatar
Martin Reinecke committed
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
            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)