sph_source module

Provides the generic SPHSource class.

Enables working with any SPH or similar simulation as input.

class martini.sources.sph_source.SPHSource(*, distance: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("Mpc")], vpeculiar: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("km / s")] = <Quantity 0. km / s>, rotation: ~scipy.spatial.transform._rotation.Rotation | None = None, L_coords: ~martini.L_coords.L_coords | None = None, ra: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("deg")] = <Quantity 0. deg>, dec: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("deg")] = <Quantity 0. deg>, h: float = 0.7, T_g: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("K")] | None = None, mHI_g: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("solMass")], xyz_g: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("kpc")], vxyz_g: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("km / s")], hsm_g: ~typing.Annotated[~astropy.units.quantity.Quantity, Unit("kpc")] | None = None, coordinate_axis: int | None = None, coordinate_frame: BaseRADecFrame = <ICRS Frame>)[source]

Bases: object

Class abstracting HI emission sources consisting of SPH simulation particles.

This class constructs an HI emission source from arrays of SPH particle properties: mass, smoothing length, temperature, position, and velocity.

Parameters:
  • distance (Quantity) – Quantity, with dimensions of length. Source distance, also used to set the velocity offset via Hubble’s law.

  • vpeculiar (Quantity, optional) – Quantity, with dimensions of velocity. Source peculiar velocity along the direction to the source centre.

  • rotation (Rotation, optional) – A rotation to apply to the source particles, specified using the Rotation class. That class supports many ways to specify a rotation (Euler angle, rotation matrices, quaternions, etc.). Refer to the scipy documentation for details. Note that the y-z plane will be the one eventually placed in the plane of the “sky”. Cannot be used at the same time as L_coords.

  • L_coords (L_coords, optional) – A named tuple specifying 3 angles. Import it as from martini import L_coords. The angles are used to orient the galaxy relative to its angular momentum vector, “L”. The routine will first identify a preferred plane based on the angular momenta of the central 1/3 of HI gas. This plane will then be rotated to lie in the plane of the “sky” (y-z plane), rotated by an angle az_rot around the angular momentum vector (rotation around x), then inclined by incl towards or away from the line of sight (rotation around y) and finally rotated on the sky to set the position angle pa (second rotation around x). All rotations are extrinsic. The position angle refers to the receding side of the galaxy measured East of North. The angles should be specified using syntax like: L_coords=L_coords(incl=0 * U.deg, pa=270 * U.deg, az_rot=0 * U.deg). These example values are the defaults. Cannot be used at the same time as rotation.

  • ra (Quantity, optional) – Quantity, with dimensions of angle. Right ascension for the source centroid.

  • dec (Quantity, optional) – Quantity, with dimensions of angle. Declination for the source centroid.

  • h (float, optional) – Dimensionless hubble constant, \(H_0 = h (100\\,\\mathrm{km}\\,\\mathrm{s}^{-1}\\,\\mathrm{Mpc}^{-1})\).

  • T_g (Quantity, optional) – Quantity, with dimensions of temperature. Particle temperature.

  • mHI_g (Quantity) – Quantity, with dimensions of mass. Particle HI mass.

  • xyz_g (Quantity) – Quantity , with dimensions of length. Particle position offset from source centroid. Note that the ‘y-z’ plane is that eventually placed in the plane of the “sky”; ‘x’ is the axis corresponding to the “line of sight”.

  • vxyz_g (Quantity) – Quantity, with dimensions of velocity. Particle velocity offset from source centroid. Note that the ‘y-z’ plane is that eventually placed in the plane of the “sky”; ‘x’ is the axis corresponding to the “line of sight”.

  • hsm_g (Quantity, optional) – Quantity, with dimensions of length. Particle SPH smoothing lengths, defined as the FWHM of the smoothing kernel. Smoothing lengths are variously defined in the literature as the radius where the kernel amplitude reaches 0, or some rational fraction of this radius (and other definitions may well exist). The FWHM requested here is not a standard choice (with the exception of SWIFT snapshots!), but has the advantage of avoiding ambiguity in the definition.

  • coordinate_axis (int, optional) – Rank of axis corresponding to position or velocity of a single particle. I.e. coordinate_axis=0 if shape is (3, N), or 1 if (N, 3). Usually prefer to omit this as it can be determined automatically, but is ambiguous for sources with exactly 3 particles.

  • coordinate_frame (BaseRADecFrame) – Optional. The coordinate frame assumed in converting particle coordinates to RA and Dec, and for transforming coordinates and velocities to the data cube frame. The frame needs to have a well-defined velocity as well as spatial origin. Recommended frames are GCRS, ICRS, HCRS, LSRK, LSRD or LSR. The frame should be passed initialized, e.g. ICRS() (not just ICRS).

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.

Note that the “line of sight” is along the ‘x’ axis.

Parameters:

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

property current_rotation: ndarray

Current rotation matrix of the source.

Returns:

The rotation matrix taking the coordinate originally passed in to the source to the current orientation.

Return type:

ndarray

See also

save_current_rotation, save_current_affine_transformations, load_current_affine_transformations

load_affine_transformations(fname: str) None[source]

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

These matrices encode an arbitrary combination of translations and rotations. When loaded they are applied to the source to transform it from its current orientation. The state is only the one that the source had when the files were saved if it has not been transformed since being loaded into martini. The rotation part of these matrices is dimensionless, while the translation part assumes dimensions of Mpc for positions or km/s for velocities (but is encoded without units).

An affine transformation matrix looks like:

\[\begin{split}A = \left[ {\begin{array}{cccc} R_{00} & R_{01} & R_{02} & t_{0} \\ R_{10} & R_{11} & R_{12} & t_{1} \\ R_{20} & R_{21} & R_{22} & t_{2} \\ 0 & 0 & 0 & 1 \\ \end{array} } \right]\end{split}\]

Where \(R\) is a rotation matrix and \(t\) is a translation vector (either in position or velocity). The affine transformation applied to a vector \(x\) is equivalent to \(Rx+t\). The loaded file should contain two affine transformation matrices (in order: one for position, one for velocity) as an array with shape (2, 4, 4). The rotational part of the two matrices should be identical. The translation parts should assume implicit units of Mpc and km/s. This is the format saved by save_current_affine_transformations().

Paramters

fnamestr

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

raises ValueError:

If the rotation part of the two matrices is not identical: this would result in an inconsistent source state.

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 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 angle. The coordinate axes extend from -lim to lim. If unspecified, the maximum angular separation between the source centre defined by its ra and dec and the particle sky positions 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.

The arguments correspond to different rotation types. Multiple types cannot be given in a single function call.

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’).

Returns:

A Rotation specifying the actual rotation done in this function call.

Return type:

Rotation

save_current_affine_transformations(fname: str) None[source]

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

These matrices encode the arbitrary combination of translations and rotations that takes the source from its initial orientation when loaded into martini to its current orientation. The rotation part is dimensionless, while the translation part assumes dimensions of Mpc for positions or km/s for velocities (but is encoded without units).

An affine transformation matrix looks like:

\[\begin{split}A = \left[ {\begin{array}{cccc} R_{00} & R_{01} & R_{02} & t_{0} \\ R_{10} & R_{11} & R_{12} & t_{1} \\ R_{20} & R_{21} & R_{22} & t_{2} \\ 0 & 0 & 0 & 1 \\ \end{array} } \right]\end{split}\]

Where \(R\) is a rotation matrix and \(t\) is a translation vector (either in position or velocity). The affine transformation applied to a vector \(x\) is equivalent to \(Rx+t\). The saved file contains two affine transformation matrices (in order: one for position, one for velocity) as an array with shape (2, 4, 4). The rotational part of the two matrices is identical.

Parameters:

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

save_current_rotation(fname: str) None[source]

Output current rotation matrix to file.

The rotation matrix can be applied to astropy coordinates (e.g. a CartesianRepresentation) as coordinates.transform(np.loadtxt(fname)). If you want to save the full coordinate transformation state of a martini source (and optionally re-load it later), see save_current_affine_transformations() and load_affine_transformations().

Parameters:

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

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

Translate the source.

Note that the “line of sight” is along the ‘x’ axis.

Parameters:

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