SunpyPlotter#
- class sunkit_pyvista.SunpyPlotter(*args, **kwargs)[source]#
-
A plotter for 3D data.
This class inherits
pyvista.Plotterso we can provide coordinate-aware plotting. For now, all coordinates are converted to a specific frame (HeliocentricInertialby default), and distance units are such that \(R_{sun} = 1\).- Parameters:
coordinate_frame (
astropy.coordinates.BaseCoordinateFrame) – Coordinate frame of the plot. The x, y, z axes of the pyvista plotter will be the x, y, z axes in this coordinate system.obstime (
astropy.time.Time) – The obstime to use for the default coordinate frame ifcoordinate_frame=is not specified. Must not be specified ifcoordinate_frameis given.kwargs (dict) – All other keyword arguments are passed through to
pyvista.Plotter.
Attributes Summary
Coordinate frame of the plot.
Methods Summary
coordinates_to_polydata(coords)Convert a set of coordinates in a
SkyCoordto apyvista.PolyDatamesh.load(filepath)Load saved meshes into this plotter.
plot_coordinates(coords, *[, radius])Plot a sphere if a single coordinate is passed and plots a line if multiple coordinates are passed.
plot_current_sheet(pfss_out, **kwargs)Plot current sheet, where \(B_r=0\).
plot_field_lines(field_lines, *[, color_func])Plot magnetic field lines.
plot_limb(m, *[, radius])Draws the solar limb as seen by the map's observer.
plot_map(m, *[, clip_interval, ...])Plot a sunpy map.
plot_quadrangle(bottom_left, *[, top_right, ...])Plot a quadrangle.
plot_solar_axis(*[, length, arrow_kwargs])Plot the solar rotation axis as an arrow.
save(filepath, *[, overwrite])Save all the meshes.
set_camera_coordinate(coord)Set the camera position.
set_camera_focus(coord)Set the camera focus.
set_view_angle(angle)Set the camera view angle.
Attributes Documentation
- coordinate_frame#
Coordinate frame of the plot.
Methods Documentation
- coordinates_to_polydata(coords)[source]#
Convert a set of coordinates in a
SkyCoordto apyvista.PolyDatamesh.- Parameters:
coords (
astropy.coordinates.SkyCoord) – Coordinates to convert.- Return type:
- load(filepath)[source]#
Load saved meshes into this plotter.
- Parameters:
filepath (
strorpathlib.Path) – Name of the file to load, should have vtm or vtmb as an extension.
- plot_coordinates(coords, *, radius=0.05, **kwargs)[source]#
Plot a sphere if a single coordinate is passed and plots a line if multiple coordinates are passed.
- Parameters:
coords (
astropy.coordinates.SkyCoord) – Coordinate(s) to plot as a center of sphere or line.radius (
int, optional) – Radius of the sphere times the radius of the sun to be plotted when a single coordinate is passed. Defaults to0.05times the radius of the sun.**kwargs – Keyword arguments are passed to
pyvista.Plotter.add_mesh.
Notes
This plots a
pyvista.Sphereobject if a single coordinate is passed and plots apyvista.Splineobject if multiple coordinates are passed.radiusis only considered when a sphere is plotted.
- plot_current_sheet(pfss_out, **kwargs)[source]#
Plot current sheet, where \(B_r=0\).
- Parameters:
pfss_out (
sunkit_magex.pfss.output.Output) – Magnetic field calculated bysunkit_magex.**kwargs – Keyword arguments are passed to
pyvista.Plotter.add_mesh.
- plot_field_lines(field_lines, *, color_func=None, **kwargs)[source]#
Plot magnetic field lines.
- Parameters:
field_lines (
sunkit_magex.pfss.fieldline.FieldLines) – Field lines to be plotted.color_func (function) –
Function to get the color for each field line. If not given, defaults to showing closed field lines in black, and open field lines in blue (positive polarity) or red (negative polarity). The function must have the signature:
def color_func(field_line: sunkit_magex.pfss.fieldline.FieldLine) -> color: Where ``color`` is any color that `pyvista` recognises (e.g. a RGBA tuple, a RGB tuple, a color string)
**kwargs – Keyword arguments are handed to
pyvista.Plotter.add_mesh.
- plot_limb(m, *, radius=0.02, **kwargs)[source]#
Draws the solar limb as seen by the map’s observer.
- Parameters:
m (
sunpy.map.Map) – Map’s limb to be plotted.radius (
float) – Radius of thepyvista.Splineused to create the limb. Defaults to0.02times the radius of the sun.**kwargs (Keyword arguments are handed to
pyvista.Plotter.add_mesh.)
- plot_map(
- m,
- *,
- clip_interval: Unit('%') = None,
- assume_spherical_screen=True,
- **kwargs,
Plot a sunpy map.
- Parameters:
m (
sunpy.map.Map) – Map to be plotted.clip_interval (two-element
Quantity, optional) – If provided, the data will be clipped to the percentile interval bounded by the two numbers.assume_spherical_screen (bool, optional) – If
True(default) then off-limb pixels are plotted usingsunpy.coordinates.screens.SphericalScreen. IfFalse, off-limb pixels are not plotted.**kwargs – Keyword arguments are handed to
pyvista.Plotter.add_mesh.
- plot_quadrangle(
- bottom_left,
- *,
- top_right=None,
- width: Unit('deg') = None,
- height: Unit('deg') = None,
- radius=0.01,
- **kwargs,
Plot a quadrangle.
This draws a quadrangle that has corners at
(bottom_left, top_right), ifwidthandheightare specified, they are respectively added to the longitude and latitude of thebottom_leftcoordinate to calculate atop_rightcoordinate.- Parameters:
bottom_left (
SkyCoord) – The bottom-left coordinate of the quadrangle. It can have shape(2,)to simultaneously definetop_right.top_right (
SkyCoord) – The top-right coordinate of the quadrangle.width (
astropy.units.Quantity, optional) – The width of the quadrangle. Required iftop_rightis omitted.height (
astropy.units.Quantity) – The height of the quadrangle. Required iftop_rightis omitted.radius (
float) – Radius of thepyvista.Splineused to create the quadrangle. Defaults to0.01times the radius of the sun.**kwargs (Keyword arguments are handed to
pyvista.Plotter.add_mesh.)
- plot_solar_axis(*, length=2.5, arrow_kwargs=None, **kwargs)[source]#
Plot the solar rotation axis as an arrow.
- Parameters:
length (float) – Length of the arrow in multiples of solar radii.
arrow_kwargs (dict) – Keyword arguments to be handed to
pyvista.Arrow.start,direction, andscalecannot be manually specified, as they are automatically set.**kwargs – Keyword arguments are handed to
pyvista.Plotter.add_mesh.
- save(filepath, *, overwrite=False)[source]#
Save all the meshes.
This saves the rendered plot as a vtm extended file as well as a directory of the individual meshes with the specified name.
- Parameters:
filepath (
strorpathlib.Path) – Name of the file to save as, should have vtm or vtmb as an extension.
Examples
>>> from sunkit_pyvista import SunpyPlotter >>> plotter = SunpyPlotter() >>> plotter.plot_solar_axis() >>> plotter.save("./filename.vtm")
- set_camera_coordinate(coord)[source]#
Set the camera position.
- Parameters:
coords (
astropy.coordinates.SkyCoord) – Camera coordinate.
- set_camera_focus(coord)[source]#
Set the camera focus.
- Parameters:
coords (
astropy.coordinates.SkyCoord) – Camera coordinate.
- set_view_angle(angle: Unit('deg'))[source]#
Set the camera view angle.
- Parameters:
angle (
astropy.units.Quantity) – The viewing angle.