combined_source module

Provides a CombinedSource class.

This can be used to join multiple individual sources together for use in a single mock observation.

class martini.sources.combined_source.CombinedSource(sources: list[SPHSource])[source]

Bases: SPHSource

Can combine multiple sources into a single source for a mock observation.

It is often convenient to load source objects using separate source classes, for example perhaps multiple galaxies each loaded with an EAGLESource. However, it would be inefficient to run martini’s source insertion routine multiple times, so it is better to combined the multiple sources into one, first, using this class.

Each source will retain its respective RA, Dec, distance, peculiar velocity, etc. These should all be set when creating the individual sources.

Parameters:

sources (list) – List of source objects to be combined (of type SPHSource or derived classes).

Examples

from astropy import units as U
from martini.souces import SPHSource, CombinedSource
from scipy.spatial.transform import Rotation

n = 1000
r = np.random.rand(n) * 5 * U.kpc
v = 100 * U.km / U.s
t = np.random.rand(n) * 2 * np.pi * U.rad

source1 = SPHSource(
    distance=10 * U.Mpc,
    vpeculiar=30 * U.km / U.s,
    rotation=Rotation.from_euler("y", (30 * U.deg).to_value(U.rad)),
    ra=49.95 * U.deg,
    dec=29.95 * U.deg,
    T_g=np.ones(n) * U.K,
    mHI_g=np.ones(n) * 1e7 * U.Msun,
    xyz_g=U.Quantity(
        [
            np.zeros(n) * U.kpc,
            r * np.cos(t),
            r * np.sin(t),
        ]
    ),
    vxyz_g=U.Quantity(
        [
            np.zeros(n) * U.km / U.s,
            -v * np.sin(t),
            v * np.cos(t),
        ]
    ),
    hsm_g=np.ones(n) * 0.3 * U.kpc,
)
source2 = SPHSource(
    distance=12 * U.Mpc,
    vpeculiar=-50 * U.km / U.s,
    rotation=Rotation.from_euler("z", (60 * U.deg).to_value(U.rad)),
    ra=50.05 * U.deg,
    dec=30.05 * U.deg,
    T_g=np.ones(n) * U.K,
    mHI_g=np.ones(n) * 1e7 * U.Msun,
    xyz_g=U.Quantity(
        [
            np.zeros(n) * U.kpc,
            r * np.cos(t),
            r * np.sin(t),
        ]
    ),
    vxyz_g=U.Quantity(
        [
            np.zeros(n) * U.km / U.s,
            -v * np.sin(t),
            v * np.cos(t),
        ]
    ),
    hsm_g=np.ones(n) * 0.3 * U.kpc,
)
source = CombinedSource([source1, source2])
datacube = DataCube(
    n_px_x=128,
    n_px_y=128,
    n_channels=128,
    px_size=5 * U.arcsec,
    channel_width=3 * U.km / U.s,
    spectral_centre=760 * U.km / U.s,
    ra=50 * U.deg,
    dec=30 * U.deg,
)
beam = GaussianBeam(
    bmaj=15 * U.arcsec,
    bmin=15 * U.arcsec,
)
sph_kernel = GaussianKernel()
spectral_model = GaussianSpectrum(sigma="thermal")
m = Martini(
    source=source,
    datacube=datacube,
    beam=beam,
    noise=None,
    sph_kernel=sph_kernel,
    spectral_model=spectral_model,
)
m.insert_source_in_cube()
m.convolve_beam()
# can write to file as usual
property T_g: Annotated[Quantity, Unit('K')]

Get the combined particle temperatures from all of the combined sources.

Returns:

Quantity, with dimensions of temperature. Particle temperatures.

Return type:

Quantity

apply_mask(mask: ndarray) None[source]

Remove particles from source arrays according to a mask.

Parameters:

mask (ndarray) – Boolean mask. Remove particles with indices corresponding to False values from the source arrays.

boost(boost_vector: Annotated[Quantity, Unit('km / s')]) None[source]

Apply an offset to the source velocity.

Not available for CombinedSource, boost individual sources instead.

Parameters:

boost_vector (Quantity) – Quantity with shape (3, ), with dimensions of velocity. Vector by which to offset the source particle velocities.

Raises:

NotImplementedError – Always raises.

property coordinates_g: CartesianRepresentation

Get the combined particle cartesian coordinates from all of the combined sources.

Returns:

The combined CartesianRepresentation object (including differentials).

Return type:

CartesianRepresentation

property current_rotation: ndarray

Current rotation matrix of the source.

Not available for CombinedSource.

Raises:

NotImplementedError – Always raises.

property distance: Annotated[Quantity, Unit('Mpc')]

The approximate distance of the source.

A CombinedSource has no well-defined distance. This property estimates a distance as the HI mass-weighted mean distance of the combined sources. This has no influence on the output produced by martini, but can influence values that are printed in messages.

Returns:

Quantity, with dimensions of length. The HI mass-weighted mean distance of the combined sources.

Return type:

Quantity

property hsm_g: Annotated[Quantity, Unit('kpc')]

Get the combined particle smoothing lengths from all of the combined sources.

Returns:

Quantity, with dimensions of length. Particle smoothing lengths.

Return type:

Quantity

load_affine_transformations(fname: str) None[source]

Load a set of affine transformation matrices (position and velocity) from a file.

This is not supported for CombinedSource.

Parameters:

fname (str) – File from which to load affine transformation matrices (in *.npy format).

Raises:

NotImplementedError – Always raises.

property mHI_g: Annotated[Quantity, Unit('solMass')]

Get the combined particle HI masses from all of the combined sources.

Returns:

Quantity, with dimensions of mass. Particle HI masses.

Return type:

Quantity

property pixcoords: Annotated[Quantity, Unit('pix')]

Get the combined pixel coordinates from all of the combined sources.

Returns:

Quantity, with dimensions of pixels. 2D coordinates of particles in pixel coordinates.

Return type:

Quantity

preview(max_points: int = 5000, fig: int | Figure = 1, lim: Annotated[Quantity, Unit('deg')] | None = None, vlim: Annotated[Quantity, Unit('km / s')] | None = None, point_scaling: str = 'auto', title: str = '', save: str | None = None) Figure[source]

Produce a figure showing the source particle coordinates and velocities.

Makes a 3-panel figure showing the projection of the combined source as it will appear in the mock observation. The first panel shows the particles in the y-z plane, coloured by the x-component of velocity (MARTINI projects the source along the x-axis). The second and third panels are position-velocity diagrams showing the x-component of velocity against the y and z coordinates, respectively.

Parameters:
  • max_points (int, optional) – Maximum number of points to draw per panel, the particles will be randomly subsampled if the source has more.

  • fig (int or Figure, optional) – Number of the figure in matplotlib, it will be created as plt.figure(fig). Or, an existing figure can be provided.

  • lim (Quantity, optional) – Quantity, with dimensions of length. The coordinate axes extent from -lim to lim. If unspecified, the maximum absolute coordinate of particles in the source is used.

  • vlim (Quantity, optional) – Quantity, with dimensions of speed. The velocity axes and colour bar extend from -vlim to vlim. If unspecified, the maximum absolute velocity of particles in the source is used.

  • point_scaling (str, optional) – By default points are scaled in size and transparency according to their HI mass and the smoothing length (loosely proportional to their surface densities, but with different scaling to achieve a visually useful plot). For some sources the automatic scaling may not give useful results, using point_scaling="fixed" will plot points of constant size without opacity.

  • title (str, optional) – A title for the figure can be provided.

  • save (str, optional) – If provided, the figure is saved using plt.savefig(save). A .png or .pdf suffix is recommended.

Returns:

The preview figure.

Return type:

Figure

rotate(rotation: Rotation | None = None, *, L_coords: L_coords | None = None) Rotation[source]

Rotate the source.

Not available for CombinedSource, rotate individual sources instead.

Parameters:
  • rotation (Rotation, optional) – A Rotation specifying the rotation. This type of object can be initialized from many ways of specifying rotations: rotation matrices, Euler angles, quaternions, etc. Refer to the scipy documentation for details.

  • L_coords (L_coords, optional) – First element containing an inclination, second element an azimuthal angle (both Quantity instances with dimensions of angle). The routine will first attempt to identify a preferred plane based on the angular momenta of the central 1/3 of particles in the source. This plane will then be rotated to lie in the ‘y-z’ plane, followed by a rotation by the azimuthal angle about its angular momentum pole (rotation about ‘x’), and then inclined (rotation about ‘y’). By default the position angle on the sky is 270 degrees, but if a third element is provided it sets the position angle (second rotation about ‘x’).

Raises:

NotImplementedError – Always raises.

save_current_affine_transformations(fname: str) None[source]

Output current affine transformation matrices (position and velocity) to a file.

This is not supported for CombinedSource.

Parameters:

fname (str) – File in which to save affine transformation matrices (in *.npy format).

Raises:

NotImplementedError – Always raises.

save_current_rotation(fname: str) None[source]

Output current rotation matrix to file.

Parameters:

fname (str) – File in which to save rotation matrix (as a text file).

Raises:

NotImplementedError – Always raises.

property skycoords: SkyCoord

Get the combined sky coordinates from all of the combined sources.

Returns:

The combined SkyCoord object.

Return type:

SkyCoord

property spectralcoords: SpectralCoord

Get the combined spectral coordinates from all of the combined sources.

Returns:

The combined SpectralCoord object.

Return type:

SpectralCoord

translate(translation_vector: Annotated[Quantity, Unit('kpc')]) None[source]

Translate the source.

Not available for CombinedSource, translate individual sources instead.

Parameters:

translation_vector (Quantity) – Quantity with shape (3, ), with dimensions of length. Vector by which to offset the source particle coordinates.

Raises:

NotImplementedError – Always raises.

property vsys: Annotated[Quantity, Unit('km / s')]

The approximate systemic velocity of the source.

A CombinedSource has no well-defined systemic velocity. This property estimates a systemic velocity as the HI mass-weighted systemic velocity of the combined sources. This has no influence on the output produced by martini, but can influence values that are printed in messages.

Returns:

Quantity, with dimensions of speed. The HI mass-weighted mean systemic velocity of the combined sources.

Return type:

Quantity