vlsvfile

VlsvFile is responsibe for implementing read and write functionality for the vlsv format, as well as defining data reducers and reduction pipelines. Some of this will be likely moved around to be more modular in some indefinite point in the future.

The vlsvfile module has all the classes related to reading and writing vlsv files

# Example:
import analysator as pt
pt.vlsvfile.
#press [tab] -> get the functions
class analysator.vlsvfile.VlasiatorReader(file_name, fsGridDecomposition=None)

Class for reading VLSV files with support for Vlasiator velocity space and structures

class analysator.vlsvfile.VlsvParticles(file_name)

Class for reading VLSV files of particle pusher output

list()

Print out a description of the content of the file. Useful for interactive usage

read_particles()

Read particle pusher data from the open vlsv file. This routine returns only active particles.

read_particles_all()

Read particle pusher data from the open vlsv file. This routine returns all particles, including inactive ones.

class analysator.vlsvfile.VlsvReader(file_name, fsGridDecomposition=None)

Class for reading VLSV files

build_cell_vertices(cid, prune_unique=False)

Builds, caches and returns the vertices that lie on the surfaces of CellIDs cid.

Parameters:
  • cid – numpy array of CellIDs

  • prune_unique – bool [False], if you suspect you might be calling the function many times with the same CellID in the list, it might be beneficial to enable this and not repeat the operation for duplicate entries.

Returns:

Dictionary of cell c (int) : set of vertex indices (3-tuple) that touch the cell c.

cellid_has_vdf(cid, pop='proton') bool

Returns whether the cid in question has a vdf or not :param coords: the cellid to test for :returns: bool

check_parameter(name)

Checks if a given parameter is in the vlsv reader

Parameters:

name – Name of the parameter

Returns:

True if the parameter is in the vlsv file, false if not

Note

This should be used for checking if a parameter exists, e.g. for different Vlasiator versions and time output

# Example usage:
vlsvReader = pt.vlsvfile.VlsvReader("test.vlsv")
if vlsvReader.check_parameter( "time" ):
   time = vlsvReader.read_parameter("time")
elif vlsvReader.check_parameter( "t" ):
   time = vlsvReader.read_parameter("t")
else:
   time = None
check_population(popname)

Checks if a given population is in the vlsv file

Parameters:

name – Name of the population

Returns:

True if the population is in the vlsv file, false if not

# Example usage:
vlsvReader = pt.vlsvfile.VlsvReader("test.vlsv")
if vlsvReader.check_population( "avgs" ):
   plot_population('avgs')
else:
   if vlsvReader.check_population( "proton" ):
      # File is newer with proton population
      plot_population('proton')
check_variable(name)

Checks if a given variable is in the vlsv reader

Parameters:

name – Name of the variable

Returns:

True if the variable is in the vlsv file, false if not

Note

This should be used for checking if a variable exists in case a function behaves differently for ex. if B vector is in the vlsv and if not

# Example usage:
vlsvReader = pt.vlsvfile.VlsvReader("test.vlsv")
if vlsvReader.check_variable( "B" ):
   # Variable can be plotted
   plot_B()
else:
   # Variable not in the vlsv file
   plot_B_vol()
construct_velocity_cell_coordinates(blocks)

Returns velocity cell coordinates in given blocks

Parameters:

blocks – list of block ids

Returns:

a numpy array containing the velocity cell ids e.g. np.array([4,2,56,44,522, ..])

construct_velocity_cell_nodes(blocks, pop='proton')

Returns velocity cell nodes in given blocks

Parameters:

blocks – list of block ids

Returns:

a numpy array containing velocity cell nodes and the keys for velocity cells

Note

This is used for constructing velocity space inside the mayavi module

See also

grid

construct_velocity_cells(blocks)

Returns velocity cells in given blocks

Parameters:

blocks – list of block ids

Returns:

a numpy array containing the velocity cell ids e.g. np.array([4,2,56,44,522, ..])

downsample_fsgrid_subarray(cellid, array)

Returns a mean value of fsgrid values underlying the SpatialGrid cellid.

fsgrid_array_to_vg(array)

Downsample, via averaging, an fsgrid array to the Vlasov grid of this reader.

Parameters:

array – array with first three dimensions corresponding to the dimensions of the fsgrid associated with this reader.

Returns:

Vlasov grid data (in file order) of array averaged to Vlasov Grid.

get_all_variables()

Returns all variables in the vlsv reader and the data reducer :returns: List of variable is in the vlsv file .. code-block:: python

# Example usage: vlsvReader = pt.vlsvfile.VlsvReader(“test.vlsv”) vars = vlsvReader.get_variables()

get_amr_level(cellid)

Returns the AMR level of a given cell defined by its cellid

Parameters:

cellid – The cell’s cellid

Returns:

The cell’s refinement level in the AMR

get_bbox_fsgrid_slicemap(low, up)

Returns a slice tuple of fsgrid indices that are contained in the (low, up) bounding box.

get_bbox_fsgrid_subarray(low, up, array)

Returns a subarray of the fsgrid array, corresponding to the (low, up) bounding box.

get_cell_bbox(cellid)

Returns the bounding box of a given cell defined by its cellid

Parameters:

cellid – The cell’s cellid

Returns:

The cell’s bbox [xmin,ymin,zmin],[xmax,ymax,zmax]

get_cell_coordinates(cellids)

Returns a given cell’s coordinates as a numpy array

Parameters:

cellids – The array of cell IDs

Returns:

a numpy array with the coordinates

See also

get_cellid()

Note

The cell ids go from 1 .. max not from 0

get_cell_corner_vertices(cids)

Builds, caches and returns the vertices that lie on the corners of CellIDs cid. :parameter cid: numpy array of CellIDs

Returns:

Dictionary of cell c (int) : 8-tuple of vertex indices (3-tuples).

get_cell_dx(cellid)

Returns the dx of a given cell defined by its cellid

Parameters:

cellid – The cell’s cellid

Returns:

The cell’s size [dx, dy, dz]

get_cell_fsgrid(cellid)

Returns a slice tuple of fsgrid indices that are contained in the SpatialGrid cell.

get_cell_fsgrid_slicemap(cellid)

Returns a slice tuple of fsgrid indices that are contained in the SpatialGrid cell.

get_cell_fsgrid_subarray(cellid, array)

Returns a subarray of the fsgrid array, corresponding to the fsgrid covered by the SpatialGrid cellid.

get_cell_indices(cellids, reflevels=None)

Returns a given cell’s indices as a numpy array

Parameters:
  • cellid – The cell’s ID, numpy array

  • reflevel – The cell’s refinement level in the AMR

Returns:

a numpy array with the coordinates

See also

get_cellid()

Note

The cell ids go from 1 .. max not from 0

get_cell_neighbor(cellidss, offsetss, periodic, prune_uniques=False)

Returns a given cells neighbor at offset (in indices)

Parameters:
  • cellids – The cell’s ID

  • offsets – The offset to the neighbor in indices

  • periodic – For each dimension, is the system periodic

Returns:

the cellid of the neighbor

Note

Returns 0 if the offset is out of bounds!

get_cellid(coords)

Returns the cell ids at given coordinates

Parameters:

coords – The cells’ coordinates

Returns:

the cell ids

Note

Returns 0 if the cellid is out of bounds!

get_cellid_locations()

Returns a dictionary with cell id as the key and the index of the cell id as the value. The index is used to locate the cell id’s values in the arrays that this reader returns

get_cellid_with_vdf(coords, pop='proton')

Returns the cell ids nearest to test points, that contain VDFs

Parameters:

coords – Test coordinates [meters] of N_in points in ND-dimensional space array with shape [N_in, ND] or [ND]

Example: cellid = vlsvReader.get_cellid_with_vdf(np.array([1e8, 0, 0])) :returns: the cell ids

get_config()

Gets config information from VLSV file

:returns a nested dictionary of dictionaries,

where keys (str) are config file group headings (appearing in ‘[]’) and values are dictionaries which contain (lists of) strings

If the same heading/parameter pair appears >once in the config file, the different values are appended to the list .

EXAMPLE: if the config contains these lines:

[proton_precipitation] nChannels = 9

then the following returns [‘9’]: vlsvReader.get_config()[‘proton_precipitation’][‘nChannels’]

get_config_string()

Gets config information from VLSV file. TAG is hardcoded to CONFIG

:returns configuration file string if config is found otherwise returns None

get_dual(pts, cellids=None)

Find the duals that contain the coordinate points pts. This will call the iterative find_ksi function to see if the resulting interpolation weights for the coordinate are in the range [0,1]; if not, it will iterate through neighbouring duals until a dual is found. :parameter pts: numpy array of coordinates (N,3)

Returns:

duals (numpy array of N 3-tuples), ksis (numpy array of interpolation weights (N, 8))

get_duals(cids)
Get the union of dual cells that cover each of CellIDs in cids.

Assumes all required duals are defined! TODO handling of missing duals, do not call separately.

Returns:

Dict of vertex-indices v (3-tuple) : 8-tuple of cellids (corners of dual cells indexed by v)

get_fsgrid_cell_size()

Read fsgrid cell size

Returns:

Maximum and minimum coordinates of the mesh, [dx, dy, dz]

get_fsgrid_coordinates(ri)

Returns real-space center coordinates of the fsgrid 3-index.

get_fsgrid_indices(coords)

Convert spatial coordinates coords to an index array [xi, yi, zi] for fsgrid

:returns 3-tuple of integers [xi, yi, zi] corresponding to fsgrid cell containing coords (low-inclusive) Example: ii = f.get_fsgrid_mesh_extent(coords) fsvar_at_coords = fsvar_array.item(ii)

get_fsgrid_mesh_extent()

Read fsgrid mesh extent

Returns:

Maximum and minimum coordinates of the mesh, [xmin, ymin, zmin, xmax, ymax, zmax]

get_fsgrid_mesh_size()

Read fsgrid mesh size

Returns:

Size of mesh in number of cells, array with three elements

get_fsgrid_slice_indices(lower, upper, eps=0.001)

Get indices for a subarray of an fsgrid variable, in the cuboid from lower to upper. This is meant for mapping a set of fsgrid cells to a given SpatialGrid cell. Shifts the corners (lower, upper) by dx_fsgrid*eps inward, if direct low-inclusive behaviour is required, set kword eps = 0.

:returns two 3-tuples of integers. Example: ii = f.get_fsgrid_mesh_extent(coords) fsvar_at_coords = fsvar_array.item(ii)

get_ionosphere_element_coords()

Read coordinates of ionosphere elements (triangle barycenters)

Returns:

[x, y, z] array of ionosphere element barycenter coordinates (in meters).

get_ionosphere_element_corners()

Read ionosphere mesh element corners

Returns:

[c1,c2,c3] array of ionosphere mesh node indices (starting from 0)

get_ionosphere_latlon_coords()

Read ionosphere node coordinates (in magnetic longitude / latitude)

Returns:

[lat,lon] array of ionosphere node coordinates

get_ionosphere_mesh_area()

Read areas of ionosphere elements (triangular mesh)

Returns:

1D array, areas of the triangular elements [m^2]

get_ionosphere_mesh_size()

Read size of the ionosphere mesh, if there is one.

Returns:

Size of the mesh in number of nodes and elements, array with two elements

get_ionosphere_node_coords()

Read ionosphere node coordinates (in cartesian GSM coordinate system).

Returns:

[x,y,z] array of ionosphere node coordinates (in meters)

get_max_refinement_level()

Returns the maximum refinement level of the AMR

get_precipitation_centre_energy(pop='proton')

Read precipitation energy bins

Returns:

Array of centre energies

get_spatial_block_size()

Read spatial mesh block size

Returns:

Size of block in number of cells, array with three elements

get_spatial_mesh_extent()

Read spatial mesh extent

Returns:

Maximum and minimum coordinates of the mesh, [xmin, ymin, zmin, xmax, ymax, zmax]

get_spatial_mesh_size()

Read spatial mesh size

Returns:

Size of mesh in number of blocks, array with three elements

get_unique_cellids(coords)
Returns a list of cellids containing all the coordinates in coords,

with no duplicate cellids. Relative order of elements is conserved.

Parameters:

coords – A list of coordinates

Returns:

a list of unique cell ids

get_velocity_block_coordinates(blocks, pop='proton')

Returns the block coordinates of the given blocks in a numpy array

Parameters:

blocks – list of block ids

Returns:

a numpy array containing the block coordinates e.g. np.array([np.array([2,1,3]), np.array([5,6,6]), ..])

get_velocity_block_indices(blocks, pop='proton')

Returns the block indices of the given blocks in a numpy array

Parameters:

blocks – list of block ids

Returns:

a numpy array containing the block indices e.g. np.array([np.array([2,1,3]), np.array([5,6,6]), ..])

get_velocity_block_size(pop='proton')

Read velocity mesh block size

Returns:

Size of block in number of cells, array with three elements

get_velocity_blocks(blockCoordinates, pop='proton')

Returns the block ids of the given block coordinates in a numpy array form

Parameters:

blockcoordinates – list of block coordinates e.g. np.array([np.array([2,1,3]), np.array([5,6,6]), ..])

Returns:

a numpy array containing the block ids e.g. np.array([4,2,56,44,2, ..])

get_velocity_cell_coordinates(vcellids, pop='proton')

Returns a given velocity cell’s coordinates as a numpy array

Arguments: :param vcellids: The velocity cell’s ID :returns: a numpy array with the coordinates

get_velocity_cell_ids(vcellcoord, pop='proton')

Returns velocity cell ids of given coordinate

Arguments: :param vcellcoords: One 3d coordinate :returns: Velocity cell id

get_velocity_mesh_dv(pop='proton')

Read velocity mesh cell size

Returns:

Velocity mesh cell size, array with three elements [dvx, dvy, dvz]

get_velocity_mesh_extent(pop='proton')

Read velocity mesh extent

Returns:

Maximum and minimum coordinates of the mesh, [vxmin, vymin, vzmin, vxmax, vymax, vzmax]

get_velocity_mesh_size(pop='proton')

Read velocity mesh size

Returns:

Size of mesh in number of blocks, array with three elements

get_vertex_coordinates_from_indices(indices)

Convert vertex indices to physical coordinates. :param indices: numpy array of vertex indices, either (3,) or (N,3)

Returns:

numpy array of coordinates, with matching shape to indices

get_vertex_indices(coordinates)

Get dual grid vertex indices for all coordinates.

Works by truncation to integer indices at fsgrid resolution, for cell low-corners. :param coordinates: np.array of coordinates, shaped either (3,) or (N,3)

list(parameter=True, variable=True, mesh=False, datareducer=False, operator=False, other=False)

Print out a description of the content of the file. Useful for interactive usage. Default is to list parameters and variables, query selection can be adjusted with keywords:

Default and supported keywords:

parameter=True variable=True mesh=False datareducer=False operator=False other=False

optimize_clear_fileindex_for_cellid()

Clears a private variable containing cell ids and their locations

Note

This should only be used for optimization purposes.

optimize_clear_fileindex_for_cellid_blocks()

Clears a private variable containing number of blocks and offsets for particular cell ids

Note

This should only be used for optimization purposes.

optimize_close_file()

Closes the vlsv file Files are opened and closed automatically upon reading and in the case of reading multiple times it will help to keep the file open with this command

Note

This should only be used for optimization purposes.

optimize_open_file()

Opens the vlsv file for reading Files are opened and closed automatically upon reading and in the case of reading multiple times it will help to keep the file open with this command

Note

This should only be used for optimization purposes.

print_config()

Prints config information from VLSV file.

:returns True if config is found otherwise returns False

print_version()

Prints version information from VLSV file. TAG is hardcoded to VERSION

:returns True if version is found otherwise returns False

read(name='', tag='', mesh='', operator='pass', cellids=-1)

Read data from the open vlsv file.

Parameters:
  • name – Name of the data array

  • tag – Tag of the data array.

  • mesh – Mesh for the data array

  • operator – Datareduction operator. “pass” does no operation on data.

  • cellids – If -1 then all data is read. If nonzero then only the vector for the specified cell id or cellids is read

Returns:

numpy array with the data

read_attribute(name='', mesh='', attribute='', tag='')

Read data from the open vlsv file.

Parameters:
  • name – Name of the data array

  • tag – Tag of the data array.

  • mesh – Mesh for the data array

  • operator – Datareduction operator. “pass” does no operation on data.

  • cellids – If -1 then all data is read. If nonzero then only the vector for the specified cell id or cellids is read

Returns:

numpy array with the data

read_blocks(cellid, pop='proton')

Read raw block data from the open file and return the data along with block ids

Parameters:

cellid – Cell ID of the cell whose velocity blocks are read

Returns:

A numpy array with block ids and data eg [array([2, 5, 6, 234, 21]), array([1.0e-8, 2.1e-8, 2.1e-8, 0, 4.0e-8])]

read_fsgrid_variable(name, operator='pass')

Reads fsgrid variables from the open vlsv file. Arguments: :param name: Name of the variable :param operator: Datareduction operator. “pass” does no operation on data :returns: ordered numpy array with the data

… seealso:: read_variable()

read_fsgrid_variable_cellid(name, cellids=-1, operator='pass')

Reads fsgrid variables from the open vlsv file. Arguments: :param name: Name of the variable :param cellids: SpatialGrid cellids for which to fetch data. Default: return full fsgrid data :param operator: Datareduction operator. “pass” does no operation on data :returns: ordered list of numpy arrays with the data

… seealso:: read_fsgrid_variable()

read_interpolated_fsgrid_variable(name, coordinates, operator='pass', periodic=[True, True, True], method='linear')

Read a linearly interpolated FSgrid variable value from the open vlsv file. Feel free to vectorize! Note that this does not account for varying centerings of fsgrid data. Arguments: :param name: Name of the (FSgrid) variable :param coords: Coordinates from which to read data :param periodic: Periodicity of the system. Default is periodic in all dimension :param operator: Datareduction operator. “pass” does no operation on data :returns: numpy array with the data

read_interpolated_ionosphere_variable(name, coordinates, operator='pass', method='linear')

Read a linearly interpolated ionosphere variable value from the open vlsv file. Arguments: :param name: Name of the (ionosphere) variable :param coords: Coordinates (x,y,z) from which to read data :param operator: Datareduction operator. “pass” does no operation on data :param method: Interpolation method. Not implemented; barycentric interp would fall under linear. :returns: numpy array with the data

read_interpolated_variable(name, coords, operator='pass', periodic=[True, True, True], method='linear')

Read a linearly interpolated variable value from the open vlsv file. Arguments: :param name: Name of the variable :param coords: Coordinates from which to read data :param periodic: Periodicity of the system. Default is periodic in all dimension :param operator: Datareduction operator. “pass” does no operation on data :param method: Interpolation method, default “linear”, options: [“nearest”, “linear”]

Returns:

numpy array with the data

read_interpolated_variable_irregular(name, coords, operator='pass', periodic=[True, True, True], method='linear', methodargs={'delaunay': {'qhull_options': 'QJ'}, 'rbf': {'neighbors': 64}})

Read a linearly interpolated variable value from the open vlsv file. Arguments:

Parameters:
  • name – Name of the variable

  • coords – Coordinates from which to read data

  • periodic – Periodicity of the system. Default is periodic in all dimension

  • operator – Datareduction operator. “pass” does no operation on data

  • method – Method for interpolation, default “linear” (“nearest”, “rbf, “delaunay”)

  • methodargs

    Dict of dicts to pass kwargs to interpolators. Default values for “rbf”, “delaunay”;

    see scipy.interpolate.RBFInterpolator for rbf and scipy.interpolate.LinearNDInterpolator for delaunay

Returns:

numpy array with the data

read_ionosphere_node_variable_at_elements(varname)

Interpolate an ionospheric node variavle to the element barycenters

Parameters:

varname – string, specifying variable (or data reducer) defined at ionosphere nodes

Returns:

specified variable interpolated to the elements’ barycenters

note: linear barycentric interpolation is just the sum of 3 corner values divided by 3!

TODO: check behavior for var.ndim>1

read_ionosphere_variable(name, operator='pass')

Reads fsgrid variables from the open vlsv file. Arguments: :param name: Name of the variable :param operator: Datareduction operator. “pass” does no operation on data :returns: numpy array with the data in node order

… seealso:: read_variable()

read_metadata(name='', tag='', mesh='')

Read variable metadata from the open vlsv file.

Parameters:
  • name – Name of the data array

  • tag – Tag of the data array.

  • mesh – Mesh for the data array

Returns:

four strings: the unit of the variable as a regular string the unit of the variable as a LaTeX-formatted string the description of the variable as a LaTeX-formatted string the conversion factor to SI units as a string

read_parameter(name)

Read a parameter from the vlsv file

Parameters:

name – Name of the parameter

Returns:

The parameter value

read_variable(name, cellids=-1, operator='pass')

Read variables from the open vlsv file. Arguments: :param name: Name of the variable :param cellids: a value of -1 reads all data :param operator: Datareduction operator. “pass” does no operation on data :returns: numpy array with the data

read_variable_from_cache(name, cellids, operator)

Read variable from cache instead of the vlsv file. :param name: Name of the variable :param cellids: a value of -1 reads all data :param operator: Datareduction operator. “pass” does no operation on data :returns: numpy array with the data, same format as read_variable

See also

read_variable()

read_variable_info(name, cellids=-1, operator='pass')

Read variables from the open vlsv file and input the data into VariableInfo

Parameters:
  • name – Name of the variable

  • cellids – a value of -1 reads all data

  • operator – Datareduction operator. “pass” does no operation on data

Returns:

numpy array with the data

See also

read_variable()

read_variable_to_cache(name, operator='pass')

Read variable from vlsv file to cache, for the whole grid and after applying operator. :param name: Name of the variable (or datareducer) :param operator: Datareduction operator. “pass” does no operation on data.

read_velocity_cells(cellid, pop='proton')

Read velocity cells from a spatial cell

Parameters:

cellid – Cell ID of the cell whose velocity cells the function will read

Returns:

Map of velocity cell ids (unique for every velocity cell) and corresponding value

#Example:

example_cellid = 1111

velocity_cell_map = vlsvReader.read_velocity_cells(example_cellid) velocity_cell_ids = velocity_cell_map.keys() velocity_cell_values = velocity_cell_map.values()

random_index = 4 # Just some index random_velocity_cell_id = velocity_cell_ids[random_index]

print (“Velocity cell value at velocity cell id “ + str(random_velocity_cell_id) + “: “ + str(velocity_cell_map[random_velocity_cell_id]))

# Getting the corresponding coordinates might be more useful than having the velocity cell id so: velocity_cell_coordinates = vlsvReader.get_velocity_cell_coordinates(velocity_cell_ids) # Get velocity cell coordinates corresponding to each velocity cell id

random_velocity_cell_coordinates = velocity_cell_ids[random_index] print(“Velocity cell value at velocity cell id “ + str(random_velocity_cell_id) + “and coordinates “ + str(random_velocity_cell_coordinates) + “: “ + str(velocity_cell_map[random_velocity_cell_id]))

See also

read_blocks()

upsample_fsgrid_subarray(cellid, var, array)

Set the elements of the fsgrid array to the value of corresponding SpatialGrid cellid. Mutator for array.

class analysator.vlsvfile.VlsvVtkReader
RequestData(request, inInfo, outInfo)

Overwritten by subclass to execute the algorithm.

RequestInformation(request, inInfo, outInfo)

Overwritten by subclass to provide meta-data to downstream pipeline.

class analysator.vlsvfile.VlsvWriter(vlsvReader, file_name, copy_meshes=None, clone=False)

Class for reading VLSV files

clone_file(vlsvReader, dst)

Simply copies overs the file in vlsvReader to a fresh new file :param vlsvReader: Some open vlsv file :param dst: Name of output file

copy_variables(vlsvReader, varlist=None)

Copies variables from vlsv reader to the file. varlist = None: list of variables to copy; if no varlist is provided, copy all variables (default)

copy_variables_list(vlsvReader, vars)

Copies variables in the list vars from vlsv reader to the file

write(data, name, tag, mesh, extra_attribs={})

Writes an array into the vlsv file

Parameters:
  • name – Name of the data array

  • tag – Tag of the data array.

  • mesh – Mesh for the data array

  • extra_attribs – Dictionary with whatever xml attributes that should be defined in the array that aren’t name, tag, or mesh

Returns:

True if the data was written successfully

write_variable(data, name, mesh, units, latex, latexunits, unitConversion, extra_attribs={})

Writes an array into the vlsv file as a variable; requires input of metadata required by VlsvReader :param data: The variable data (array) :param name: Name of the data array :param mesh: Mesh for the data array :param latex: LaTeX string representation of the variable name :param units: plaintext string representation of the unit :param latexunits: LaTeX string representation of the unit :param unitConversion: string representation of the unit conversion to get to SI :param extra_attribs: Dictionary with whatever xml attributes that should be defined in the array that aren’t name, tag, or mesh.

Returns:

True if the data was written successfully

write_variable_info(varinfo, mesh, unitConversion, extra_attribs={})

Writes an array into the vlsv file as a variable; requires input of metadata required by VlsvReader

Parameters:
  • varinfo – VariableInfo object containing -data: The variable data (array) -name: Name of the data array -latex: LaTeX string representation of the variable name -units: plaintext string representation of the unit -latexunits: LaTeX string representation of the unit

  • mesh – Mesh for the data array

  • unitConversion – string representation of the unit conversion to get to SI

  • extra_attribs – Dictionary with whatever xml attributes that should be defined in the array that aren’t name, tag, or mesh, or contained in varinfo. Can be used to overwrite varinfo values besids name.

Returns:

True if the data was written successfully

write_velocity_space(vlsvReader, cellid, blocks_and_values)

Writes given velocity space into vlsv file

Parameters:
  • vlsvReader – Some open vlsv reader file with velocity space in the given cell id

  • cellid – Given cellid

  • blocks_and_values – Blocks and values in list format e.g. [[block1,block2,..], [block1_values, block2_values,..]] where block1_values are velocity block values (list length 64)