/
transform.py
481 lines (388 loc) · 15.8 KB
/
transform.py
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
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
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
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
133
134
135
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
161
162
163
164
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
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
import numpy
import dask.array as da
from six import string_types
from nbodykit.utils import deprecate
from nbodykit import _global_options
def StackColumns(*cols):
"""
Stack the input dask arrays vertically, column by column.
This uses :func:`dask.array.vstack`.
Parameters
----------
*cols : :class:`dask.array.Array`
the dask arrays to stack vertically together
Returns
-------
:class:`dask.array.Array` :
the dask array where columns correspond to the input arrays
Raises
------
TypeError
If the input columns are not dask arrays
"""
cols = da.broadcast_arrays(*cols)
return da.vstack(cols).T
def ConcatenateSources(*sources, **kwargs):
"""
Concatenate CatalogSource objects together, optionally including only
certain columns in the returned source.
.. note::
The returned catalog object carries the meta-data from only
the first catalog supplied to this function (in the ``attrs`` dict).
Parameters
----------
*sources : subclass of :class:`~nbodykit.base.catalog.CatalogSource`
the catalog source objects to concatenate together
columns : str, list of str, optional
the columns to include in the concatenated catalog
Returns
-------
CatalogSource :
the concatenated catalog source object
Examples
--------
>>> from nbodykit.lab import *
>>> source1 = UniformCatalog(nbar=100, BoxSize=1.0)
>>> source2 = UniformCatalog(nbar=100, BoxSize=1.0)
>>> print(source1.csize, source2.csize)
>>> combined = transform.ConcatenateSources(source1, source2, columns=['Position', 'Velocity'])
>>> print(combined.csize)
"""
from nbodykit.base.catalog import CatalogSource
columns = kwargs.get('columns', None)
if isinstance(columns, string_types):
columns = [columns]
# concatenate all columns, if none provided
if columns is None or columns == []:
columns = sources[0].columns
# check comms
if not all(src.comm == sources[0].comm for src in sources):
raise ValueError("cannot concatenate sources: comm mismatch")
# check all columns are there
for source in sources:
if not all(col in source for col in columns):
raise ValueError(("cannot concatenate sources: columns are missing "
"from some sources"))
# the total size
size = numpy.sum([src.size for src in sources], dtype='intp')
data = {}
for col in columns:
data[col] = da.concatenate([src[col] for src in sources], axis=0)
toret = CatalogSource._from_columns(size, sources[0].comm, **data)
toret.attrs.update(sources[0].attrs)
return toret
def ConstantArray(value, size, chunks=100000):
"""
Return a dask array of the specified ``size`` holding a single value.
This uses numpy's "stride tricks" to avoid replicating
the data in memory for each element of the array.
Parameters
----------
value : float
the scalar value to fill the array with
size : int
the length of the returned dask array
chunks : int, optional
the size of the dask array chunks
"""
ele = numpy.array(value)
toret = numpy.lib.stride_tricks.as_strided(ele, [size] + list(ele.shape), [0] + list(ele.strides))
return da.from_array(toret, chunks=chunks)
def CartesianToEquatorial(pos, observer=[0,0,0], frame='icrs'):
"""
Convert Cartesian position coordinates to equatorial right ascension
and declination, using the specified observer location.
.. note::
RA and DEC will be returned in degrees, with RA in the range [0,360]
and DEC in the range [-90, 90].
Parameters
----------
pos : array_like
a N x 3 array holding the Cartesian position coordinates
observer : array_like
a length 3 array holding the observer location
frame : string
A string, 'icrs' or 'galactic'. The frame of the input position.
Use 'icrs' if the cartesian position is already in Equatorial.
Returns
-------
ra, dec : array_like
the right ascension and declination coordinates, in degrees. RA
will be in the range [0,360] and DEC in the range [-90, 90]
"""
if isinstance(pos, da.Array):
pos = da.rechunk(pos, chunks=_global_options['dask_chunk_size'])
pos, observer = da.broadcast_arrays(pos, observer)
# recenter based on observer
pos = pos - observer
if frame == 'icrs':
# from equatorial to equatorial
s = da.hypot(pos[:,0], pos[:,1])
lon = da.arctan2(pos[:,1], pos[:,0])
lat = da.arctan2(pos[:,2], s)
# convert to degrees
lon = da.rad2deg(lon)
lat = da.rad2deg(lat)
# wrap lon to [0,360]
lon = da.mod(lon-360., 360.)
ra, dec = lon, lat
else:
from astropy.coordinates import SkyCoord
def convert_coord(pos):
x, y, z = pos.T
try:
sc = SkyCoord(x, y, z, representation_type='cartesian', frame=frame)
except:
sc = SkyCoord(x, y, z, representation='cartesian', frame=frame)
scg = sc.transform_to(frame='icrs')
scg.representation = 'unitspherical'
ra, dec = scg.ra.value, scg.dec.value
# must preserve the shape.
ang = numpy.stack([ra, dec], axis=1)
return ang
ang = da.map_blocks(convert_coord, pos, dtype=pos.dtype, chunks=(pos.chunks[0], (2,)))
ra, dec = ang.T
return da.stack((ra, dec), axis=0)
def CartesianToSky(pos, cosmo, velocity=None, observer=[0,0,0], zmax=100., frame='icrs'):
r"""
Convert Cartesian position coordinates to RA/Dec and redshift,
using the specified cosmology to convert radial distances from
the origin into redshift.
If velocity is supplied, the returned redshift accounts for the
additional peculiar velocity shift.
Users should ensure that ``zmax`` is larger than the largest possible
redshift being considered to avoid an interpolation exception.
.. note::
Cartesian coordinates should be in units of Mpc/h and velocity
should be in units of km/s.
Parameters
----------
pos : dask array
a N x 3 array holding the Cartesian position coordinates in Mpc/h
cosmo : :class:`~nbodykit.cosmology.cosmology.Cosmology`
the cosmology used to meausre the comoving distance from ``redshift``
velocity : array_like
a N x 3 array holding velocity in km/s
observer : array_like, optional
a length 3 array holding the observer location
zmax : float, optional
the maximum possible redshift, should be set to a reasonably large
value to avoid interpolation failure going from comoving distance
to redshift
frame : string ('icrs' or 'galactic')
speciefies which frame the Cartesian coordinates is.
Returns
-------
ra, dec, z : dask array
the right ascension (in degrees), declination (in degrees), and
redshift coordinates. RA will be in the range [0,360] and DEC in the
range [-90, 90]
Notes
-----
If velocity is provided, redshift-space distortions are added to the
real-space redshift :math:`z_\mathrm{real}`, via:
.. math::
z_\mathrm{redshift} = ( v_\mathrm{pec} / c ) (1 + z_\mathrm{reals})
Raises
------
TypeError
If the input columns are not dask arrays
"""
from astropy.constants import c
from scipy.interpolate import interp1d
if isinstance(pos, da.Array):
pos = da.rechunk(pos, chunks=_global_options['dask_chunk_size'])
pos, observer = da.broadcast_arrays(pos, observer)
# recenter position
pos = pos - observer
# RA,dec coordinates (in degrees)
ra, dec = CartesianToEquatorial(pos, frame=frame)
# the distance from the origin
r = da.linalg.norm(pos, axis=-1)
def z_from_comoving_distance(x):
zgrid = numpy.logspace(-8, numpy.log10(zmax), 1024)
zgrid = numpy.concatenate([[0.], zgrid])
rgrid = cosmo.comoving_distance(zgrid)
return interp1d(rgrid, zgrid)(x)
# invert distance - redshift relation
z = r.map_blocks(z_from_comoving_distance)
# add in velocity offsets?
if velocity is not None:
vpec = (pos*velocity).sum(axis=-1) / r
z += vpec / c.to('km/s').value * (1 + z)
return da.stack((ra, dec, z), axis=0)
def SkyToUnitSphere(ra, dec, degrees=True, frame='icrs'):
"""
Convert sky coordinates (``ra``, ``dec``) to Cartesian coordinates on
the unit sphere.
Parameters
----------
ra : :class:`dask.array.Array`; shape: (N,)
the right ascension angular coordinate
dec : :class:`dask.array.Array`; ; shape: (N,)
the declination angular coordinate
degrees : bool, optional
specifies whether ``ra`` and ``dec`` are in degrees or radians
frame : string ('icrs' or 'galactic')
speciefies which frame the Cartesian coordinates is.
Returns
-------
pos : :class:`dask.array.Array`; shape: (N,3)
the cartesian position coordinates, where columns represent
``x``, ``y``, and ``z``
Raises
------
TypeError
If the input columns are not dask arrays
"""
ra, dec = da.broadcast_arrays(ra, dec)
if frame == 'icrs':
# no frame transformation
# put into radians from degrees
if degrees:
ra = da.deg2rad(ra)
dec = da.deg2rad(dec)
# cartesian coordinates
x = da.cos( dec ) * da.cos( ra )
y = da.cos( dec ) * da.sin( ra )
z = da.sin( dec )
return da.vstack([x,y,z]).T
else:
from astropy.coordinates import SkyCoord
if degrees:
ra = da.deg2rad(ra)
dec = da.deg2rad(dec)
def convert_coord(ra, dec):
try:
sc = SkyCoord(ra[:, 0], dec[:, 0], unit='rad', representation_type='unitspherical', frame='icrs')
except:
sc = SkyCoord(ra[:, 0], dec[:, 0], unit='rad', representation='unitspherical', frame='icrs')
scg = sc.transform_to(frame=frame)
scg = scg.cartesian
x, y, z = scg.x.value, scg.y.value, scg.z.value
return numpy.stack([x, y, z], axis=1)
return da.map_blocks(convert_coord, ra[:, None], dec[:, None], dtype=ra.dtype, chunks=(ra.chunks[0], (3,)))
def SkyToCartesian(ra, dec, redshift, cosmo, degrees=True, frame='icrs'):
"""
Convert sky coordinates (``ra``, ``dec``, ``redshift``) to a
Cartesian ``Position`` column.
.. warning::
The returned Cartesian position is in units of Mpc/h.
Parameters
-----------
ra : :class:`dask.array.Array`; shape: (N,)
the right ascension angular coordinate
dec : :class:`dask.array.Array`; shape: (N,)
the declination angular coordinate
redshift : :class:`dask.array.Array`; shape: (N,)
the redshift coordinate
cosmo : :class:`~nbodykit.cosmology.cosmology.Cosmology`
the cosmology used to meausre the comoving distance from ``redshift``
degrees : bool, optional
specifies whether ``ra`` and ``dec`` are in degrees
frame : string ('icrs' or 'galactic')
speciefies which frame the Cartesian coordinates is.
Returns
-------
pos : :class:`dask.array.Array`; shape: (N,3)
the cartesian position coordinates, where columns represent
``x``, ``y``, and ``z`` in units of Mpc/h
Raises
------
TypeError
If the input columns are not dask arrays
"""
ra, dec, redshift = da.broadcast_arrays(ra, dec, redshift)
# pos on the unit sphere
pos = SkyToUnitSphere(ra, dec, degrees=degrees, frame=frame)
# multiply by the comoving distance in Mpc/h
r = redshift.map_blocks(lambda z: cosmo.comoving_distance(z), dtype=redshift.dtype)
return r[:,None] * pos
def HaloConcentration(mass, cosmo, redshift, mdef='vir'):
"""
Return halo concentration from halo mass, based on the analytic fitting
formulas presented in
`Dutton and Maccio 2014 <https://arxiv.org/abs/1402.7073>`_.
.. note::
The units of the input mass are assumed to be :math:`M_{\odot}/h`
Parameters
----------
mass : array_like
either a numpy or dask array specifying the halo mass; units
assumed to be :math:`M_{\odot}/h`
cosmo : :class:`~nbodykit.cosmology.cosmology.Cosmology`
the cosmology instance used in the analytic formula
redshift : float
compute the c(M) relation at this redshift
mdef : str, optional
string specifying the halo mass definition to use; should be
'vir' or 'XXXc' or 'XXXm' where 'XXX' is an int specifying the
overdensity
Returns
-------
concen : :class:`dask.array.Array`
a dask array holding the analytic concentration values
References
----------
Dutton and Maccio, "Cold dark matter haloes in the Planck era: evolution
of structural parameters for Einasto and NFW profiles", 2014, arxiv:1402.7073
"""
from halotools.empirical_models import NFWProfile
mass, redshift = da.broadcast_arrays(mass, redshift)
kws = {'cosmology':cosmo.to_astropy(), 'conc_mass_model':'dutton_maccio14', 'mdef':mdef}
def get_nfw_conc(mass, redshift):
kw1 = {}
kw1.update(kws)
kw1['redshift'] = redshift
model = NFWProfile(**kw1)
return model.conc_NFWmodel(prim_haloprop=mass)
return da.map_blocks(get_nfw_conc, mass, redshift, dtype=mass.dtype)
def HaloVelocityDispersion(mass, cosmo, redshift, mdef='vir'):
""" Compute the velocity dispersion of halo from Mass.
This is a simple model suggested by Martin White.
See http://adsabs.harvard.edu/abs/2008ApJ...672..122E
"""
mass, redshift = da.broadcast_arrays(mass, redshift)
def compute_vdisp(mass, redshift):
h = cosmo.efunc(redshift)
return 1100. * (h * mass / 1e15) ** 0.33333
return da.map_blocks(compute_vdisp, mass, redshift, dtype=mass.dtype)
def HaloRadius(mass, cosmo, redshift, mdef='vir'):
r"""
Return halo radius from halo mass, based on the specified mass definition.
This is independent of halo profile, and simply returns
.. math::
R = \left [ 3 M /(4\pi\Delta) \right]^{1/3}
where :math:`\Delta` is the density threshold, which depends on cosmology,
redshift, and mass definition
.. note::
The units of the input mass are assumed to be :math:`M_{\odot}/h`
Parameters
----------
mass : array_like
either a numpy or dask array specifying the halo mass; units
assumed to be :math:`M_{\odot}/h`
cosmo : :class:`~nbodykit.cosmology.cosmology.Cosmology`
the cosmology instance
redshift : float
compute the density threshold which determines the R(M) relation
at this redshift
mdef : str, optional
string specifying the halo mass definition to use; should be
'vir' or 'XXXc' or 'XXXm' where 'XXX' is an int specifying the
overdensity
Returns
-------
radius : :class:`dask.array.Array`
a dask array holding the halo radius
"""
from halotools.empirical_models import halo_mass_to_halo_radius
mass, redshift = da.broadcast_arrays(mass, redshift)
kws = {'cosmology':cosmo.to_astropy(), 'mdef':mdef}
def mass_to_radius(mass, redshift):
return halo_mass_to_halo_radius(mass=mass, redshift=redshift, **kws)
return da.map_blocks(mass_to_radius, mass, redshift, dtype=mass.dtype)
# deprecated functions
vstack = deprecate("nbodykit.transform.vstack", StackColumns, "nbodykit.transform.StackColumns")
concatenate = deprecate("nbodykit.transform.concatenate", ConcatenateSources, "nbodykit.transform.ConcatenateSources")
SkyToCartesion = deprecate("nbodykit.transform.SkyToCartesion", SkyToCartesian, "nbodykit.transform.SkyToCartesian")