After a plot is generated using the standard tools (e.g. SlicePlot,
ProjectionPlot, etc.), it can be annotated with any number of callbacks
before being saved to disk. These callbacks can modify the plots by adding
lines, text, markers, streamlines, velocity vectors, contours, and more.
Callbacks can be applied to plots created with
SlicePlot
,
ProjectionPlot
,
OffAxisSlicePlot
, or
OffAxisProjectionPlot
by calling
one of the annotate_
methods that hang off of the plot object.
The annotate_
methods are dynamically generated based on the list
of available callbacks. For example:
slc = SlicePlot(ds,0,'density')
slc.annotate_title('This is a Density plot')
would add the TitleCallback()
to
the plot object. All of the callbacks listed below are available via
similar annotate_
functions.
To clear one or more annotations from an existing plot, see the annotate_clear() function.
For a brief demonstration of a few of these callbacks in action together, see the cookbook recipe: Annotating Plots to Include Lines, Text, Shapes, etc..
Many of the callbacks (e.g.
TextLabelCallback
) are specified
to occur at user-defined coordinate locations (like where to place a marker
or text on the plot). There are several different coordinate systems used
to identify these locations. These coordinate systems can be specified with
the coord_system
keyword in the relevant callback, which is by default
set to data
. The valid coordinate systems are:
data
– the 3D dataset coordinates
plot
– the 2D coordinates defined by the actual plot limits
axis
– the MPL axis coordinates: (0,0) is lower left; (1,1) is upper right
figure
– the MPL figure coordinates: (0,0) is lower left, (1,1) is upper right
Here we will demonstrate these different coordinate systems for an projection of the x-plane (i.e. with axes in the y and z directions):
The underlying functions are more thoroughly documented in Callback List.
annotate_clear
(index=None)¶This function will clear previous annotations (callbacks) in the plot. If no index is provided, it will clear all annotations to the plot. If an index is provided, it will clear only the Nth annotation to the plot. Note that the index goes from 0..N, and you can specify the index of the last added annotation as -1.
annotate_arrow
(self, pos, length=0.03, coord_system='data', plot_args=None)¶(This is a proxy for
ArrowCallback
.)
Overplot an arrow pointing at a position for highlighting a specific feature. Arrow points from lower left to the designated position with arrow length “length”.
annotate_clumps
(self, clumps, plot_args=None)¶(This is a proxy for
ClumpContourCallback
.)
Take a list of clumps
and plot them as a set of
contours.
annotate_contour
(self, field, ncont=5, factor=4, take_log=False, clim=None, plot_args=None, label=False, text_args=None, data_source=None)¶(This is a proxy for
ContourCallback
.)
Add contours in field
to the plot. ncont
governs the number of
contours generated, factor
governs the number of points used in the
interpolation, take_log
governs how it is contoured and clim
gives
the (upper, lower) limits for contouring.
annotate_quiver
(self, field_x, field_y, factor, scale=None, scale_units=None, normalize=False)¶(This is a proxy for
QuiverCallback
.)
Adds a ‘quiver’ plot to any plot, using the field_x
and field_y
from
the associated data, skipping every factor
datapoints scale
is the
data units per arrow length unit using scale_units
(see
matplotlib.axes.Axes.quiver for more info)
annotate_cquiver
(self, field_x, field_y, factor)¶(This is a proxy for
CuttingQuiverCallback
.)
Get a quiver plot on top of a cutting plane, using field_x
and
field_y
, skipping every factor
datapoint in the discretization.
annotate_grids
(self, alpha=0.7, min_pix=1, min_pix_ids=20, draw_ids=False, periodic=True, min_level=None, max_level=None, cmap='B-W Linear_r', edgecolors=None, linewidth=1.0)¶(This is a proxy for
GridBoundaryCallback
.)
Adds grid boundaries to a plot, optionally with alpha-blending via the
alpha
keyword. Cuttoff for display is at min_pix
wide. draw_ids
puts the grid id in the corner of the grid. (Not so great in projections...)
annotate_halos
(self, halo_catalog, circle_args=None, width=None, annotate_field=None, text_args=None, factor=1.0)¶(This is a proxy for
HaloCatalogCallback
.)
Accepts a HaloCatalog
and plots a circle at the location of each halo with the radius of the
circle corresponding to the virial radius of the halo. If width
is set
to None (default) all halos are plotted, otherwise it accepts a tuple in
the form (1.0, ‘Mpc’) to only display halos that fall within a slab with
width width
centered on the center of the plot data. The appearance of
the circles can be changed with the circle_kwargs dictionary, which is
supplied to the Matplotlib patch Circle. One can label each of the halos
with the annotate_field, which accepts 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, whereas
'particle_identifier'
shows the halo number). font_kwargs contains the
arguments controlling the text appearance of the annotated field.
Factor is the 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.
annotate_line
(self, p1, p2, coord_system='data', plot_args=None)¶(This is a proxy for
LinePlotCallback
.)
Overplot a line with endpoints at p1 and p2. p1 and p2 should be 2D or 3D coordinates consistent with the coordinate system denoted in the “coord_system” keyword.
annotate_magnetic_field
(self, factor=16, scale=None, scale_units=None, normalize=False)¶(This is a proxy for
MagFieldCallback
.)
Adds a ‘quiver’ plot of magnetic field to the plot, skipping all but every
factor
datapoint. scale
is the data units per arrow length unit using
scale_units
(see matplotlib.axes.Axes.quiver for more info). if
normalize
is True
, the magnetic fields will be scaled by their local
(in-plane) length, allowing morphological features to be more clearly seen
for fields with substantial variation in field strength.
annotate_marker
(self, pos, marker='x', coord_system='data', plot_args=None)¶(This is a proxy for
MarkerAnnotateCallback
.)
Overplot a marker on a position for highlighting specific features.
annotate_particles
(self, width, p_size=1.0, col='k', marker='o', stride=1.0, ptype=None, minimum_mass=None, alpha=1.0)¶(This is a proxy for
ParticleCallback
.)
Adds particle positions, based on a thick slab along axis
with a
width
along the line of sight. p_size
controls the number of pixels
per particle, and col
governs the color. ptype
will restrict plotted
particles to only those that are of a given type. minimum_mass
will
require that the particles be of a given mass minimum mass in solar units.
annotate_sphere
(self, center, radius, circle_args=None, coord_system='data', text=None, text_args=None)¶(This is a proxy for
SphereCallback
.)
Overplot a circle with designated center and radius with optional text.
annotate_streamlines
(self, field_x, field_y, factor=16, density = 1, plot_args=None)¶(This is a proxy for
StreamlineCallback
.)
Add streamlines to any plot, using the field_x
and field_y
from the
associated data, using nx
and ny
starting points that are bounded by
xstart
and ystart
. To begin streamlines from the left edge of the
plot, set start_at_xedge
to True
; for the bottom edge, use
start_at_yedge
. A line with the qmean vector magnitude will cover
1.0/factor
of the image.
annotate_text
(self, pos, text, coord_system='data', text_args=None, inset_box_args=None)¶(This is a proxy for
TextLabelCallback
.)
Overplot text on the plot at a specified position. If you desire an inset box around your text, set one with the inset_box_args dictionary keyword.
annotate_title
(self, title='Plot')¶(This is a proxy for
TitleCallback
.)
Accepts a title
and adds it to the plot.
annotate_velocity
(self, factor=16, scale=None, scale_units=None, normalize=False)¶(This is a proxy for
VelocityCallback
.)
Adds a ‘quiver’ plot of velocity to the plot, skipping all but every
factor
datapoint. scale
is the data units per arrow length unit using
scale_units
(see matplotlib.axes.Axes.quiver for more info). if
normalize
is True
, the velocity fields will be scaled by their local
(in-plane) length, allowing morphological features to be more clearly seen
for fields with substantial variation in field strength (normalize is not
implemented and thus ignored for Cutting Planes).
annotate_timestamp
(x_pos=None, y_pos=None, corner='lower_left', time=True, redshift=False, time_format='t = {time:.0f} {units}', time_unit=None, redshift_format='z = {redshift:.2f}', use_inset_box=False, text_args=None, inset_box_args=None)¶(This is a proxy for
TimestampCallback
.)
Annotates the timestamp and/or redshift of the data output at a specified location in the image (either in a present corner, or by specifying (x,y) image coordinates with the x_pos, y_pos arguments. If no time_units are specified, it will automatically choose appropriate units. It allows for custom formatting of the time and redshift information, as well as the specification of an inset box around the text.
annotate_scale
(corner='lower_right', coeff=None, unit=None, pos=None, max_frac=0.2, min_frac=0.018, text_args=None, inset_box_args=None)¶(This is a proxy for
ScaleCallback
.)
Annotates the scale of the plot at a specified location in the image (either in a preset corner, or by specifying (x,y) image coordinates with the pos argument. Coeff and units (e.g. 1 Mpc) refer to the distance scale you desire to show on the plot. If no coeff and units are specified, an appropriate pair will be determined such that your scale bar is never smaller than min_frac or greater than max_frac of your plottable axis length. For additional text and plot arguments for the text and line, include them as dictionaries to pass to text_args and plot_args.
annotate_triangle_facets
(triangle_vertices, plot_args=None)¶(This is a proxy for
TriangleFacetsCallback
.)
This add a line collection of a SlicePlot’s plane-intersection with the triangles to the plot. This callback is ideal for a dataset representing a geometric model of triangular facets.
annotate_ray
(ray, plot_args=None)¶(This is a proxy for
RayCallback
.)
Adds a line representing the projected path of a ray across the plot. The ray can be either a YTOrthoRayBase, YTRayBase, or a LightRay object. annotate_ray() will properly account for periodic rays across the volume.