NEI

class plasmapy_nei.nei.NEI(inputs, abundances: ~typing.Dict | str | None = None, T_e: ~typing.Callable | ~astropy.units.quantity.Quantity | None = None, n: ~typing.Callable | ~astropy.units.quantity.Quantity | None = None, time_input: ~astropy.units.quantity.Quantity | None = None, time_start: ~astropy.units.quantity.Quantity | None = None, time_max: ~astropy.units.quantity.Quantity | None = None, max_steps: int | ~numpy.integer = 10000, tol: int | float = 1e-15, dt: ~astropy.units.quantity.Quantity | None = None, dt_max: ~astropy.units.quantity.Quantity = <Quantity inf s>, dt_min: ~astropy.units.quantity.Quantity = <Quantity 0. s>, adapt_dt: bool | None = None, safety_factor: int | float = 1, verbose: bool = False)

Bases: object

Perform and analyze a non-equilibrium ionization simulation.

Parameters:
  • inputs

  • T_e (astropy.units.Quantity or callable) – The electron temperature, which may be a constant, an array of temperatures corresponding to the times in time_input, or a function that yields the temperature as a function of time.

  • n (astropy.units.Quantity or callable) – The number density multiplicative factor. The number density of each element will be n times the abundance given in abundances. For example, if abundance['H'] = 1, then this will correspond to the number density of hydrogen (including neutral hydrogen and protons). This factor may be a constant, an array of number densities over time, or a function that yields a number density as a function of time.

  • time_input (astropy.units.Quantity, optional) – An array containing the times associated with n and T_e in units of time.

  • time_start (astropy.units.Quantity, optional) – The start time for the simulation. If density and/or temperature are given by arrays, then this argument must be greater than time_input[0]. If this argument is not supplied, then time_start defaults to time_input[0] (if given) and zero seconds otherwise.

  • time_max (astropy.units.Quantity) – The maximum time for the simulation. If density and/or temperature are given by arrays, then this argument must be less than time_input[-1].

  • max_steps (int) – The maximum number of time steps to be taken during a simulation.

  • dt (astropy.units.Quantity) – The time step. If adapt_dt is False, then dt is the time step for the whole simulation.

  • dt_max (astropy.units.Quantity) – The maximum time step to be used with an adaptive time step.

  • dt_min (astropy.units.Quantity) – The minimum time step to be used with an adaptive time step.

  • adapt_dt (bool) – If True, change the time step based on the characteristic ionization and recombination time scales and change in temperature. Not yet implemented.

  • safety_factor (float or int) – A multiplicative factor to multiply by the time step when adapt_dt is True. Lower values improve accuracy, whereas higher values reduce computational time. Not yet implemented.

  • tol (float) – The absolute tolerance to be used in comparing ionic fractions.

  • verbose (bool, optional) – A flag stating whether or not to print out information for every time step. Setting verbose to True is useful for testing. Defaults to False.

  • abundances (dict) –

Examples

>>> import numpy as np
>>> import astropy.units as u
>>> inputs = {'H': [0.9, 0.1], 'He': [0.9, 0.099, 0.001]}
>>> abund = {'H': 1, 'He': 0.085}
>>> n = u.Quantity([1e9, 1e8], u.cm**-3)
>>> T_e = np.array([10000, 40000]) * u.K
>>> time = np.array([0, 300]) * u.s
>>> dt = 0.25 * u.s

The initial conditions can be accessed using the initial attribute.

>>> sim = NEI(inputs=inputs, abundances=abund, n=n, T_e=T_e, time_input=time, adapt_dt=False, dt=dt)

After having inputted all of the necessary information, we can run the simulation.

>>> results = sim.simulate()

The initial results are stored in the initial attribute.

>>> sim.initial.ionic_fractions['H']
array([0.9, 0.1])

The final results can be access with the final attribute.

>>> sim.final.ionic_fractions['H']
array([0.16665179, 0.83334821])
>>> sim.final.ionic_fractions['He']
array([0.88685261, 0.11218358, 0.00096381])
>>> sim.final.T_e
<Quantity 40000. K>

Both initial and final are instances of the IonizationStateCollection class.

Notes

The ionization and recombination rates are from Chianti version 8.7. These rates include radiative and dielectronic recombination. Photoionization is not included.

Attributes Summary

T_e_input

The temperature input.

abundances

Return the abundances.

adapt_dt

Return True if the time step is set to be adaptive, False if the time step is set to not be adapted, and None if this attribute was not set.

dt_input

Return the inputted time step.

dt_max

dt_min

The minimum time step.

eigen_data_dict

Return a dict containing EigenData instances for each element.

elements

A list of the elements.

final

Return the ionization states of the plasma at the end of the simulation.

initial

Return the ionization states of the plasma at the beginning of the simulation.

max_steps

The maximum number of steps that a simulation will be allowed to take.

n_input

The number density factor input.

results

Return the SimulationResults class instance that corresponds to the simulation results.

safety_factor

The multiplicative factor that the time step is to be multiplied by when using an adaptive time step.

time_input

time_max

The maximum time allowed for the simulation.

time_start

The start time of the simulation.

tol

The tolerance for comparisons between different ionization states.

verbose

Return True if verbose output during a simulation is requested, and False otherwise.

Methods Summary

electron_temperature(time)

equil_ionic_fractions([T_e, time])

Return the equilibrium ionic fractions for a temperature or at a given time.

hydrogen_number_density(time)

in_time_interval(time[, buffer])

Return True if the time is between time_start - buffer and time_max + buffer , and False otherwise.

index_to_time(index)

Returns the time value or array given the index/indices

save([filename])

Save the NEI instance to an HDF5 file.

set_timestep([dt])

Set the time step for the next non-equilibrium ionization time advance.

simulate()

Perform a non-equilibrium ionization simulation.

time_advance()

Advance the simulation by one time step.

time_to_index(time)

Returns the closest index value or array for the given time(s)

Attributes Documentation

T_e_input

The temperature input.

abundances

Return the abundances.

adapt_dt

Return True if the time step is set to be adaptive, False if the time step is set to not be adapted, and None if this attribute was not set.

dt_input

Return the inputted time step.

dt_max
dt_min

The minimum time step.

eigen_data_dict

Return a dict containing EigenData instances for each element.

elements

A list of the elements.

final

Return the ionization states of the plasma at the end of the simulation.

initial

Return the ionization states of the plasma at the beginning of the simulation.

max_steps

The maximum number of steps that a simulation will be allowed to take.

n_input

The number density factor input.

results

Return the SimulationResults class instance that corresponds to the simulation results.

safety_factor

The multiplicative factor that the time step is to be multiplied by when using an adaptive time step.

time_input
time_max

The maximum time allowed for the simulation.

time_start

The start time of the simulation.

tol

The tolerance for comparisons between different ionization states.

verbose

Return True if verbose output during a simulation is requested, and False otherwise.

Methods Documentation

electron_temperature(time: astropy.units.Quantity) astropy.units.Quantity
equil_ionic_fractions(T_e: astropy.units.Quantity | None = None, time: astropy.units.Quantity | None = None) Dict[str, ndarray]

Return the equilibrium ionic fractions for a temperature or at a given time.

Parameters:
Returns:

equil_ionfracs – The equilibrium ionic fractions for the elements contained within this class

Return type:

dict

Notes

Only one of T_e and time may be included as an argument. If neither T_e or time is provided and the temperature for the simulation is given by a constant, the this method will assume that T_e is the temperature of the simulation.

hydrogen_number_density(time: astropy.units.Quantity) astropy.units.Quantity
in_time_interval(time: Unit("s"), buffer: Unit("s") = <Quantity 1.e-09 s>)

Return True if the time is between time_start - buffer and time_max + buffer , and False otherwise.

Raises:
index_to_time(index)

Returns the time value or array given the index/indices

Parameters:

index (array-like) – A value or array of values representing the index of the time array created by the simulation

Returns:

get_time – The time value associated with index input(s)

Return type:

astropy.units.Quantity

save(filename: str = 'nei.h5')

Save the NEI instance to an HDF5 file. Not implemented.

set_timestep(dt: astropy.units.Quantity | None = None)

Set the time step for the next non-equilibrium ionization time advance.

Parameters:

dt (astropy.units.Quantity, optional) – The time step to be used for the next time advance.

Notes

If dt is not None, then the time step will be set to dt.

If dt is not set and the adapt_dt attribute of an NEI instance is True, then this method will calculate the time step corresponding to how long it will be until the temperature rises or drops into the next temperature bin. If this time step is between dtmin and dtmax, then

If dt is not set and the adapt_dt attribute is False, then this method will set the time step as what was inputted to the NEI class upon instantiation in the dt argument or through the NEI class’s dt_input attribute.

Raises:

NEIError – If the time step cannot be set, for example if the dt argument is invalid or the time step cannot be adapted.

simulate() SimulationResults

Perform a non-equilibrium ionization simulation.

Returns:

results – The results from the simulation (which are also stored in the results attribute of the NEI instance this method was called from.

Return type:

Simulation

time_advance()

Advance the simulation by one time step.

time_to_index(time)

Returns the closest index value or array for the given time(s)

Parameters:

time (array-like,) – A value or array of values representing the values of the time array created by the simulation

Returns:

index – The index value associated with the time input(s)

Return type:

int or array-like,