API Reference

This page is a detailed API reference guide, where you can see what each function does, and how its api is set up.

swarmpyfac.fac module

The fac module is used for calculated field alligned currents.

This module contains functions used to calculated field alligned currents, and some of the components used in the calculations. This module is dependent on the swarm_util module, which contains functions for utilities and some of the more general components needed.

This module currently only allows caluculation of the single satelite version of the FAC product, and it uses data available from viresclient’s python client or another data source with similar format. The major difference between this and the official ESA product, is that the input from models are their effects at each measurement point, since this is the format delivered by viresclient.

This means that there is a small algorithmic difference compared to what is used in the official ESA products, though experimental comparison have found the difference to be less than 1.2*10^(-4) micro A/m^2, and the average difference to be less than 1% of this.

Running this module as a script will (outside of running some tests) perform the default experiment.

Overview

This section lists each of the functions in this module, and give a short description of them. You can use this as an easier way to look up what you might need.

split_models

Takes measurements and up to several models, and returns the combined model and the residual on the measurements.

radial_current

Calculate the radial currents for time-like series.

single_sat_fac

Calculate the field alligned currents for a single satelinte on time-like series.

single_sat_fac_full

Like single_sat_fac, but it handles input from a dictionary, and can figure out which models are pressent itself.

fetch_data

Fetch the data needed, and can figure out itself whether it needs to download it first.

fac_from_file

Do the entire single satelite field alligned currents calculations, but takes a possible file as input the same way as fetch_data.

Examples

>>> from pyfac.fac import fac_from_file  
>>> import datetime as date
>>> results, input_data = fac_from_file(
...     start=date.datetime(2017, 4, 11),
...     end=date.datetime(2017, 4, 12),
...     force_download=True)  
>>> lat, fac = results[1][:, 0], results[3]  
swarmpyfac.fac.fac_from_file(**options)[source]

Compute FAC for a single satelite from data in an outside source.

This function is inteded as the topmost function of this module. It is therefore suggested to just strait up use this function, when one wants to use this as a module instead of as a script.

The function will read a cdf-file or download it, if such a file is not already cached. It will then calculate FAC (field aligned currents), and other related measures, and return those results.

Parameters

**options – options to be passed to other lower level function(s), mainly fetch_data.

Returns

  • The return is a tuple consiting of the following

  • time (ndarray) – An ndarray of scalars describing the time each point in the calculations corrisond to. Note that it is off compared to the parameter time.

  • position (ndarray) – An ndarray of spherical coordinates describing the postion of each point in the time series where the calculations are suppose to apply at. The coordinates are arranged as theta, phi, r. Note tht this is off compared to the parameter positions.

  • irc (ndarray) – An ndarray of vectors of the radial currents (the radial components, of currents) found when solving one of Maxwells equations on the magnetic field. Note that this is given in the VSC (spacecraft velocity) frame.

  • fac (ndarray) – An ndarray of vectors of the field alligned currents. Note that this is given in the VSC (spacecraft velocity) frame.

  • V (ndarray) – An ndarray of vectors of the velocity in the NEC frame.

  • V_VSC (ndarray) – An ndarray of vectors of the velocity in the VSC frame.

  • pos_ltl (ndarray) – An ndarray of scalars of positions where local time is acounted for. Each position reference the same position as positions output. Note tht this is off compared to the parameter positions.

  • inclination – An ndarray of scalars of the inclination of the aproximation of the magnetic field, where the calculations are supposed to be applied at.

Examples

>>> from pyfac.fac import *  
>>> fac_from_file()  
(array([...]), array([[...]]), array([...]), array([...]), array([[...]]), array([[...]]), array([[...]]), array([...]))
swarmpyfac.fac.fetch_data(force_download=False, temp_file='tempdata.cdf', **options)[source]

Fetch cdf form data from a file or by downloading it.

This function is made to hide the complexity of getting data from a cdf file that you might need to download first. The idea is a dowloaded file will be saved as some temporary file temp_file and then opened. If the file have already been downloaded then it will be assumed that the file at temp_file is the one wanted, unless force_download is set to true to indicate that this cache should be flushed and overwritten by a new download.

Parameters
  • force_download (boolean, optional) – A flag on whether to force downloading of the data and thereby overwrite the temp_file if it already exists. If set to False the program will only download the data if it cannot find the temp_file. Default is False.

  • temp_file (string) – A string of the file path for where the downloaded file should be. If the file already exists then this may take priority depending on the force_download flag.

Returns

A dictionary map of _name_pairings, where the values have been replaced with the data entry of the previous value.

Return type

dict

Examples

>>> from pyfac.fac import *  
>>> fetch_data()  
{...}
swarmpyfac.fac.radial_current(delta_B, velocity, delta_time, MU_0=1.2566370614359173e-06)[source]

Calculate the radial curents (IRC) over a time series.

Parameters
  • delta_B (ndarray) – The change in the B (magnetic) field in a time series in an ndarray of vectors. Standard use requires this to be in the VSC frame.

  • velocity (ndarray) – An ndarray of vectors of the spartial velocities at each point in the timeseries.

  • delta_time (ndarray) – An ndarray of scalars of the time steps between each point in the timeseries of delta_B and velocity.

  • MU_0 (float, optional) – The constant of permubility of vacum, in the given units. Defaults to swarm_util.MU_0, which assumes velocity and delta_time are in base SI units, and delta_B in nT.

Returns

An ndarray of vectors, each of which describes the radial current at the respective point in the timseries.

Return type

ndarray

Examples

>>> from pyfac.fac import *  
>>> vecs = sw.delta(np.arange(15).reshape(5,3))
>>> radial_current(vecs,vecs,np.ones(len(vecs)))
array([-0., -0., -0., -0.])
>>> radial_current(vecs, sw.pack_3d(np.ones((4,)),-1.5 * np.ones((4,)),np.zeros((4,))),np.ones(len(vecs)))
array([-1989.43678865, -1989.43678865, -1989.43678865, -1989.43678865])
swarmpyfac.fac.single_sat_fac(time, positions, B_res, B_model)[source]

Calculate the Field Aligned Currents on a single satelites path.

This function runs through the same logic as needed for the single satelite version of the official FAC product. This function is not designed to deal with data flaged in some way as bad. There are several extra outputs usefull for debuging and comparisons.

Parameters
  • time (ndarray) – An array of scalars of time stamps in seconds.

  • poisitions (ndarray) – An array of vectors (theta, phi,r) of the position of the satelite in spherical coordinates.

  • B_res (ndarray) – An array of vectors of the residual of the magnetic field measurements after the model have been substracted. The vectors are assumed to be in the NEC coordinate system.

  • B_model (ndarray) – An array of vectors of the model of the magnetic field at the point of measurements. The vectors are assumed to be in the NEC coordinate system.

Returns

  • The return is a tuple consiting of the following

  • time (ndarray) – An ndarray of scalars describing the time each point in the calculations corrisond to. Note that it is off compared to the parameter time.

  • position (ndarray) – An ndarray of spherical coordinates describing the postion of each point in the time series where the calculations are suppose to apply at. The coordinates are arranged as theta, phi, r. Note tht this is off compared to the parameter positions.

  • irc (ndarray) – An ndarray of vectors of the radial currents (the radial components, of currents) found when solving one of Maxwells equations on the magnetic field. Note that this is given in the VSC (spacecraft velocity) frame.

  • fac (ndarray) – An ndarray of vectors of the field alligned currents. Note that this is given in the VSC (spacecraft velocity) frame.

  • V (ndarray) – An ndarray of vectors of the velocity in the NEC frame.

  • V_VSC (ndarray) – An ndarray of vectors of the velocity in the VSC frame.

  • pos_ltl (ndarray) – An ndarray of scalars of positions where local time is acounted for. Each position reference the same position as positions output. Note tht this is off compared to the parameter positions.

  • inclination – An ndarray of scalars of the inclination of the aproximation of the magnetic field, where the calculations are supposed to be applied at.

Examples

>>> from pyfac.fac import *  
>>> time = np.arange(0,10,0.1)
>>> positions = sw.pack_3d(np.sin(0.01 * time)*90,(0.01 * time % np.pi)/sw.TO_RADIANS ,np.ones(len(time)))
>>> B_res = np.random.randn(len(time),3)* sw.as_3d(10 *positions[:,1])
>>> B_model = sw.as_3d(positions[:,2])*sw.pack_3d(np.sin(positions[:,0]),np.zeros(len(time)),2*np.cos(positions[:,0]))
>>> single_sat_fac(time,positions,B_res,B_model)  
(array([...]), array([[...]]), array([...]), array([...]), array([[...]]), array([[...]]), array([[...]]), array([...]))
swarmpyfac.fac.single_sat_fac_full(input_data)[source]

Calculate the Field Aligned Currents on a single satelites path.

This is mainly designed as a wrapper for single_sat_fac, which can extract the relevant parameters from a dictionary

This function runs through the same logic as needed for the single satelite version of the official FAC product. This function is not designed to deal with data flaged in some way as bad. There are several extra outputs usefull for debuging and comparisons.

Parameters
  • is only one parameter (There) –

  • dictionary functions effectively as a superset of (This) –

  • relevant parameters. The relevant parameters are given below (the) –

  • time (ndarray) – An array of scalars of time stamps in seconds.

  • Theta (ndarray) – An array of scalars representing the polar angle component of the satelite in spherical coordinates.

  • phi (ndarray) – An array of scalars representing the azimuth angle component of the satelite in spherical coordinates.

  • r (ndarray) – An array of scalars representing the radial angle component of the satelite in spherical coordinates.

  • B_base (ndarray) – An array of vectors of the magnetic field measurements. The vectors are assumed to be in the NEC coordinate system.

  • B_core (ndarray, optional) – An array of vectors of the magnetic field component derived from a model of Earth’s core.

  • B_lithosphere (ndarray, optional) – An array of vectors of the magnetic field component derived from a model of the Earthø’s lithosphere (crust).

  • B_magnetosphere_primary (ndarray, optional) – An array of vectors of the magnetic field component derived from a model of the primary magnetosphere field.

  • B_magnetosphere_secondary (ndarray, optional) – An array of vectors of the magnetic field component derived from a model of the induced field from the magnetosphere.

Returns

  • The return is a tuple consiting of the following

  • time (ndarray) – An ndarray of scalars describing the time each point in the calculations corrisond to. Note that it is off compared to the parameter time.

  • position (ndarray) – An ndarray of spherical coordinates describing the postion of each point in the time series where the calculations are suppose to apply at. The coordinates are arranged as theta, phi, r. Note tht this is off compared to the parameter positions.

  • irc (ndarray) – An ndarray of vectors of the radial currents (the radial components, of currents) found when solving one of Maxwells equations on the magnetic field. Note that this is given in the VSC (spacecraft velocity) frame.

  • fac (ndarray) – An ndarray of vectors of the field alligned currents. Note that this is given in the VSC (spacecraft velocity) frame.

  • V (ndarray) – An ndarray of vectors of the velocity in the NEC frame.

  • V_VSC (ndarray) – An ndarray of vectors of the velocity in the VSC frame.

  • pos_ltl (ndarray) – An ndarray of scalars of positions where local time is acounted for. Each position reference the same position as positions output. Note tht this is off compared to the parameter positions.

  • inclination – An ndarray of scalars of the inclination of the aproximation of the magnetic field, where the calculations are supposed to be applied at.

Examples

>>> from pyfac.fac import *  
>>> time = np.arange(0,10,0.1)
>>> positions = sw.pack_3d(np.sin(0.01 * time)*90,(0.01 * time % np.pi)/sw.TO_RADIANS ,np.ones(len(time)))
>>> B_res = np.random.randn(len(time),3)* sw.as_3d(10 *positions[:,1])
>>> B_core = sw.as_3d(positions[:,2])*sw.pack_3d(np.sin(positions[:,0]),np.zeros(len(time)),2*np.cos(positions[:,0]))
>>> single_sat_fac_full({'time':time,'Theta':positions[:,0],'phi':positions[:,1],'r':positions[:,2],'B_base':B_res+B_core,'B_core':B_core})  
(array([...]), array([[...]]), array([...]), array([...]), array([[...]]), array([[...]]), array([[...]]), array([...]))
swarmpyfac.fac.split_models(baseline, *models)[source]

Seperate the effect of models from the baseline.

Seperate the effect of model and risidual on the baseline, base on precalculated effects of each model. It is assumed that the models combine linearly with each other and a residual (meaning modelA + modelB + residual = baseline, but not modelA * modelB + residual = baseline or other such variations).

Parameters
  • baseline (ndarray) – An array of the base measurements.

  • *models (ndarray) – Arrays of the effects of each model to be removed. all models must have the same length as baseline.

Returns

A tuple, with the first component being an ndarray of the residual and the second component being an ndarray of the combined model. The result statisfies that the sum of the components is equal to the baseline.

Return type

tuple of ndarray

Examples

>>> from pyfac.fac import *  
>>> modelA = np.arange(15).reshape(5,3)**2
>>> modelB = sw.as_3d(np.arange(7,2,-1))
>>> data = (np.random.randn(5,3)*0.1 + modelA + modelB
...         + sw.pack_3d(np.sin(np.arange(5)*0.5),
...                      np.cos(np.arange(5)*0.5),
...                      np.tan(np.arange(5)*0.5)))
>>> split_models(data,modelA,modelB) 
(array([[...,...,...],
       [...,...,...],
       [...,...,...],
       [...,...,...],
       [...,...,...]]), array([[  7.,   8.,  11.],
       [ 15.,  22.,  31.],
       [ 41.,  54.,  69.],
       [ 85., 104., 125.],
       [147., 172., 199.]]))

swarmpyfac.utils module

The swarm_util module provides a host of utilities function.

This module stores a vararity of different functions and constants, which is deemed general usefull when working on swarm data, specifically when such data is loaded from the viresclient python client.

There is a mix of general constants, such as the permuability of vacum, and different types of functions. The functions can generally be categorised into twi types: Functions for working on arrays of scalars or vectors, and functions to work with viresclient. These categories have their own sections, with a short descriptions on each.

Array Utilities

These functions are used to work on numpy arrays of scalars or vectors. They are mostly geared toward dealign with time-like series.

pack_3d

Take 3 equaly long arrays of scalars and make an array of 3d vectors.

as_3d

Like pack_3d, but uses dublicates of the same array of scalars. Usefull for aggregating scalar-vector operations over an array.

map_3d

A map function over arrays of 3d vectors, where the function maped takes a 3d vector to a 3d vector.

NEC_to_VSC

Build a transformation from the NEC frame to the VSC frame. The return is a function, which will transfrom its input from the NEC frame to the VSC frame.

delta

Calculate the changes of a time-like series.

means

Calculate the intermediate value (by linear interpolation) for a time-like series.

spherical_delta

Calculate the changes in spherical coordinates as NEC frame vectors.

curl

Calculate a single component of curl for an array of vectors.

inclination

Calculate the inclination for an awway of vectors.

Viresclient

These are the functions related to fetching and loading data from the viresclient python client.

request_data

Builds and send a request for a CDF file to vires. Most usefull for defaulting all the relevant parameters.

read_cdf

A function to read a cdf file based on a dictionary pairing of attributes you want to read from the cdf file, and what you want to refere to them as afterwards.

swarmpyfac.utils.NEC_to_VSC(velocities)[source]

Construct a function to transform from the NEC frame to the VSC frame.

This function generates another function, which is a transform of vectors from the NEC frame to the VSC frame. This specific transformation may be different for each data point.

Parameters

velocities (ndarray) – The velocities describe thw forward direction in the VSC frame, and they can be different for each data point.

Returns

A function that takes a list of vectors of the same length as velocities and transform them from the NEC frame to the VSC frame.

Return type

function

Examples

>>> from pyfac.utils import *  
>>> v = np.arange(15).reshape(5,3)
>>> trans = NEC_to_VSC(v)
>>> trans(v) 
array([[ 0.7...,  0.7...,  2. ...],
       [ 3.5...,  3.5...,  5. ...],
       [ 6.5...,  6.5...,  8. ...],
       [ 9.5...,  9.5..., 11. ...],
       [12.5..., 12.5..., 14. ...]])
>>> trans((v[::-1]-2.5) * (v + 3)) 
array([[ 49.8...,   9.5...,  57.5...],
       [ 46.0...,  46.4...,  68. ...],
       [ 34.8...,  42.4...,  60.5...],
       [  7.0...,  19.1...,  35. ...],
       [-38.4..., -22.4...,  -8.5...]])
swarmpyfac.utils.as_3d(scalars)[source]

As pack_3d, but with the same scalars in on all components.

Parameters

scalars (array-like) – An array-like of scalars.

Returns

An array of 3 dimensional vectors with the same length as xs.

Return type

ndarray

Examples

>>> from pyfac.utils import *  
>>> vec = np.arange(15).reshape(5,3)
>>> as_3d(range(5)) * vec
array([[ 0.,  0.,  0.],
       [ 3.,  4.,  5.],
       [12., 14., 16.],
       [27., 30., 33.],
       [48., 52., 56.]])
swarmpyfac.utils.curl(delta_x, delta_field, target_index=2)[source]

Compute a single component of curl for an array of vectors.

This function is used to calculate a finite difference apporximation of curl. Its starting point is similar to a finite difference approximation of a partial derivative, where you use both the recoded change in the function and what you differentiate with respect to. In practice this calculate an array of in independent curl operations, as given by a single step.

Parameters
  • delta_x (ndarray) – An array of the finite difference in the considered step.

  • delta_field (ndarray) – An array of the finite difference of the field under the relevant considered step.

  • target_index (int, optional) – The index of the dimension of the curl that should be given as output.

Returns

An array of the curl-like quantity, but only target_index dimension of it.

Return type

ndarray

Examples

>>> from pyfac.utils import *  
>>> vecs = delta(np.arange(15).reshape(5,3))
>>> curl(vecs,vecs)
array([0., 0., 0., 0.])
>>> curl(pack_3d(np.ones((4,)),-1.5 * np.ones((4,)),np.zeros((4,))),vecs)
array([5., 5., 5., 5.])
swarmpyfac.utils.delta(vectors)[source]

Computes the finite difference on an arralike.

Parameters

vectors (ndarray) – An array of vectors or scalars.

Returns

An array of vectors or scalars of shape (n-1,…), when vectors input have shape (n,…).

Return type

ndarray

Examples

>>> from pyfac.utils import *  
>>> delta(np.arange(5))
array([1, 1, 1, 1])
>>> delta(np.array([4, 3, 7, 4, 2, 92]))
array([-1,  4, -3, -2, 90])
>>> delta(np.array([[1, 0, 0], [0, 1, 3], [7, -13, 4]]))
array([[ -1,   1,   3],
       [  7, -14,   1]])
swarmpyfac.utils.inclination(vectors)[source]

Calculates the inclination of vectors.

Calculate the inclination for each individual vector in vectors. This means this is equivalent to the tilt of the vectors toward the third dimention.

Parameters

vectors (ndarray) – An array of 3d vectors.

Returns

An array of scalars, where each is the inclination of their respective vector in vectors.

Return type

ndarray

Examples

>>> from pyfac.utils import *  
>>> inclination(np.arange(15).reshape(5,3))
array([1.10714872, 0.78539816, 0.71469295, 0.68539502, 0.66942998])
>>> np.sin(_)
array([0.89442719, 0.70710678, 0.65538554, 0.63297887, 0.62053909])
swarmpyfac.utils.map_3d(function, vectors)[source]

map_3d maps a function for 1d data to each dimension in 3d data.

Parameters
  • function (ndarray -> ndarray) – A function that takes a scalar array and maps it to another scalar array.

  • vectors (ndarray) – An array of 3 dimensional vectors.

Returns

An array of 3 dimensional vectors.

Return type

ndarray

Examples

>>> from pyfac.utils import *  
>>> vec = np.arange(15).reshape(5,3)
>>> map_3d(lambda xs: [x > 2 and x%4 != 0 for x in xs],vec)
array([[0., 0., 0.],
       [1., 0., 1.],
       [1., 1., 0.],
       [1., 1., 1.],
       [0., 1., 1.]])
swarmpyfac.utils.means(vectors)[source]

Computes the means between this and the next point.

Parameters

vectors (ndarray) – An array of vectors or scalars.

Returns

An array of vectors or scalars of shape (n-1,…), when vectors input have shape (n,…).

Return type

ndarray

Examples

>>> from pyfac.utils import *  
>>> means(np.array(list(range(5))))
array([0.5, 1.5, 2.5, 3.5])
>>> means(np.array([4, 3, 7, 4, 2, 92]))
array([ 3.5,  5. ,  5.5,  3. , 47. ])
>>> means(np.array([[1, 0, 0], [0, 1, 3], [7, -13, 4]]))
array([[ 0.5,  0.5,  1.5],
       [ 3.5, -6. ,  3.5]])
swarmpyfac.utils.pack_3d(xs, ys, zs)[source]

Pack 3 array-likes into an array of 3d vectors.

Parameters
  • xs (array-like) – The first (xs) components of the final vector. Must be equivalent to an array of scalars, generators are not accepted (though ranges are). Must have a reasoable response to len(xs), which must be the same as len(ys) and len(zs).

  • ys (array-like) – The second (ys) components of the final vector. Must be equivalent to an array of scalars, generators are not accepted. Must have a reasoable response to len(ys), which must be the same as len(xs) and len(zs).

  • zs (array-like) – The third (zs) components of the final vector. Must be equivalent to an array of scalars, generators are not accepted. Must have a reasoable response to len(zs), which must be the same as len(ys) and len(xs).

Returns

An array of vectors each with 3 dimensions, and the same length as xs.

Return type

ndarray

Examples

>>> from pyfac.utils import *  
>>> pack_3d(range(5),np.arange(5)*2-3,[np.cos(xs*np.pi) for xs in range(5)])
array([[ 0., -3.,  1.],
       [ 1., -1., -1.],
       [ 2.,  1.,  1.],
       [ 3.,  3., -1.],
       [ 4.,  5.,  1.]])
swarmpyfac.utils.read_cdf(file, **name_pairings)[source]

Read data from a cdf file.

Read a cdf file and construct a dictionary based on name_parings. This may open and read a file as a side-effect.

Parameters
  • file (str or cdfread) – The file to read the data from. Accepts both a string filepath to the file and an opened cdfread file.

  • **name_pairings – Pairs of names used to link a name in the output with its corrisponding name in the cdf file.

Returns

The keys are the keywords used in name_parings, while the values are the variables loaded from the cdf file with the name of the corrisponding argument.

Return type

dict

Examples

We can directly ask for variables:

>>> from pyfac.utils import *  
>>> read_cdf('tempdata.cdf') 
{}
>>> read_cdf('tempdata.cdf', time='Timestamp') 
{'time': array([...])}

We can also make use of predefined options:

>>> options = {'time': 'Timestamp', 'B_base': 'B_NEC'}
>>> read_cdf('tempdata.cdf', **options) 
{'time': array([...]), 'B_base': array([[...]])}

Note

These examples are skiped by doctest due to dependency on the existance of a ‘tempdata.cdf’ file.

swarmpyfac.utils.request_data(start=datetime.datetime(2016, 1, 1, 0, 0), end=datetime.datetime(2016, 1, 2, 0, 0), target_file='tempdata.cdf', filters=[{'parameter': 'Latitude', 'minimum': - 90.0, 'maximum': 90.0}], url='https://vires.services/ows', collection='SW_OPER_MAGA_LR_1B', product_options={'auxiliaries': ['QDLat', 'QDLon']}, sampling_step='PT1S', models=['MCO_SHA_2C', 'MLI_SHA_2C', 'MMA_SHA_2C-Primary', 'MMA_SHA_2C-Secondary'], measurements=['F', 'B_NEC'], to_file=True, **credentials)[source]

Request data from a vires server.

This function sets up a request to a vires server. It is a wrapper on viresclient functionality, where default arguments have been provided.

The data will be saved to as a cdf file.

Parameters
  • start (datatime, optional) – The inclusive starting time of the requested data. Defaults to 1st of Janurary 2016

  • end (datetime, optional) – The exclusive end time of the requested data. Defaults to 2nd of Janurary 2016

  • target_file (str, optional) – Filepath for the output cdf file. Must include name and extension. Defaults to ‘tempdata.cdf’

  • filters (list(dict), optional) – A list of dictionaries, where each dictionary is the options for a filter to be applied to the data.

  • target_url (str, optional) – The full url for the server to request the data at. Defaults to ‘https://staging.viresdisc.vires.services/openows

  • collection (str, optional) – See viresclient set_collection. Defaults to ‘SW_OPER_MAGA_LR_1B’.

  • product_options (dict, optional) – Extra options to viresclient set_products. Defaults to {‘auxiliaries’: [‘QDLat’, ‘QDLon’]}

  • sampling_step (str, optional) – Describes the sampling frequency, since vires may otherwise downsample the data. The default is 1Hz by ‘PT1S’.

  • models (list(str), optional) – The models calculated at the data points. See viresclient set_products models for details. Defaults to [ ‘MCO_SHA_2C’, ‘MLI_SHA_2C’, ‘MMA_SHA_2C-Primary’, ‘MMA_SHA_2C-Secondary’]

  • measurements (list(str), optional) – Measured quantities to be included for each data point. Defaults to [‘F’, ‘B_NEC’]

  • to_file (boolean, optional) – If set to true, the the result will be saved to a file as cdf, otherwise the datastructure will be returned as is. Defaults to True.

  • **credentials (**str, optional) –

    extra parameters to pass to SwarmRequest. Intented for passing credentails (though other options are also possible). It is backwards compatible with calls using parameter names for credentials. The specific credentials are as below: token : str, optional

    Token for authenticating with the target_url. You need either a token or a username/password combination. The default options for the website will be used if all are set to None.

    usernamestr, optional

    username for authenticating with the target_url. You need either a token or a username/password combination. The default options for the website will be used if all are set to None.

    passwordstr, optional

    username for authenticating with the target_url. You need either a token or a username/password combination. The default options for the website will be used if all are set to None.

Examples

>>> from pyfac.utils import *  
>>> request_data() 
>>> request_data(date.datetime(2017, 7, 5),
...              date.datetime(2017, 7, 6))

Note

These examples are skiped by doctest due to security reasons, and unintentional side-effects.

swarmpyfac.utils.spherical_delta(positions)[source]

Computes the change of sherical positions as a NEC vector.

Parameters

positions (ndarray) – An array of positions (3d vectors) described in spherical coordinates.

Returns

ndarray – An array of vectors in the NEC frame.

Return type

ndarray

Examples

>>> from pyfac.utils import *  
>>> vecs = np.arange(15).reshape(5,3)
>>> spherical_delta(vecs)
array([[ 0.18317585,  0.18311308, -3.        ],
       [ 0.34018372,  0.33913504, -3.        ],
       [ 0.49719158,  0.49293804, -3.        ],
       [ 0.65419945,  0.64324482, -3.        ]])