Magnetosphere regions and bow shock: How to find
Note: The analysator methods here are prototypes and should be adapted and modified as needed.
Bow shock
Plasma properties for estimating bow shock position:
- plasma compression:
\(n_p > 2n_{p, sw}\) [Battarbee_et_al_2020] (Vlasiator)
- solar wind core heating:
\(T_{core} > 4T_{sw}\) [Battarbee_et_al_2020] (Vlasiator)
\(T_{core} = 3T_{sw}\) [Suni_et_al_2021] (Vlasiator)
- magnetosonic Mach number:
\(M_{ms} < 1\) [Battarbee_et_al_2020] (Vlasiator)
In analysator:
regions has an option to find the bow shock. Default method uses 1.5*solar wind density as limit.
Usage example:
datafile = "vlsvbulkfile.vlsv"
outfilen = "bowshock.vlsv"
RegionFlags(datafile, outfilen, regions=["bowshock"])
Magnetosheath
properties:
- density:
\(8 cm^{-3}\) [Hudges_Introduction_to_space_physics_Ch_9]
- temperature:
ion: \(150 eV\) [Hudges_Introduction_to_space_physics_Ch_9]
electron: \(25 eV\) [Hudges_Introduction_to_space_physics_Ch_9]
- magnetic field:
- plasma \(\beta\):
In analysator:
regions has an option to find the magnetosheath using bow shock and magnetopause:
Usage example:
datafile = "vlsvbulkfile.vlsv"
outfilen = "magnetosheath.vlsv"
RegionFlags(datafile, outfilen, regions=["magnetosheath"])
Polar cusps
Properties:
- plasma density: high density in comparison to solar wind
Ion density \(\geq\) solar wind ion density [Pitout_et_al_2006] (Cluster spacecraft data)
- ion energy:
mean ion energy \(~2-3 keV\) [Pitout_et_al_2006] /[Stenuit_et_al_2001]
- energy flux:
energy flux
In analysator:
regions has an option to find cusps using convex hull of the magnetosphere.
Note that the default parameters may catch some dayside cells near the magnetopause and some non-cusp areas inside dayside in addition to cusps.
Areas near magnetopause can be excluded for example by only including cells where the SDF is less than e.g. -1 Re.
Usage example:
datafile = "vlsvbulkfile.vlsv"
outfilen = "cusps.vlsv"
RegionFlags(datafile, outfilen, regions=["cusps"])
Tail lobes
- plasma density: low
below \(0.03 cm^{-3}\) [Grison_et_al_2025] (Cluster spacecraft data)
\(0.01 cm^{-3}\) [Koskinen_Space_Storms] p.38
less than \(0.1 cm^{-3}\) [Wolf_Introduction_to_space_physics_Ch_10] p.291
- plasma \(\beta\): low
typically around \(0.05\) [Grison_et_al_2025] (Cluster spacecraft data)
\(3e-3\) [Koskinen_Space_Storms] p.38
- temperature:
ion temperature \(300 eV\) [Koskinen_Space_Storms] p.38
electron temperature \(50 eV\) [Koskinen_Space_Storms] p.38
- magnetic field:
\(20 nT\) [Koskinen_Space_Storms] p.38
open magnetic field lines [Wolf_Introduction_to_space_physics_Ch_10] p.291
strong and stable magnetic field towards the Earth (northern lobe) and away from the Earth (southern lobe) [Coxon_et_al_2016]
Separated from the plasma sheet by the plasma sheet boundary layer (PSBL)
In analysator:
regions has an option to find tail lobes.
Note that the default parameters to find lobes may include some areas from the dayside.
Usage example:
datafile = "vlsvbulkfile.vlsv"
outfilen = "lobes.vlsv"
RegionFlags(datafile, outfilen, regions=["lobes"])
Low-latitude boundary layer (LLBL)
Properties:
- density:
ion number densities between those of magnetosphere and magnetosheath [Hudges_Introduction_to_space_physics_Ch_9] p.267
- temperature
ion temperatures between those of magnetosphere and magnetosheath [Hudges_Introduction_to_space_physics_Ch_9] p.267
unknown field line configuration, probably a mix of open and closed field lines [Hudges_Introduction_to_space_physics_Ch_9] p.262
High-latitude boundary layer (HLBL)
Includes the plasma mantle on the tail side and the entry layer on the dayside
Properties:
open magnetic field lines [Hudges_Introduction_to_space_physics_Ch_9] p.261
Plasma sheet boundary layer (PSBL)
The plasma sheet boundary layer is a very thin boundary layer separating the tail lobes from the tail plasma sheet [Koskinen_Johdatus]
Properties:
- density:
\(0.1 cm^{-3}\) [Koskinen_Space_Storms] p.38
- temperature:
ion temperature \(1000 eV\) [Koskinen_Space_Storms] p.38
electron temperature \(150 eV\) [Koskinen_Space_Storms] p.38
- magnetic field:
\(20 nT\) [Koskinen_Space_Storms] p.38
- plasma \(\beta\) :
\(0.1\) [Koskinen_Space_Storms] p.38
probably closed magnetic field lines [Wolf_Introduction_to_space_physics_Ch_10] p.291
Central plasma sheet
Properties:
- density:
\(0.3 cm^{-3}\) [Koskinen_Space_Storms] p.38
\(0.1-1 cm^{-3}\) [Wolf_Introduction_to_space_physics_Ch_10] p.291
- temperature: hot
ion temperature \(4200 eV\) [Koskinen_Space_Storms] p.38
electron temperature \(600 eV\) [Koskinen_Space_Storms] p.38
- magnetic field:
\(10 nT\) [Koskinen_Space_Storms] p.38, [Hudges_Introduction_to_space_physics_Ch_9]
- plasma \(\beta\): high
\(6\) [Koskinen_Space_Storms] p.38
Mostly closed magnetic field lines [Wolf_Introduction_to_space_physics_Ch_10]
Inner plasma sheet: unusually low plasma beta may exist (e.g., cold tenuous plasma near the neutral sheet after long periods of northward IMF) [Boakes_et_al_2014], (Cluster spacecraft data)
In analysator:
regions has an option to find the central plasma sheet.
Note that the default parameters may catch some dayside cells in addition to plasma sheet.
Usage example:
datafile = "vlsvbulkfile.vlsv"
outfilen = "CPS.vlsv"
RegionFlags(datafile, outfilen, regions=["central_plasma_sheet"])
Script and functions for creating sidecar files with SDF/region/boundary tags of plasma regions.
Usage example where bow shock conditions are given to replace default conditions:
datafile = "bulkfile.vlsv"
outfilen = "regions.vlsv"
RegionFlags(datafile, outfilen, regions=["all"],
region_conditions = {"bowshock": {"density": [2e6, None]}}
- regions.RegionFlags(datafile, outfilen, regions=['all'], ignore_boundaries=True, region_flag_type='01', magnetopause_kwargs={}, region_conditions={})[source]
Creates a sidecar .vlsv file with flagged cells for regions and boundaries in near-Earth plasma environment. Region flags (named flag_region, flags are fractions of filled conditions or 1/0): magnetosheath, magnetosphere, cusps, lobe_N, lobe_S, central_plasma_sheet Boundary signed distance flags (named “SDF_boundary”, flags are signed distances to boundary in m with inside being negative distance): magnetopause, bowshock
- possible regions: “all”, “boundaries” (magnetopause, bow shock), “large_areas” (boundaries + upstream, magnetosheath, magnetosphere), “magnetosphere”, “bowshock”,
“cusps”, “lobes”, “central_plasma_sheet”
Note that different runs may need different tresholds for region parameters and region accuracy should be verified visually Beta* convex hull is most likely the best magnetopause method for regions, for just magnetopause with different options use magnetopause.py
- Parameters:
datafile – .vlsv bulk file name (and path)
outfilen – sidecar .vlsv file name (and path)
- Kword ignore_boundaries:
True: do not take cells in the inner/outer boundaries of the simulation into account when looking for regions (applicable for cusps and CPS for now)
- Kword region_flag_type:
“01” or “fraction”, whether flags are binary (all conditions must be satisfied) or fractions of how many given conditions are met
- Kword magnetopause_method:
default “beta_star”, other options: see magnetopause.py, use alpha=None
- Kword region_conditions:
optional dict where keys are str region names and values are condition dictionaries, for setting own conditions for bow shock, cusps, lobes or central plasma sheet
- regions.bowshock_SDF(f, variable_dict, query_points, own_condition_dict=None)[source]
Finds the bow shock by making a convex hull of either default variables/values based on upstream values or a user-defined variable/value dictionary. Returns signed distances from the convex hull to given query_points.
- Parameters:
f – a VlsvReader
variable_dict – dictionary containing names and variable arrays needed
query_points – xyz-coordinates of all points where the SDF will be calculated ([x1, y1, z1], [x2, y2, z2], …)
- Kword own_condition_dict:
optional, dictionary with string variable names as keys and tresholds as values to pass to treshold_mask()-function
- regions.box_mask(f, coordpoints=None, cells=None, marginal=[150000000.0, 50000000.0, 50000000.0, 50000000.0, 50000000.0, 50000000.0])[source]
Crops simulation box for calculations, output flags outside of cropped box will be 0
- Parameters:
f – a VlsvReader
- Kword coordpoints:
coordinate points of cells to be masked
- Kword cells:
cellIDs to be masked
- Kword marginal:
6-length list of wanted marginal lengths from mesh edges in meters [negx, negy, negz, posx, posy, posz]
- Returns:
Boolean mask in input order of coordinates/cells, -1 if no coorpoints or cells were given and mask could not be done
- regions.make_region_flags(variable_dict, condition_dict, flag_type='01', mask=None)[source]
Makes region flags
example:
make_region_flags(variable_dict = {"density": f.read_variable(name="proton/vg_rho", cellids=-1)}, condition_dict = {"density": [None, 1e7]})
results in flag array in shape of read_variable output where cells where “proton/vg_rho” is less or equal to 1e7 are 1 and others are 0
- Parameters:
variable_dict – dictionary containing names and data arrays of variables, names must match those in condition_dict
condition_dict – dictionary containing names and conditions for variable data in variable_dict
- Kword flag_type:
default “01”, optionally “fraction”. “01” gives flags 1: all variables fill all conditions, 0: everything else; “fraction” flags are rounded fractions of varibles that fulfill conditions (e.g. flag 0.5 means one of two conditions was filled)
- Kword mask:
deafault None, optionally boolean array to use for variable dict arrays to search for a region from subset of cells (cells False in mask will be automatically flagged as 0)
- regions.treshold_mask(data_array, value)[source]
Chooses and flags cells where variable values match those given. If value is given as float, relative tolerance can be given
- Parameters:
data_array – array to mask, e.g. output of f.read_variable(name=”proton/vg_rho”, cellids=-1)
variable – str, variable name
value – value/values to use for masking; a float or int for exact match, a (value, relative tolerance) tuple, or [min value, max value] list pair where either can be None for less than or eq./more than or eq. value
- Returns:
0/1 mask in same order as data_array, 1: variable value in array inside treshold values, 0: outside
- regions.vtkDelaunay3d_SDF(query_points, coordinates, alpha=None)[source]
Gives a signed distance to a convex hull or alpha shape surface created from given coordinates. Note: if using alpha, SDF might not work, especially if the point cloud used for Delaunay does not fully cover e.g. lobes
- Parameters:
all_points – points ([x, y, z] coordinates in m) for which a signed distance to surface will be calculated
coordinates – coordinates (array of [x, y, z]:s in m) that are used to make a surface.
- Kword alpha:
alpha to be given to vtkDelaunay3D (e.g. R_E), removes surface edges that have length more than alpha so that the resulting surface is not convex. None -> convex hull
- Returns:
vtkDataSetSurfaceFilter() object, array of signed distances (negative sign: inside the surface, positive sign: outside the surface)
- regions.write_flags(writer, flags, flag_name, mask=None)[source]
Writes flags into .vlsv-file with writer
- Parameters:
writer – a vlsvwriter
flags – an array of flags in order of cellids
flag_name – string, name of the flag variable
- Kword mask:
if flags are made from cropped cells, give mask used for cropping
References
Battarbee, M., Ganse, U., Pfau-Kempf, Y., Turc, L., Brito, T., Grandin, M., Koskela, T., and Palmroth, M.: Non-locality of Earth’s quasi-parallel bow shock: injection of thermal protons in a hybrid-Vlasov simulation, Ann. Geophys., 38, 625-643, https://doi.org/10.5194/angeo-38-625-2020, 2020
Suni, J., Palmroth, M., Turc, L., Battarbee, M., Johlander, A., Tarvus, V., et al. (2021). Connection between foreshock structures and the generation of magnetosheath jets: Vlasiator results. Geophysical Research Letters, 48, e2021GL095655. https://doi. org/10.1029/2021GL095655
Grison, B., Darrouzet, F., Maggiolo, R. et al. Localization of the Cluster satellites in the geospace environment. Sci Data 12, 327 (2025). https://doi.org/10.1038/s41597-025-04639-z
Koskinen, H. E. J. (2011). Johdatus plasmafysiikkaan ja sen avaruussovellutuksiin. Limes ry.
Koskinen, H. E. J. (2011). Physics of Space Storms: From the Solar Surface to the Earth. Springer-Verlag. https://doi.org/10.1007/978-3-642-00319-6
Pitout, F., Escoubet, C. P., Klecker, B., and Rème, H.: Cluster survey of the mid-altitude cusp: 1. size, location, and dynamics, Ann. Geophys., 24, 3011–3026, https://doi.org/10.5194/angeo-24-3011-2006, 2006.
Coxon,J.C.,C.M.Jackman, M. P. Freeman, C. Forsyth, and I. J. Rae (2016), Identifying the magnetotail lobes with Cluster magnetometer data, J. Geophys. Res. Space Physics, 121, 1436–1446, doi:10.1002/2015JA022020.
Hudges, W. J. (1995) The magnetopause, magnetotail and magnetic reconnection. In Kivelson, M. G., & Russell, C. T. (Eds.), Introduction to space physics (pp.227-287). Cambridge University Press.
Wolf, R. A. (1995) Magnetospheric configuration. In Kivelson, M. G., & Russell, C. T. (Eds.), Introduction to space physics (pp.288-329). Cambridge University Press.
Boakes, P. D., Nakamura, R., Volwerk, M., and Milan, S. E. (2014). ECLAT Cluster Spacecraft Magnetotail Plasma Region Identifications (2001–2009). Dataset Papers in Science, 2014(1):684305. eprint: https://onlinelibrary.wiley.com/doi/pdf/10.1155/2014/684305
Stenuit, H., Sauvaud, J.-A., Delcourt, D. C., Mukai, T., Kokubun, S., Fujimoto, M., Buzulukova, N. Y., Kovrazhkin, R. A., Lin, R. P., and Lepping, R. P. (2001). A study of ion injections at the dawn and dusk polar edges of the auroral oval. Journal of Geophysical Research: Space Physics, 106(A12):29619–29631. eprint: https://agupubs.onlinelibrary.wiley.com/doi/pdf/10.1029/2001JA900060.