yt.visualization.streamlines module

class yt.visualization.streamlines.Streamlines(ds, positions, xfield='velocity_x', yfield='velocity_x', zfield='velocity_x', volume=None, dx=None, length=None, direction=1, get_magnitude=False)[source]

Bases: ParallelAnalysisInterface

A collection of streamlines that flow through the volume

The Streamlines object contains a collection of streamlines defined as paths that are parallel to a specified vector field.

Parameters:
  • ds (Dataset) – This is the dataset to streamline

  • positions (array_like) – An array of initial starting positions of the streamlines.

  • xfield (str or tuple of str, optional) – The x component of the vector field to be streamlined. Default:’velocity_x’

  • yfield (str or tuple of str, optional) – The y component of the vector field to be streamlined. Default:’velocity_y’

  • zfield (str or tuple of str, optional) – The z component of the vector field to be streamlined. Default:’velocity_z’

  • volume (yt.extensions.volume_rendering.AMRKDTree, optional) – The volume to be streamlined. Can be specified for finer-grained control, but otherwise will be automatically generated. At this point it must use the AMRKDTree. Default: None

  • dx (float, optional) – Optionally specify the step size during the integration. Default: minimum dx

  • length (float, optional) – Optionally specify the length of integration. Default: np.max(self.ds.domain_right_edge-self.ds.domain_left_edge)

  • direction (real, optional) – Specifies the direction of integration. The magnitude of this value has no effect, only the sign.

  • get_magnitude (bool, optional) – Specifies whether the Streamlines.magnitudes array should be filled with the magnitude of the vector field at each point in the streamline. This seems to be a ~10% hit to performance. Default: False

Examples

>>> import matplotlib.pyplot as plt
... import numpy as np
... from mpl_toolkits.mplot3d import Axes3D
... import yt
... from yt.visualization.api import Streamlines
>>> # Load the dataset and set some parameters
>>> ds = load("IsolatedGalaxy/galaxy0030/galaxy0030")
>>> c = np.array([0.5] * 3)
>>> N = 100
>>> scale = 1.0
>>> pos_dx = np.random.random((N, 3)) * scale - scale / 2.0
>>> pos = c + pos_dx
>>> # Define and construct streamlines
>>> streamlines = Streamlines(
...     ds, pos, "velocity_x", "velocity_y", "velocity_z", length=1.0
... )
>>> streamlines.integrate_through_volume()
>>> # Make a 3D plot of the streamlines and save it to disk
>>> fig = plt.figure()
>>> ax = Axes3D(fig)
>>> for stream in streamlines.streamlines:
...     stream = stream[np.all(stream != 0.0, axis=1)]
...     ax.plot3D(stream[:, 0], stream[:, 1], stream[:, 2], alpha=0.1)
>>> plt.savefig("streamlines.png")
clean_streamlines()[source]
comm = None
get_dependencies(fields)
integrate_through_volume()[source]
partition_index_2d(axis)
partition_index_3d(ds, padding=0.0, rank_ratio=1)
partition_index_3d_bisection_list()

Returns an array that is used to drive _partition_index_3d_bisection, below.

partition_region_3d(left_edge, right_edge, padding=0.0, rank_ratio=1)

Given a region, it subdivides it into smaller regions for parallel analysis.

path(streamline_id)[source]

Returns an YTSelectionContainer1D object defined by a streamline.

Parameters:

streamline_id (int) – This defines which streamline from the Streamlines object that will define the YTSelectionContainer1D object.

Return type:

An YTStreamline YTSelectionContainer1D object

Examples

>>> from yt.visualization.api import Streamlines
>>> streamlines = Streamlines(ds, [0.5] * 3)
>>> streamlines.integrate_through_volume()
>>> stream = streamlines.path(0)
>>> fig, ax = plt.subplots()
>>> ax.set_yscale("log")
>>> ax.plot(stream["t"], stream[("gas", "density")], "-x")
yt.visualization.streamlines.sanitize_length(length, ds)[source]