yt_astro_analysis.halo_analysis.halo_catalog.plot_modifications.HaloCatalogCallback¶
- class yt_astro_analysis.halo_analysis.halo_catalog.plot_modifications.HaloCatalogCallback(halo_catalog, circle_args=None, circle_kwargs=None, width=None, annotate_field=None, radius_field='virial_radius', center_field_prefix='particle_position', text_args=None, font_kwargs=None, factor=1.0)¶
Plots circles at the locations of all the halos in a halo catalog with radii corresponding to the virial radius of each halo.
- Parameters:
halo_catalog (
Dataset
,) –YTDataContainer
, orHaloCatalog
The object containing halos to be overplotted. This can be a HaloCatalog object, a loaded halo catalog dataset, or a data container from a halo catalog dataset.circle_args (list) – Contains the arguments controlling the appearance of the circles, supplied to the Matplotlib patch Circle.
width (tuple) – The width over which to select halos to plot, useful when overplotting to a slice plot. Accepts a tuple in the form (1.0, ‘Mpc’).
annotate_field (str) – A field contained in the halo catalog to add text to the plot near the halo. Example: annotate_field = ‘particle_mass’ will write the halo mass next to each halo.
radius_field (str) – A field contained in the halo catalog to set the radius of the circle which will surround each halo. Default: ‘virial_radius’.
center_field_prefix (str) – Accepts a field prefix which will be used to find the fields containing the coordinates of the center of each halo. Ex: ‘particle_position’ will result in the fields ‘particle_position_x’ for x ‘particle_position_y’ for y, and ‘particle_position_z’ for z. Default: ‘particle_position’.
text_args (dict) – Contains the arguments controlling the text appearance of the annotated field.
factor (float) – A number the virial radius is multiplied by for plotting the circles. Ex: factor = 2.0 will plot circles with twice the radius of each halo virial radius.
Examples
>>> import yt >>> dds = yt.load("Enzo_64/DD0043/data0043") >>> hds = yt.load("rockstar_halos/halos_0.0.bin") >>> p = yt.ProjectionPlot( ... dds, "x", ("gas", "density"), weight_field=("gas", "density") ... ) >>> p.annotate_halos(hds) >>> p.save()
>>> # plot a subset of all halos >>> import yt >>> dds = yt.load("Enzo_64/DD0043/data0043") >>> hds = yt.load("rockstar_halos/halos_0.0.bin") >>> # make a region half the width of the box >>> dregion = dds.box( ... dds.domain_center - 0.25 * dds.domain_width, ... dds.domain_center + 0.25 * dds.domain_width, ... ) >>> hregion = hds.box( ... hds.domain_center - 0.25 * hds.domain_width, ... hds.domain_center + 0.25 * hds.domain_width, ... ) >>> p = yt.ProjectionPlot( ... dds, ... "x", ... ("gas", "density"), ... weight_field=("gas", "density"), ... data_source=dregion, ... width=0.5, ... ) >>> p.annotate_halos(hregion) >>> p.save()
>>> # plot halos from a HaloCatalog >>> import yt >>> from yt.extensions.astro_analysis.halo_analysis import HaloCatalog >>> dds = yt.load("Enzo_64/DD0043/data0043") >>> hds = yt.load("rockstar_halos/halos_0.0.bin") >>> hc = HaloCatalog(data_ds=dds, halos_ds=hds) >>> p = yt.ProjectionPlot( ... dds, "x", ("gas", "density"), weight_field=("gas", "density") ... ) >>> p.annotate_halos(hc) >>> p.save()
|
Plots circles at the locations of all the halos in a halo catalog with radii corresponding to the virial radius of each halo. |