Magnetopause: how to find
Note: The analysator methods here are prototypes and should be adapted and modified as needed.
Magnetopause
Plasma beta, beta*
Modified plasma beta
[Xu_et_al_2016], [Brenner_et_al_2021]
The beta* gets values below 1 inside the magnetosphere near magnetopause and can be used to create a magnetopause surface.
Caveats: magnetotail current sheet has beta* \(>\) 1
in analysator:
see magnetopause.magnetopause() “method” keyword options “beta_star”, “beta_star_with_connectivity”
datareducer: beta_star, vg_beta_star
Xu, S., M. W. Liemohn, C. Dong, D. L. Mitchell, S. W. Bougher, and Y. Ma (2016), Pressure and ion composition boundaries at Mars, J. Geophys. Res. Space Physics, 121, 6417–6429, doi:10.1002/2016JA022644.
Brenner A, Pulkkinen TI, Al Shidi Q and Toth G (2021) Stormtime Energetics: Energy Transport Across the Magnetopause in a Global MHD Simulation. Front. Astron. Space Sci. 8:756732. doi: 10.3389/fspas.2021.756732
Field line connectivity
Bulkfiles from newer Vlasiator runs include connection (whether magnetic field lines are closed-closed, closed-open, open/closed or open-open) and coordinates for field line start and end points that can be used in finding the mangetopause location.
Solar wind flow
Method used in e.g. [Palmroth_et_al_2003]
Streamlines of velocity field that are traced from outside the bow shock curve around the magnetopause.
Caveats: sometimes some streamlines can curve into the magnetotail or dayside magnetoshpere
In analysator:
see
calculations.find_magnetopause_sw_streamline_2d()
calculations.find_magnetopause_sw_streamline_3d()
Streamlines are traced from outside the bow shock towards Earth. A subsolar point for the magnetopause is chosen to be where streamlines get closest to Earth in x-axis [y/z~0].
From subsolar point towards Earth, space is divided radially by spherical coordiante (theta from x-axis) angles theta and phi, and magnetopause is located by looking at streamline point distances from origo and marked to be at middle of the sector
From the Earth towards negative x the space is divided into yz-planes and streamlines are interpolated to certain x-points. Each yz-plane is then divided into 2d sectors and magnetopause is marked to be in the middle of the sector with the radius of the n:th closest streamline to the x-axis.
For subsolar point, radial dayside, and -x planes the closest streamline point can be changed to be n:th closest by setting keyword ignore, in which case ignore closest streamline points are not taken into account.
2d: Note: no radial dayside, uses yt-package instead of analysator’s fieldtracer
3d: After the magnetopause points are chosen, they are made into a surface by setting connnection lines between vertices (magnetopause points) to make surface triangles. This surface is then made into a vtkPolyData and returned as vtkDataSetSurfaceFilter that can be written into e.g. .vtp file.
Palmroth, M., T. I. Pulkkinen, P. Janhunen, and C.-C. Wu (2003), Stormtime energy transfer in global MHD simulation, J. Geophys. Res., 108, 1048, doi:10.1029/2002JA009446, A1.
Shue et al. (1997)
The Shue et al. (1997) [Shue_et_al_1997] magnetopause model:
where \(D_p\) is the dynamic pressure of the solar wind and \(B_z\) the magnetic field z-component magnitude. \(r_0\) is the magnetopause standoff distance and \(\alpha\) is the level of tail flaring.
The magnetopause radius as a function of angle \(\theta\) is
In analysator:
see shue
Shue, J.-H., Chao, J. K., Fu, H. C., Russell, C. T., Song, P., Khurana, K. K., and Singer, H. J. (1997). A new functional form to study the solar wind control of the magnetopause size and shape. Journal of Geophysical Research: Space Physics, 102(A5):9497–9511. eprint: https://agupubs.onlinelibrary.wiley.com/doi/pdf/10.1029/97JA00196.
In analysator:
magnetopause in scripts for 3d runs
Constructs the magnetopause surface with vtk’s vtkDelaunay3d triangulation with optional alpha to make the surface non-convex.
Uses regions.py functions.
Important: SDF of non-convex surface might not always work as expected inside the magnetosphere if e.g. the lobes are hollow.
options (magnetopause() method keyword) and some notes:
- solar wind flow (“streamlines”)
uses magnetopause_sw_streamline_3d.py from pyCalculations
if streamlines make a turn so that the velocity points sunwards (to pos. x), the streamline is ignored from that point onwards
- streamlines that enter an area where beta* is below 0.4 are also stopped
- this can affect the far magnetotail as streamlines may turn inwards into lobes at some point and may result in magnetopause not forming
due to not enough streamline points where magnetopause is supposed to be
very dependent on how solar wind streamlines behave
streamline seeds and other options can greatly affect the resulting magnetopause
- beta* (“beta_star”)
beta* treshold might need tweaking as sometimes there are small low beta* areas in the magnetosheath that get taken in distorting the magnetopause shape at nose
convex hull (Delaunay_alpha=None) usually makes a nice rough magnetopause but goes over any inward dips (like polar cusps)
alpha shape (Delaunay_alpha= e.g. 1*R_E) does a better job at cusps and delicate shapes like vortices but might fail at SDF inside the magnetosphere
Delaynay3d has an easier time if the treshold is something like [0.4, 0.5] and not [0.0, 0.5], but this can affect the alpha shape if alpha is used
- beta* with magnetic field line connectivity (“beta_star_with_connectivity”)
includes closed-closed magnetic field line areas if available, otherwise like “beta_star”
can help with nose shape as beta* can be set lower to exclude magnetosheath low beta* areas while still getting the full dayside from field lines
- Shue et al. 1997 (“shue”)
uses shue.py from scripts
a rough theoretical magnetopause using Shue et al. 1997 method based on B_z, solar wind density, and solar wind velocity
- user-defined parameter tresholds (“dict”)
creates a magnetopause (or other area) using the Delaunay3d triangulation of some area where user-defined tresholds given as dictionary
dictionary key is data name in datafile and value is treshold used, if dictionary has multiple conditions, they all need to be fulfilled
dictionary example: {“vg_rho”: [None, 1e5]} makes a magnetopause using cells where density is less than 1e5
Functions for finding the magnetopause from .vlsv data
Example use with own arguments passed to streamline function:
datafile = "bulkfile.vlsv"
vtpoutfilen = "magnetopause.vtp"
vlsvoutfilen = "magnetopause.vlsv"
SW_args = {"seeds_n":300, "seeds_x0":150e6, "seeds_range":[-5*6371000, 5*6371000],
"dl":5e5, "iterations":int(1500), "end_x":-50*6371000, "x_point_n":100, "sector_n":36*2}
surface, SDF = magnetopause(datafile,
method="streamlines",
return_SDF=True,
return_surface=True,
method_args=SW_args)
write_vtk_surface_to_file(surface, vtpoutfilen)
write_SDF_to_file(SDF, datafile, vlsvoutfilen)
Example use for 0.0-0.4 beta* magnetopause with tetrahedrons with sides longer than 1Re excluded (aplha = 1Re):
datafile = "bulkfile.vlsv"
vtpoutfilen = "magnetopause.vtp"
surface, __ = magnetopause(datafile,
method="beta_star_with_connectivity",
beta_star_range=[0.0, 0.4],
Delaunay_alpha=1*6371000,
return_SDF=False,
return_surface=True)
write_vtk_surface_to_file(surface, vtpoutfilen)
- magnetopause.magnetopause(datafilen, method='beta_star_with_connectivity', own_tresholds=None, return_surface=True, return_SDF=True, SDF_points=None, Delaunay_alpha=None, beta_star_range=[0.4, 0.5], method_args={})[source]
- Finds the magnetopause using the specified method. Surface is constructed using vtk’s Delaunay3d triangulation which results in a convex hull if no Delaunay_alpha is given.
Returns vtk.vtkDataSetSurfaceFilter object and/or signed distances (negative -> inside magnetopause) (=SDF) to all cells Note that using alpha for Delaunay might make SDF different from expected inside the magnetosphere, especially if surface is constructed with points not everywhere in the magnetosphere (e.g. beta* 0.4-0.5) or if simulation grid size is larger than alpha
- Parameters:
datafilen – a .vlsv bulk file name (and path)
method – str, default “beta_star_with_connectivity”, other options “beta_star”, “streamlines”, “shue”, “dict”
own_tresholds – if using method “dict”, a dictionary with conditions for the magnetopause/magnetosphere must be given where key is data name in datafile and value is treshold used (see treshold_mask())
return_surface – True/False, return vtkDataSetSurfaceFilter object
return_SDF – True/False, return array of distances in m to SDF_points in point input order, negative distance inside the surface
SDF_points – optionally give array of own points to calculate signed distances to. If not given, distances will be to cell centres in the order of f.read_variable(“CellID”) output
Delaunay_alpha – alpha (float) to give to vtkDelaunay3d, None -> convex hull, alpha=__: surface egdes longer than __ will be excluded
beta_star_range – [min, max] treshold rage to use with methods “beta_star” and “beta_star_with_connectivity”
method_args – dict of keyword arguments to be passed down to external functions (for streamlines and shue)
- Returns:
vtkDataSetSurfaceFilter object of convex hull or alpha shape if return_surface=True, signed distance field of convex hull or alpha shape of magnetopause if return_SDF=True