Note
Go to the end to download the full example code.
Multiple TracerMP Instances#
Perform simple tracing using the TracerMP class.
This example demonstrates how to use the TracerMP
class to perform forward tracing of magnetic field lines from a set of default starting points
across both Coronal and Heliospheric magnetic field models.
import matplotlib.pyplot as plt
import numpy as np
from mapflpy.tracer import TracerMP
from mapflpy.utils import plot_traces
from mapflpy.data import fetch_cor_magfiles, fetch_hel_magfiles
Load in the magnetic field files
The fetch_cor_magfiles() function returns a tuple of file paths
corresponding to the radial, theta, and phi components of the magnetic field data.
cor_magnetic_field = fetch_cor_magfiles()
hel_magnetic_field = fetch_hel_magfiles()
The TracerMP class can be instantiated directly with the magnetic
field file paths (along with any additional “mapfl params” as keyword arguments).
For now, the default mapfl configuration is used i.e.
DEFAULT_PARAMS.
Note
The TracerMP class utilizes multiprocessing to instantiate a
new mapflpy_fortran object in a worker process, and communicates with it via
inter-process communication. The Pipe() protocol requires that
the multiprocessing.connection.Connection is properly opened and closed to
avoid resource leaks. Therefore, it is recommended to use the TracerMP
class within a context manager (the with statement) to ensure that resources are
properly managed.
Alternatively, the connect() and
disconnect() methods can be used to manually manage the
connection lifecycle.
tracer_cor = TracerMP(*cor_magnetic_field, context=CONTEXT)
tracer_hel = TracerMP(*hel_magnetic_field, context=CONTEXT)
tracer_cor.connect()
tracer_hel.connect()
traces_cor = tracer_cor.trace_fwd()
Forward trace field lines from the coronal model, then use the end points at r=30 Rs
as launch points for the heliospheric model. When tracing with the heliospheric model,
we need to set the minimum domain radius to 30 Rs using the domain_r_min_ parameter.
outerboundary = np.isclose(traces_cor.end_pos[0,...], 30)
hel_launch_points = traces_cor.end_pos[:, outerboundary]
tracer_hel['domain_r_min_'] = 30.0
traces_hel = tracer_hel.trace_fwd(launch_points=hel_launch_points)
tracer_cor.disconnect()
tracer_hel.disconnect()
Traces can be plotted using the plot_traces() utility function
ax = plt.figure().add_subplot(projection='3d')
plot_traces(traces_cor, ax=ax)
plot_traces(traces_hel, ax=ax, color='orange')
FOV = 100 # Rsun
for dim in 'xyz':
getattr(ax, f'set_{dim}lim3d')((-FOV, FOV))
ax.set_box_aspect([1, 1, 1])
plt.show()

Total running time of the script: (0 minutes 3.168 seconds)