Lightweaver API Documentation

Submodules

lightweaver.atmosphere module

class lightweaver.atmosphere.Atmosphere(structure: Layout, temperature: ndarray, vturb: ndarray, ne: ndarray, nHTot: ndarray, B: ndarray | None = None, gammaB: ndarray | None = None, chiB: ndarray | None = None)[source]

Bases: object

Storage for all atmospheric data. These arrays will be shared directly with the backend, so a modification here also modifies the data seen by the backend. Be careful to modify these arrays in place, as their data is shared by direct memory reference. Use the class methods to construct atmospheres of different dimensionality.

structure

A layout structure holding the atmospheric stratification, and velocity description.

Type:

Layout

temperature

The atmospheric temperature structure.

Type:

np.ndarray

vturb

The atmospheric microturbulent velocity structure.

Type:

np.ndarray

ne

The electron density structure in the atmosphere.

Type:

np.ndarray

nHTot

The total hydrogen number density distribution throughout the atmosphere.

Type:

np.ndarray

B

The magnitude of the stratified magnetic field throughout the atmosphere (Tesla).

Type:

np.ndarray, optional

gammaB

Co-altitude (latitude) of magnetic field vector (radians) throughout the atmosphere from the local vertical.

Type:

np.ndarray, optional

chiB

Azimuth of magnetic field vector (radians) in the x-y plane, measured from the x-axis.

Type:

np.ndarray, optional

property Ndim: int

int The dimensionality (1, 2, or 3) of the atmospheric model.

Type:

Ndim

property Noutgoing: int

int The number of cells at the top of the atmosphere (that each produce a spectrum).

Type:

Noutgoing

property Nrays

int Number of rays in angular discretisation used.

Type:

Nrays

property Nspace

int Total number of points in the atmospheric spatial discretistaion.

Type:

Nspace

property Nx: int

int The number of points in the x-direction discretisation.

Type:

Nx

property Ny: int

int The number of points in the y-direction discretisation.

Type:

Ny

property Nz: int

int The number of points in the y-direction discretisation.

Type:

Nz

property cmass: ndarray

np.ndarray Column mass [kg m-2].

Type:

cmass

configure_bcs(upOnly: bool = False)[source]

Configure the required angular information for all boundary conditions on the model.

Parameters:

upOnly (bool, optional) – Whether to only configure boundary conditions for up-going rays. (default: False)

dimensioned_unit_view()[source]

Returns a view over the contents of Layout reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, and the correct astropy.units.

dimensioned_view()[source]

Returns a view over the contents of Layout reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, as these are all stored under a flat scheme.

classmethod make_1d(scale: ScaleType, depthScale: ndarray, temperature: ndarray, vlos: ndarray, vturb: ndarray, ne: ndarray | None = None, hydrogenPops: ndarray | None = None, nHTot: ndarray | None = None, B: ndarray | None = None, gammaB: ndarray | None = None, chiB: ndarray | None = None, lowerBc: BoundaryCondition | None = None, upperBc: BoundaryCondition | None = None, convertScales: bool = True, abundance: AtomicAbundance | None = None, logG: float = 2.44, Pgas: ndarray | None = None, Pe: ndarray | None = None, Ptop: float | None = None, PeTop: float | None = None, verbose: bool = False)[source]

Constructor for 1D Atmosphere objects. Optionally will use an equation of state (EOS) to estimate missing parameters.

If sufficient information is provided (i.e. all required parameters and ne and (hydrogenPops or nHTot)) then the EOS is not invoked to estimate any thermodynamic properties. If both of nHTot and hydrogenPops are omitted, then the electron pressure will be used with the Wittmann equation of state to estimate the mass density, and the hydrogen number density will be inferred from this and the abundances. If, instead, ne is omitted, then the mass density will be used with the Wittmann EOS to estimate the electron pressure. If both of these are omitted then the EOS will be used to estimate both. If:

  • Pgas is provided, then this gas pressure will define the atmospheric stratification and will be used with the EOS.

  • Pe is provided, then this electron pressure will define the atmospheric stratification and will be used with the EOS.

  • Ptop is provided, then this gas pressure at the top of the atmosphere will be used with the log gravitational acceleration logG, and the EOS to estimate the missing parameters assuming hydrostatic equilibrium.

  • PeTop is provided, then this electron pressure at the top of the atmosphere will be used with the log gravitational acceleration logG, and the EOS to estimate the missing parameters assuming hydrostatic equilibrium.

  • If all of Pgas, Pe, Ptop, PeTop are omitted then Ptop will be estimated from the gas pressure in the FALC model at the temperature at the top boundary. The hydrostatic reconstruction will then continue as usual.

convertScales will substantially slow down this function due to the slow calculation of background opacities used to compute tauRef. If an atmosphere is constructed with a Geometric stratification, and an estimate of tauRef is not required before running the main RT module, then this can be set to False. All of these parameters can be provided as astropy Quantities, and will be converted in the constructor.

Parameters:
  • scale (ScaleType) – The type of stratification used along the z-axis.

  • depthScale (np.ndarray) – The z-coordinates used along the chosen stratification. The stratification is expected to start at the top of the atmosphere (closest to the observer), and descend along the observer’s line of sight.

  • temperature (np.ndarray) – Temperature structure of the atmosphere [K].

  • vlos (np.ndarray) – Velocity structure of the atmosphere along z [m/s].

  • vturb (np.ndarray) – Microturbulent velocity structure of the atmosphere [m/s].

  • ne (np.ndarray) – Electron density structure of the atmosphere [m-3].

  • hydrogenPops (np.ndarray, optional) – Detailed (per level) hydrogen number density structure of the atmosphere [m-3], 2D array [Nlevel, Nspace].

  • nHTot (np.ndarray, optional) – Total hydrogen number density structure of the atmosphere [m-3]

  • B (np.ndarray, optional.) – Magnetic field strength [T].

  • gammaB (np.ndarray, optional) – Co-altitude of magnetic field vector [radians].

  • chiB (np.ndarray, optional) – Azimuth of magnetic field vector (in x-y plane, from x) [radians].

  • lowerBc (BoundaryCondition, optional) – Boundary condition for incoming radiation at the minimal z coordinate (default: ThermalisedRadiation).

  • upperBc (BoundaryCondition, optional) – Boundary condition for incoming radiation at the maximal z coordinate (default: ZeroRadiation).

  • convertScales (bool, optional) – Whether to automatically compute tauRef and cmass for an atmosphere given in a stratification of m (default: True).

  • abundance (AtomicAbundance, optional) – An instance of AtomicAbundance giving the abundances of each atomic species in the given atmosphere, only used if the EOS is invoked. (default: DefaultAtomicAbundance)

  • logG (float, optional) – The log10 of the magnitude of gravitational acceleration [m/s2] (default: 2.44).

  • Pgas (np.ndarray, optional) – The gas pressure stratification of the atmosphere [Pa], optionally used by the EOS.

  • Pe (np.ndarray, optional) – The electron pressure stratification of the atmosphere [Pa], optionally used by the EOS.

  • Ptop (np.ndarray, optional) – The gas pressure at the top of the atmosphere [Pa], optionally used by the EOS for a hydrostatic reconstruction.

  • Petop (np.ndarray, optional) – The electron pressure at the top of the atmosphere [Pa], optionally used by the EOS for a hydrostatic reconstruction.

  • verbose (bool, optional) – Explain decisions made with the EOS to estimate missing parameters (if invoked) through print calls (default: False).

Raises:

ValueError – if incorrect arguments or unable to construct estimate missing parameters.

classmethod make_2d(height: ndarray, x: ndarray, temperature: ndarray, vx: ndarray, vz: ndarray, vturb: ndarray, ne: ndarray | None = None, nHTot: ndarray | None = None, B: ndarray | None = None, gammaB: ndarray | None = None, chiB: ndarray | None = None, xUpperBc: BoundaryCondition | None = None, xLowerBc: BoundaryCondition | None = None, zUpperBc: BoundaryCondition | None = None, zLowerBc: BoundaryCondition | None = None, abundance: AtomicAbundance | None = None, verbose=False)[source]

Constructor for 2D Atmosphere objects.

No provision for estimating parameters using hydrostatic equilibrium is provided, but one of ne, or nHTot can be omitted and inferred by use of the Wittmann equation of state. The atmosphere must be defined on a geometric stratification. All atmospheric parameters are expected in a 2D [z, x] array.

Parameters:
  • height (np.ndarray) – The z-coordinates of the atmospheric grid. The stratification is expected to start at the top of the atmosphere (closest to the observer), and descend along the observer’s line of sight.

  • x (np.ndarray) – The (horizontal) x-coordinates of the atmospheric grid.

  • temperature (np.ndarray) – Temperature structure of the atmosphere [K].

  • vx (np.ndarray) – x-component of the atmospheric velocity [m/s].

  • vz (np.ndarray) – z-component of the atmospheric velocity [m/s].

  • vturb (np.ndarray) – Microturbulent velocity structure [m/s].

  • ne (np.ndarray) – Electron density structure of the atmosphere [m-3].

  • nHTot (np.ndarray, optional) – Total hydrogen number density structure of the atmosphere [m-3].

  • B (np.ndarray, optional.) – Magnetic field strength [T].

  • gammaB (np.ndarray, optional) – Inclination (co-altitude) of magnetic field vector to the z-axis [radians].

  • chiB (np.ndarray, optional) – Azimuth of magnetic field vector (in x-y plane, from x) [radians].

  • xLowerBc (BoundaryCondition, optional) – Boundary condition for incoming radiation at the minimal x coordinate (default: PeriodicRadiation).

  • xUpperBc (BoundaryCondition, optional) – Boundary condition for incoming radiation at the maximal x coordinate (default: PeriodicRadiation).

  • zLowerBc (BoundaryCondition, optional) – Boundary condition for incoming radiation at the minimal z coordinate (default: ThermalisedRadiation).

  • zUpperBc (BoundaryCondition, optional) – Boundary condition for incoming radiation at the maximal z coordinate (default: ZeroRadiation).

  • convertScales (bool, optional) – Whether to automatically compute tauRef and cmass for an atmosphere given in a stratification of m (default: True).

  • abundance (AtomicAbundance, optional) – An instance of AtomicAbundance giving the abundances of each atomic species in the given atmosphere, only used if the EOS is invoked. (default: DefaultAtomicAbundance)

  • verbose (bool, optional) – Explain decisions made with the EOS to estimate missing parameters (if invoked) through print calls (default: False).

Raises:

ValueError – if incorrect arguments or unable to construct estimate missing parameters.

quadrature(Nrays: int | None = None, mu: Sequence[float] | None = None, wmu: Sequence[float] | None = None)[source]

Compute the angular quadrature for solving the RTE and Kinetic Equilibrium in a given atmosphere.

Procedure varies with dimensionality.

By convention muz is always positive, as the direction on this axis is determined by the toObs term that is used internally to the formal solver.

1D:

If a number of rays is given (typically 3 or 5), then the Gauss-Legendre quadrature for this set is used. If mu and wmu are instead given then these will be validated and used.

2+D:

If the number of rays selected is in the list of near optimal quadratures for unpolarised radiation provided by Stepan et al 2020 (A&A, 646 A24), then this is used. Otherwise an exception is raised.

The available quadratures are:

Points

Order

1

3

3

7

6

9

7

11

10

13

11

15

Parameters:
  • Nrays (int, optional) – The number of rays to use in the quadrature. See notes above.

  • mu (sequence of float, optional) – The cosine of the angle made between the between each of the set of rays and the z axis, only used in 1D.

  • wmu (sequence of float, optional) – The integration weights for each mu, must be provided if mu is provided.

Raises:

ValueError – on incorrect input.

rays(muz: float | Sequence[float], mux: float | Sequence[float] | None = None, muy: float | Sequence[float] | None = None, wmu: float | Sequence[float] | None = None, upOnly: bool = False)[source]

Set up the rays on the Atmosphere for computing the intensity in a particular direction (or set of directions).

If only the z angle is set then the ray is assumed in the x-z plane. If either muz or muy is omitted then this angle is inferred by normalisation of the projection.

By convention muz is always positive, as the direction on this axis is determined by the toObs term that is used internally to the formal solver.

Parameters:
  • muz (float or sequence of float, optional) – The angular projections along the z axis.

  • mux (float or sequence of float, optional) – The angular projections along the x axis.

  • muy (float or sequence of float, optional) – The angular projections along the y axis.

  • wmu (float or sequence of float, optional) – The integration weights for the given ray if J is to be integrated for angle set.

  • upOnly (bool, optional) – Whether to only configure boundary conditions for up-only rays. (default: False)

Raises:

ValueError – if the angular projections or integration weights are incorrectly normalised.

property tauRef: ndarray

np.ndarray Reference optical depth at 500 nm.

Type:

tauRef

unit_view()[source]

Returns a view over the contents of the Layout with the correct astropy.units.

property vlos: ndarray

np.ndarray z component of plasma velocity (present for all Ndim) [m/s]. Only available when Ndim==1`.

Type:

vz

property vx: ndarray

np.ndarray x component of plasma velocity (present for Ndim >= 2) [m/s].

Type:

vx

property vy: ndarray

np.ndarray y component of plasma velocity (present for Ndim == 3) [m/s].

Type:

vy

property vz: ndarray

np.ndarray z component of plasma velocity (present for all Ndim) [m/s]. Aliased to vlos when Ndim==1

Type:

vz

property x: ndarray

np.ndarray Ordinates of grid points along the x-axis (present for Ndim >= 2) [m].

Type:

x

property xLowerBc: BoundaryCondition

BoundaryCondition Boundary condition for the plane of minimal x-coordinate.

Type:

xLowerBc

property xUpperBc: BoundaryCondition

BoundaryCondition Boundary condition for the plane of maximal x-coordinate.

Type:

xUpperBc

property y: ndarray

np.ndarray Ordinates of grid points along the y-axis (present for Ndim == 3) [m].

Type:

y

property yLowerBc: BoundaryCondition

BoundaryCondition Boundary condition for the plane of minimal y-coordinate.

Type:

yLowerBc

property yUpperBc: BoundaryCondition

BoundaryCondition Boundary condition for the plane of maximal y-coordinate.

Type:

yUpperBc

property z: ndarray

np.ndarray Ordinates of grid points along the z-axis (present for all Ndim) [m].

Type:

z

property zLowerBc: BoundaryCondition

BoundaryCondition Boundary condition for the plane of minimal z-coordinate.

Type:

zLowerBc

property zUpperBc: BoundaryCondition

BoundaryCondition Boundary condition for the plane of maximal z-coordinate.

Type:

zUpperBc

class lightweaver.atmosphere.BoundaryCondition[source]

Bases: object

Base class for boundary conditions.

Defines the interface; do not use directly.

These attributes are only available after set_required_angles has been called.
mux

The mu_x to return from compute_bc (in order).

Type:

np.ndarray

muy

The mu_y to return from compute_bc (in order).

Type:

np.ndarray

muz

The mu_z to return from compute_bc (in order).

Type:

np.ndarray

indexVector

A 2D array of integer shape (mu, toObs) - where mu is the mu index on the associated atmosphere - relating each index of the second (Nrays) axis of a pair of (mu, toObs). Used to construct and destructure this array.

Type:

np.ndarray

compute_bc(atmos: Atmosphere, spect: LwSpectrum) ndarray[source]

Called when the radiation boundary condition is needed by the backend.

Parameters:
  • atmos (Atmosphere) – The atmospheric object in which to compute the radiation.

  • spect (LwSpectrum) – The computational spectrum object provided by the Context.

Returns:

result – This function needs to return a contiguous array of shape [Nwave, Nrays, Nbc], where Nwave is the number of wavelengths in the wavelength grid, Nrays is the number of rays in the angular quadrature (also including up/down directions) ordered as specified by the mux/y/z and indexVector variables on this object, Nbc is the number of spatial positions the boundary condition needs to be defined at ordered in a flattened [Nz, Ny, Nx] fashion. (dtype: <f8)

Return type:

np.ndarray

set_required_angles(mux, muy, muz, indexVector)[source]

The angles (and their ordering) to be used for this boundary condition (in the case of a callable)

class lightweaver.atmosphere.Layout(Ndim: int, x: ndarray, y: ndarray, z: ndarray, vx: ndarray, vy: ndarray, vz: ndarray, xLowerBc: BoundaryCondition, xUpperBc: BoundaryCondition, yLowerBc: BoundaryCondition, yUpperBc: BoundaryCondition, zLowerBc: BoundaryCondition, zUpperBc: BoundaryCondition, stratifications: Stratifications | None = None)[source]

Bases: object

Storage for basic atmospheric parameters whose presence is determined by problem dimensionality, boundary conditions and optional stratifications.

Ndim

Number of dimensions in model.

Type:

int

x

Ordinates of grid points along the x-axis (present for Ndim >= 2) [m].

Type:

np.ndarray

y

Ordinates of grid points along the y-axis (present for Ndim == 3) [m].

Type:

np.ndarray

z

Ordinates of grid points along the z-axis (present for all Ndim) [m].

Type:

np.ndarray

vx

x component of plasma velocity (present for Ndim >= 2) [m/s].

Type:

np.ndarray

vy

y component of plasma velocity (present for Ndim == 3) [m/s].

Type:

np.ndarray

vz

z component of plasma velocity (present for all Ndim) [m/s]. Aliased to vlos when Ndim==1

Type:

np.ndarray

xLowerBc

Boundary condition for the plane of minimal x-coordinate.

Type:

BoundaryCondition

xUpperBc

Boundary condition for the plane of maximal x-coordinate.

Type:

BoundaryCondition

yLowerBc

Boundary condition for the plane of minimal y-coordinate.

Type:

BoundaryCondition

yUpperBc

Boundary condition for the plane of maximal y-coordinate.

Type:

BoundaryCondition

zLowerBc

Boundary condition for the plane of minimal z-coordinate.

Type:

BoundaryCondition

zUpperBc

Boundary condition for the plane of maximal z-coordinate.

Type:

BoundaryCondition

property Noutgoing: int

Number of grid points at which the outgoing radiation is computed.

property Nspace: int

Number of spatial points present in the grid.

property Nx: int

Number of grid points along the x-axis.

property Ny: int

Number of grid points along the y-axis.

property Nz: int

Number of grid points along the z-axis.

property cmass

Alias to self.stratifications.cmass, if computed.

property dimensioned_shape

Tuple defining the shape to which the arrays of atmospheric paramters can be reshaped to be indexed in a 1/2/3D fashion.

dimensioned_unit_view() Layout[source]

Returns a view over the contents of Layout reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, and the correct astropy.units.

dimensioned_view() Layout[source]

Returns a view over the contents of Layout reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, as these are all stored under a flat scheme.

classmethod make_1d(z: ndarray, vz: ndarray, lowerBc: BoundaryCondition, upperBc: BoundaryCondition, stratifications: Stratifications | None = None) Layout[source]

Construct 1D Layout.

classmethod make_2d(x: ndarray, z: ndarray, vx: ndarray, vz: ndarray, xLowerBc: BoundaryCondition, xUpperBc: BoundaryCondition, zLowerBc: BoundaryCondition, zUpperBc: BoundaryCondition, stratifications: Stratifications | None = None) Layout[source]

Construct 2D Layout.

classmethod make_3d(x: ndarray, y: ndarray, z: ndarray, vx: ndarray, vy: ndarray, vz: ndarray, xLowerBc: BoundaryCondition, xUpperBc: BoundaryCondition, yLowerBc: BoundaryCondition, yUpperBc: BoundaryCondition, zLowerBc: BoundaryCondition, zUpperBc: BoundaryCondition, stratifications: Stratifications | None = None) Layout[source]

Construct 3D Layout.

property tauRef

Alias to self.stratifications.tauRef, if computed.

unit_view() Layout[source]

Returns a view over the contents of the Layout with the correct astropy.units.

class lightweaver.atmosphere.NoBc[source]

Bases: BoundaryCondition

Indicates no boundary condition on the axis because it is invalid for the current simulation. Used only by the backend.

class lightweaver.atmosphere.PeriodicRadiation[source]

Bases: BoundaryCondition

Periodic boundary condition. Commonly used on the x-axis in 2D simulations.

class lightweaver.atmosphere.ScaleType(value)[source]

Bases: Enum

Atmospheric scales used in the definition of 1D atmospheres to allow the correct conversion to a height based system. Options:

  • Geometric

  • ColumnMass

  • Tau500

class lightweaver.atmosphere.Stratifications(cmass: ndarray, tauRef: ndarray)[source]

Bases: object

Stores the optional derived z-stratifications of an atmospheric model.

cmass

Column mass [kg m-2].

Type:

np.ndarray

tauRef

Reference optical depth at 500 nm.

Type:

np.ndarray

dimensioned_unit_view(shape) Stratifications[source]

Makes an instance of Stratifications reshaped to the provided shape with the correct astropy.units for multi-dimensional atmospheres. For internal use.

Parameters:

shape (tuple) – Shape to reform the stratifications, provided by Layout.dimensioned_shape.

Returns:

stratifications – Reshaped stratifications with units.

Return type:

Stratifications

dimensioned_view(shape) Stratifications[source]

Makes an instance of Stratifications reshaped to the provided shape for multi-dimensional atmospheres. For internal use.

Parameters:

shape (tuple) – Shape to reform the stratifications, provided by Layout.dimensioned_shape.

Returns:

stratifications – Reshaped stratifications.

Return type:

Stratifications

unit_view() Stratifications[source]

Makes an instance of Stratifications with the correct astropy.units For internal use.

Returns:

stratifications – The same data with units applied.

Return type:

Stratifications

class lightweaver.atmosphere.ThermalisedRadiation[source]

Bases: BoundaryCondition

Thermalised radiation (blackbody) boundary condition. Commonly used for photospheric situations.

class lightweaver.atmosphere.ZeroRadiation[source]

Bases: BoundaryCondition

Zero radiation boundary condition. Commonly used for coronal situations.

lightweaver.atmosphere.get_top_pressure(eos: Wittmann, temp, ne=None, rho=None)[source]

Return a pressure for the top of atmosphere. For internal use.

In order this is deduced from:
  • the electron density ne, if provided

  • the mass density rho, if provided

  • the electron pressure present in FALC

Returns:

pressure – pressure IN CGS [dyn cm-2]

Return type:

float

lightweaver.atomic_model module

class lightweaver.atomic_model.AtomicContinuum(j: int, i: int)[source]

Bases: AtomicTransition

Base class for atomic continua.

alpha(wavelength: ndarray) ndarray[source]

Returns the cross-section as a function of wavelength

Parameters:

wavelength (np.ndarray) – The wavelengths at which to compute the cross-section

Returns:

alpha – The cross-section for each wavelength

Return type:

np.ndarray

property lambda0: float

The maximum (edge) wavelength at which this transition contributes [nm].

property lambda0_m: float

The maximum (edge) wavelength at which this transition contributes [m].

property lambdaEdge: float

The maximum (edge) wavelength at which this transition contributes [nm].

property minLambda: float

The minimum wavelength at which this transition contributes.

property polarisable: bool

Returns whether this continuum is polarisable, always False.

wavelength() ndarray[source]

The wavelength grid on which this continuum’s cross section is defined.

class lightweaver.atomic_model.AtomicLevel(E: float, g: float, label: str, stage: int, J: Fraction | None = None, L: int | None = None, S: Fraction | None = None)[source]

Bases: object

Description of atomic level in model atom.

E

Energy above ground state [cm-1]

Type:

float

g

Statistical weight of level

Type:

float

label

Name for level

Type:

str

stage

Ionisation of level with 0 being neutral

Type:

int

atom

AtomicModel that holds this level, will be initialised by the atom.

Type:

AtomicModel

J

Total quantum angular momentum.

Type:

Fraction, optional

L

Orbital angular momentum.

Type:

int, optional

S

Spin.

Type:

Fraction, optional

property E_SI

Returns E in Joule.

property E_eV

Returns E in electron volt.

property lsCoupling: bool

Returns whether the L-S coupling formalism can be applied to this level.

class lightweaver.atomic_model.AtomicLine(j: int, i: int, f: float, type: LineType, quadrature: LineQuadrature, broadening: LineBroadening, gLandeEff: float | None = None)[source]

Bases: AtomicTransition

Base class for atomic lines, holding their specialised information over transitions.

f

Oscillator strength.

Type:

float

type

Should the line be treated in PRD or CRD.

Type:

LineType

quadrature

Wavelength quadrature for integrating line properties over.

Type:

LineQuadrature

broadening

Object describing the broadening processes to be used in conjunction with the quadrature to generate the line profile.

Type:

LineBroadening

gLandeEff

Optionally override LS-coupling (if available for this transition), and just directly set the effective Lande g factor (if it isn’t).

Type:

float, optional

property Aji: float

Return the Einstein A coefficient for this line.

property Bij: float

Return the Einstein B_{ij} coefficient for this line.

property Bji: float

Return the Einstein B_{ji} coefficient for this line.

compute_phi(state: LineProfileState) LineProfileResult[source]

Compute the line profile, intended to be called from the backend.

property lambda0: float

Return the line rest wavelength [nm].

property lambda0_m: float

Return the line rest wavelength [m].

property overlyingContinuumLevel: AtomicLevel

Find the first overlying continuum level.

property polarisable: bool

Return whether sufficient information is available to compute full Stokes solutions for this line.

wavelength(vMicroChar=3000.0) ndarray[source]

Returns the wavelength grid for this transition based on the LineQuadrature.

Parameters:

vMicroChar (float, optional) – Characterisitc microturbulent velocity to assume when computing the line quadrature (default 3e3 m/s).

zeeman_components() ZeemanComponents | None[source]

Returns the Zeeman components of a line, if possible or None.

class lightweaver.atomic_model.AtomicModel(element: Element, levels: Sequence[AtomicLevel], lines: Sequence[AtomicLine], continua: Sequence[AtomicContinuum], collisions: Sequence[CollisionalRates])[source]

Bases: object

Container class for the complete description of a model atom.

element

The element or ion represented by this model.

Type:

Element

levels

The levels in use in this model.

Type:

list of AtomicLevel

lines

The atomic lines present in this model.

Type:

list of AtomicLine

continua

The atomic continua present in this model.

Type:

list of AtomicContinuum

collisions

The collisional rates present in this model.

Type:

list of CollisionalRates

property transitions: Sequence[AtomicTransition]

List of all atomic transitions present on the model.

vBroad(atmos: Atmosphere) ndarray[source]

Computes the atomic broadening velocity structure for a given atmosphere from the thermal motions and microturbulent velocity.

class lightweaver.atomic_model.AtomicTransition(j: int, i: int)[source]

Bases: object

Basic storage class for atomic transitions. Both lines and continua are derived from this.

property transId: Tuple[Element, int, int]

Unique identifier (transition ID) for transition (assuming one copy of each Element), used in creating a SpectrumConfiguration etc.

class lightweaver.atomic_model.ExplicitContinuum(j: int, i: int, wavelengthGrid: Sequence[float], alphaGrid: Sequence[float])[source]

Bases: AtomicContinuum

Specific version of atomic continuum with tabulated cross-section against wavelength. Interpolated using weno4. .. attribute:: wavelengthGrid

Wavelengths at which cross-section is tabulated [nm].

type:

list of float

alphaGrid

Tabulated cross-sections [m2].

Type:

list of float

alpha(wavelength: ndarray) ndarray[source]

Computes cross-section as a function of wavelength.

Parameters:

wavelength (np.ndarray) – Wavelengths at which to compute the cross-section [nm].

Returns:

alpha – Cross-section at associated wavelength.

Return type:

np.ndarray

property minLambda: float

The minimum wavelength at which this transition contributes.

wavelength() ndarray[source]

Returns the wavelength grid at which this transition needs to be computed to be correctly integrated. Specific handling is added to ensure that it is treated properly close to the edge.

class lightweaver.atomic_model.HydrogenicContinuum(j: int, i: int, NlambdaGen: int, alpha0: float, minWavelength: float)[source]

Bases: AtomicContinuum

Specific case of a Hydrogenic continuum, approximately falling off as 1/nu**3 towards higher frequencies (additional effects from Gaunt factor).

NlambaGen

The number of points to generate for the wavelength grid.

Type:

int

alpha0

The cross-section at the edge wavelength [m2].

Type:

float

minWavelength

The minimum wavelength below which this transition is assumed to no longer contribute [nm].

Type:

float

alpha(wavelength: ndarray) ndarray[source]

Computes cross-section as a function of wavelength.

Parameters:

wavelength (np.ndarray) – Wavelengths at which to compute the cross-section [nm].

Returns:

alpha – Cross-section at associated wavelength.

Return type:

np.ndarray

property minLambda: float

The minimum wavelength at which this transition contributes.

wavelength() ndarray[source]

Returns the wavelength grid at which this transition needs to be computed to be correctly integrated.

class lightweaver.atomic_model.LineProfileResult(phi: ndarray, aDamp: ndarray, Qelast: ndarray)[source]

Bases: object

Dataclass for returning the line profile and associated data that needs to be saved (damping parameter and elastic collision rate) from the frontend to the backend.

class lightweaver.atomic_model.LineProfileState(wavelength: ndarray, vlosMu: ndarray, atmos: Atmosphere, eqPops: SpeciesStateTable, default_voigt_callback: Callable[[ndarray, ndarray], ndarray], vBroad: ndarray | None = None)[source]

Bases: object

Dataclass used to communicate line profile calculations from the backend to the frontend whilst allowing the backend to provide an overrideable optimised voigt implementation for the default case.

wavelength

Wavelengths at which to compute the line profile [nm]

Type:

np.ndarray

vlosMu

Bulk velocity projected onto each ray in the angular integration scheme [m/s] in an array of [Nmu, Nspace].

Type:

np.ndarray

atmos

The associated atmosphere.

Type:

Atmosphere

eqPops

The associated populations for each species present in the simulation.

Type:

SpeciesStateTable

default_voigt_callback

Computes the Voigt profile for the default case, takes the damping parameter aDamp and broadening velocity vBroad as arguments, and returns the line profile phi (in this case phi_num in the tech report).

Type:

callable

vBroad

Cache to avoid recomputing vBroad every time. May be None.

Type:

np.ndarray, optional

class lightweaver.atomic_model.LineQuadrature[source]

Bases: object

Describes the wavelength quadrature to be used for integrating properties associated with a line.

doppler_units(line: AtomicLine) ndarray[source]

Return the quadrature in Doppler units.

wavelength(line: AtomicLine, vMicroChar: float = 3000.0) ndarray[source]

Return the quadrature in nm.

class lightweaver.atomic_model.LineType(value)[source]

Bases: Enum

Enum to show if the line should be treated in CRD or PRD.

class lightweaver.atomic_model.LinearCoreExpWings(qCore: float, qWing: float, Nlambda: int)[source]

Bases: LineQuadrature

RH-Style line quadrature, with approximately linear core spacing and exponential wing spacing, by using a function of the form q(n) = a*(n + (exp(b*n)-1)) with n in [0, N) satisfying the following conditions:

  • q[0] = 0

  • q[(N-1)/2] = qcore

  • q[N-1] = qwing.

If qWing <= 2 * qCore, linear grid spacing will be used for this transition.

doppler_units(line: AtomicLine) ndarray[source]

Return the quadrature in Doppler units.

wavelength(line: AtomicLine, vMicroChar=3000.0) ndarray[source]

Return the quadrature in nm.

class lightweaver.atomic_model.LinearQuadrature(Nlambda: int, deltaLambda: float)[source]

Bases: LineQuadrature

Simple linearly spaced wavelength grid. Primarily provided for CRTAF interaction.

Nlambdaint

The number of wavelength points in the wavelength grid (typically odd).

deltaLambdaint

The half-width of the grid (i.e. from core to one edge) [nm].

doppler_units(line: AtomicLine) ndarray[source]

Return the quadrature in Doppler units.

wavelength(line: AtomicLine, vMicroChar: float = 3000.0) ndarray[source]

Return the quadrature in nm.

class lightweaver.atomic_model.TabulatedQuadrature(wavelengthGrid: Sequence[float])[source]

Bases: LineQuadrature

Tabulated wavelength quadrature. Primarily provided for CRTAF interaction.

wavelengthGridSequence[float]

The wavelength sample points [nm].

doppler_units(line: AtomicLine) ndarray[source]

Return the quadrature in Doppler units.

wavelength(line: AtomicLine, vMicroChar: float = 3000.0) ndarray[source]

Return the quadrature in nm.

class lightweaver.atomic_model.VoigtLine(j: int, i: int, f: float, type: LineType, quadrature: LineQuadrature, broadening: LineBroadening, gLandeEff: float | None = None)[source]

Bases: AtomicLine

Specialised line profile for the default case of a Voigt profile.

compute_phi(state: LineProfileState) LineProfileResult[source]

Computes the line profile.

In the case of a VoigtLine the line profile simply uses the default_voigt_callback from the backend.

Parameters:

state (LineProfileState) – The information from the backend

Returns:

result – The line profile, as well as the damping parameter ‘a’ and and the broadening velocity.

Return type:

LineProfileResult

damping(atmos: Atmosphere, eqPops: SpeciesStateTable, vBroad: ndarray | None = None)[source]

Computes the damping parameter and elastic collision rate.

Parameters:
  • atmos (Atmosphere) – The atmosphere to consider.

  • eqPops (SpeciesStateTable) – The populations in this atmosphere.

  • vBroad (np.ndarray, optional) – The broadening velocity, will be used if passed, or computed using atom.vBroad if not.

Returns:

  • aDamp (np.ndarray) – The Voigt damping parameter.

  • Qelast (np.ndarray) – The rate of elastic collisions broadening the line – needed for PRD.

lightweaver.atomic_model.reconfigure_atom(atom: AtomicModel)[source]

Re-perform all atomic set up after modifying parameters.

lightweaver.atomic_set module

class lightweaver.atomic_set.AtomicState(model: AtomicModel, abundance: float, nStar: ndarray, nTotal: ndarray, detailed: bool = False, pops: ndarray | None = None, radiativeRates: Dict[Tuple[int, int], ndarray] | None = None)[source]

Bases: object

Container for the state of an atomic model during a simulation.

This hold both the model, as well as the simulations properties such as abundance, populations and radiative rates.

model

The python model of the atom.

Type:

AtomicModel

abundance

The abundance of the species as a fraction of H abundance.

Type:

float

nStar

The LTE populations of the species.

Type:

np.ndarray

nTotal

The total species population at each point in the atmosphere.

Type:

np.ndarray

detailed

Whether the species has detailed populations.

Type:

bool

pops

The NLTE populations for the species, if detailed is True.

Type:

np.ndarray, optional

radiativeRates

If detailed the radiative rates for the species will be present here, stored under (i, j) and (j, i) for each transition.

Type:

Dict[(int, int), np.ndarray], optional

dimensioned_unit_view(shape)[source]

Returns a view over the contents of AtomicState reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, and the correct astropy.units.

Parameters:

shape (tuple) – The shape to reshape to, this can be obtained from Atmosphere.structure.dimensioned_shape

Returns:

state – An instance of self with the arrays reshaped to the appropriate dimensionality.

Return type:

AtomicState

dimensioned_view(shape)[source]

Returns a view over the contents of AtomicState reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, as these are all stored under a flat scheme.

Parameters:

shape (tuple) – The shape to reshape to, this can be obtained from Atmosphere.structure.dimensioned_shape

Returns:

state – An instance of self with the arrays reshaped to the appropriate dimensionality.

Return type:

AtomicState

property element: Element

The element associated with this model.

property mass: float

The mass of the element associated with this model.

property n: ndarray

The NLTE populations, if present, or the LTE populations.

property name: str

The name of the element associated with this model.

set_n_to_lte()[source]

Reset the NLTE populations to LTE.

unit_view()[source]

Returns a view over the contents of the AtomicState with the correct astropy.units.

update_nTotal(atmos: Atmosphere)[source]

Update nTotal assuming either the abundance or nHTot have changed.

class lightweaver.atomic_set.AtomicStateTable(atoms: List[AtomicState])[source]

Bases: object

Container for AtomicStates.

The __getitem__ on this class is intended to be smart, and should work correctly with ints, strings, or Elements and return the associated AtomicState. This object is not normally constructed directly by the user, but will instead be interacted with as a means of transporting information to and from the backend.

dimensioned_unit_view(shape)[source]

Returns a view over the contents of AtomicStateTable reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, and the correct astropy.units.

Parameters:

shape (tuple) – The shape to reshape to, this can be obtained from Atmosphere.structure.dimensioned_shape

Returns:

state – An instance of self with the arrays reshaped to the appropriate dimensionality.

Return type:

AtomicStateTable

dimensioned_view(shape)[source]

Returns a view over the contents of AtomicStateTable reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, as these are all stored under a flat scheme.

Parameters:

shape (tuple) – The shape to reshape to, this can be obtained from Atmosphere.structure.dimensioned_shape

Returns:

state – An instance of self with the arrays reshaped to the appropriate dimensionality.

Return type:

AtomicStateTable

unit_view()[source]

Returns a view over the contents of the AtomicStateTable with the correct astropy.units.

class lightweaver.atomic_set.RadiativeSet(atoms: ~typing.List[~lightweaver.atomic_model.AtomicModel], abundance: ~lightweaver.atomic_table.AtomicAbundance = <lightweaver.atomic_table.AtomicAbundance object>)[source]

Bases: object

Used to configure the atomic models present in the simulation and then set up the global wavelength grid and initial populations. All atoms start passive.

Parameters:
  • atoms (list of AtomicModel) – The atomic models to be used in the simulation (active, detailed, and background).

  • abundance (AtomicAbundance, optional) – The abundance to be used for each species.

abundance

The abundances in use.

Type:

AtomicAbundance

elements

The elements present in the simulation.

Type:

list of Elements

atoms

Mapping from Element to associated model.

Type:

Dict[Element, AtomicModel]

passiveSet

Set of atoms (designmated by their Elements) set to passive in the simulation.

Type:

set of Elements

detailedStaticSet

Set of atoms (designmated by their Elements) set to “detailed static” in the simulation.

Type:

set of Elements

activeSet

Set of atoms (designmated by their Elements) set to active in the simulation.

Type:

set of Elements

property activeAtoms: List[AtomicModel]

List of AtomicModels set to active.

compute_eq_pops(atmos: Atmosphere, mols: MolecularTable | None = None, nlteStartingPops: Dict[Element, ndarray] | None = None)[source]

Compute the starting populations for the simulation with all NLTE atoms in LTE or otherwise using the provided populations.

Parameters:
  • atmos (Atmosphere) – The atmosphere for which to compute the populations.

  • mols (MolecularTable, optional) – Molecules to be included in the populations (default: None)

  • nlteStartingPops (Dict[Element, np.ndarray], optional) – Starting population override for any active or detailed static species.

Returns:

eqPops – The configured initial populations.

Return type:

SpeciesStatTable

compute_wavelength_grid(extraWavelengths: ndarray | None = None, lambdaReference=500.0) SpectrumConfiguration[source]

Compute the global wavelength grid from the current configuration of the RadiativeSet.

Parameters:
  • extraWavelengths (np.ndarray, optional) – Extra wavelengths to add to the global array [nm].

  • lambdaReference (float, optional) – If a difference reference wavelength is to be used then it should be specified here to ensure it is in the global array.

Returns:

spect – The configured wavelength grids needed to set up the backend.

Return type:

SpectrumConfiguration

property detailedAtoms: List[AtomicModel]

List of AtomicModels set to detailed static.

is_active(name: int | Tuple[int, int] | str | Element) bool[source]

Check if an atom (designated by int, (int, int), str, or Element) is active.

is_detailed(name: int | Tuple[int, int] | str | Element) bool[source]

Check if an atom (designated by int, (int, int), str, or Element) is passive.

is_passive(name: int | Tuple[int, int] | str | Element) bool[source]

Check if an atom (designated by int, (int, int), str, or Element) is passive.

iterate_lte_ne_eq_pops(atmos: Atmosphere, mols: MolecularTable | None = None, nlteStartingPops: Dict[Element, ndarray] | None = None, direct: bool = False) SpeciesStateTable[source]

Compute the starting populations for the simulation with all NLTE atoms in LTE or otherwise using the provided populations. Additionally computes a self-consistent LTE electron density.

Parameters:
  • atmos (Atmosphere) – The atmosphere for which to compute the populations.

  • mols (MolecularTable, optional) – Molecules to be included in the populations (default: None)

  • nlteStartingPops (Dict[Element, np.ndarray], optional) – Starting population override for any active or detailed static species.

  • direct (bool) – Whether to use the direct electron density solver (essentially, Lambda iteration for the fixpoint), this may require many more iterations than the Newton-Krylov solver used otherwise, but may also be more robust.

Returns:

eqPops – The configured initial populations.

Return type:

SpeciesStatTable

property passiveAtoms: List[AtomicModel]

List of AtomicModels set to passive.

set_active(*args: str)[source]

Set one (or multiple) atoms active.

set_detailed_static(*args: str)[source]

Set one (or multiple) atoms to detailed static

set_passive(*args: str)[source]

Set one (or multiple) atoms passive.

class lightweaver.atomic_set.SpeciesStateTable(atmosphere: Atmosphere, abundance: AtomicAbundance, atomicPops: AtomicStateTable, molecularTable: MolecularTable, molecularPops: List[ndarray], HminPops: ndarray)[source]

Bases: object

Container for the species populations in the simulation. Similar to AtomicStateTable but also holding the molecular populations and the atmosphere object.

The __getitem__ is intended to be smart, returning in order of priority on the name match (int, str, Element), H- populations, molecular populations, NLTE atomic populations, LTE atomic populations. This object is not normally constructed directly by the user, but will instead be interacted with as a means of transporting information to and from the backend.

atmosphere

The atmosphere object.

Type:

Atmosphere

abundance

The abundance of all species present in the atmosphere.

Type:

AtomicAbundance

atomicPops

The atomic populations state container.

Type:

AtomicStateTable

molecularTable

The molecules present in the simulation.

Type:

MolecularTable

molecularPops

The populations of each molecule in the molecularTable

Type:

list of np.ndarray

HminPops

H- ion populations throughout the atmosphere.

Type:

np.ndarray

dimensioned_unit_view()[source]

Returns a view over the contents of SpeciesStateTable reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, and the correct astropy.units.

dimensioned_view()[source]

Returns a view over the contents of SpeciesStateTable reshaped so all data has the correct (1/2/3D) dimensionality for the atmospheric model, as these are all stored under a flat scheme.

unit_view()[source]

Returns a view over the contents of the SpeciesStateTable with the correct astropy.units.

update_lte_atoms_Hmin_pops(atmos: Atmosphere, conserveCharge=False, updateTotals=False, maxIter=2000, quiet=False, tol=0.001)[source]

Under the assumption that the atmosphere has changed, update the LTE atomic populations and the H- populations.

Parameters:
  • atmos (Atmosphere) – The atmosphere object.

  • conserveCharge (bool) – Whether to conserveCharge and adjust the electron density in atmos based on the change in ionisation of the non-detailed species (default: False).

  • updateTotals (bool, optional) – Whether to update the totals of each species from the abundance and total hydrogen density (default: False).

  • maxIter (int, optional) – The maximum number of iterations to take looking for a stable solution (default: 2000).

  • quiet (bool, optional) – Whether to print information about the update (default: False)

  • tol (float, optional) – The tolerance of relative change at which to consider the populations converged (default: 1e-3)

class lightweaver.atomic_set.SpectrumConfiguration(radSet: RadiativeSet, wavelength: ndarray, models: List[AtomicModel], transWavelengths: Dict[Tuple[Element, int, int], ndarray], blueIdx: Dict[Tuple[Element, int, int], int], redIdx: Dict[Tuple[Element, int, int], int], activeTrans: Dict[Tuple[Element, int, int], bool], activeWavelengths: Dict[Tuple[Element, int, int], ndarray])[source]

Bases: object

Container for the configuration of common wavelength grid and species active at each wavelength.

radSet

The set of atoms involved in the creation of this simulation.

Type:

RadiativeSet

wavelength

The common wavelength array used for this simulation.

Type:

np.ndarray

models

The models for the active and detailed static atoms present in this simulation.

Type:

list of AtomicModel

transWavelengths

The local wavelength grid for each transition stored in a dictionary by transition ID.

Type:

Dict[(Element, i, j), np.ndarray]

blueIdx

The index at which each local grid starts in the global wavelength array.

Type:

Dict[(Element, i, j), int]

redIdx

The index at which each local grid has ended in the global wavelength array (exclusive, i.e. transWavelength = globalWavelength[blueIdx:redIdx]).

Type:

Dict[(Element, i, j), int]

activeTrans

Whether this transition is ever active (contributing in either an active or detailed static sense) over the range of wavelength.

Type:

Dict[(Element, i, j), bool]

activeWavelengths

A mask of the wavelengths at which this transition is active.

Type:

Dict[(Element, i, j), np.ndarray]

Properties
----------
NprdTrans

The number of PRD transitions present on the active transitions.

Type:

int

subset_configuration(wavelengths) SpectrumConfiguration[source]

Computes a SpectrumConfiguration for a sub-region of the global wavelength array.

This is typically used for computing a final formal solution on a single ray through the atmosphere. In this situation all lines are set to contribute throughout the entire grid, to avoid situations where small jumps in intensity occur from lines being cut off.

Parameters:

wavelengths (np.ndarray) – The grid on which to produce the new SpectrumConfiguration.

Returns:

spectrumConfig – The subset spectrum configuration.

Return type:

SpectrumConfiguration

lightweaver.atomic_set.chemical_equilibrium_fixed_ne(atmos: Atmosphere, molecules: MolecularTable, atomicPops: AtomicStateTable, abundance: AtomicAbundance) SpeciesStateTable[source]

Compute the molecular populations from the current atmospheric model and atomic populations.

This method assumes that the number of electrons bound in molecules is insignificant, and neglects this. Intended for internal use.

Parameters:
Returns:

state – The combined state object of atomic and molecular populations.

Return type:

SpeciesState

lightweaver.atomic_set.hminus_pops(atmos: Atmosphere, hPops: AtomicState) ndarray[source]

Compute the H- ion populations for a given atmosphere

Parameters:
  • atmos (Atmosphere) – The atmosphere object.

  • hPops (AtomicState) – The hydrogen populations state associated with atmos.

Returns:

HminPops – The H- populations.

Return type:

np.ndarray

lightweaver.atomic_set.lte_pops(atomicModel: AtomicModel, temperature: ndarray, ne: ndarray, nTotal: ndarray, nStar=None, debye: bool = True) ndarray[source]

Compute the LTE populations for a given atomic model under given thermodynamic conditions.

Parameters:
  • atomicModel (AtomicModel) – The atomic model to consider.

  • temperature (np.ndarray) – The temperature structure in the atmosphere.

  • ne (np.ndarray) – The electron density in the atmosphere.

  • nTotal (np.ndarray) – The total population of the species at each point in the atmosphere.

  • nStar (np.ndarray, optional) – An optional array to store the result in.

  • debye (bool, optional) – Whether to consider Debye shielding (default: True).1

Returns:

ltePops – The ltePops for the species.

Return type:

np.ndarray

lightweaver.atomic_table module

class lightweaver.atomic_table.AtomicAbundance(abundanceData: dict | None = None, abundDex=True, metallicity: float = 0.0)[source]

Bases: object

Container and accessor for atomic abundance data. This can be instantiated with a subset of atomic abundances, which will be used in conjunction with the the default values for non-specified elements.

Parameters:
  • abundanceData (dict, optional) – Contains the abundance data to override. For elements this should be dict[Element] = abundance, and for isotopes dict[Isotope] = isotope fraction.

  • abundDex (bool, optional) – Whether the supplied abundance is in dex (with Hydrogen abundance of 12.0) or in relative Hydrogen abundance (default: True i.e. in dex).

  • metallicity (float, optional) – Enhance the metallic abundance by a factor of 10**metallicity, (default: 0.0).

static apply_metallicity(abunds, metallicity)[source]

Used to adjust the metallicity of the abundances, by a fraction 10**metallicity.

compute_stats()[source]

Compute the total abundance (totalAbundance), mass per H atom (massPerH), and average mass per Element (avgMass).

convert_isotopes_to_abundances()[source]

Converts the isotope fractions to relative Hydrogen abundance.

static dex_to_decimal(abunds)[source]

Used to convert from absolute abundance in dex to relative fractional abundance.

get_primary_isotope(x: Element) Isotope[source]

Returns the Isotope with the highest abundance of a particular Element.

static load_default_abundance_data() dict[source]

Load the default abundances and convert to required format for AtomicAbundance class (i.e. dict where dict[Element] = abundance, and dict[Isotope] = isotope fraction).

class lightweaver.atomic_table.Element(Z: int)[source]

Bases: object

A simple value comparable description of an element (just proton number Z), that can be quickly and easily compared, whilst allowing access to things from the periodic table.

property mass

Returns the mass of the element in AMU.

property name

Returns the name of the element as a string.

class lightweaver.atomic_table.Isotope(N, Z)[source]

Bases: Element

A simple value comparable isotope description, inheriting from Element.

property element

Returns the underlying Element of which this isotope is a family member.

property element_mass

Returns the average mass of the element.

property element_name

Returns the name of the Element as a string.

property mass

Returns the mass of the isotope in AMU.

property name

Returns the name of the isotope as a string.

class lightweaver.atomic_table.KuruczPf(element: Element, abundance: float, Tpf: ndarray, pf: ndarray, ionPot: ndarray)[source]

Bases: object

Storage and functions relating to Bob Kurucz’s partition functions. Based on the data used in RH.

element

The element associated with this partition funciton.

Type:

Element

abundance

The abundance of this element.

Type:

float

Tpf

The temperature grid on which the partition function is defined.

Type:

np.ndarray

pf

The partition function data.

Type:

np.ndarray

ionPot

The ionisation potential of each level.

Type:

np.ndarray

fj(atmos: Atmosphere) Tuple[ndarray, ndarray][source]

Compute the fractional population of the species in each ionisation stage and partial derivative wrt n_e for each location in a given atmosphere.

Parameters:

atmos (Atmosphere) – The atmosphere in which to compute the populations.

Returns:

  • fj (np.ndarray) – The fractional populations [Nstage x Nspace].

  • dfj (np.ndarray) – The derivatives of the fractional populations [Nstage x Nspace].

fjk(atmos: Atmosphere, k: int) Tuple[ndarray, ndarray][source]

Compute the fractional population of the species in each ionisation stage and partial derivative wrt n_e at one point in a given atmosphere.

Parameters:
  • atmos (Atmosphere) – The atmosphere in which to compute the populations.

  • k (int) – The spatial index at which to compute the populations

Returns:

  • fj (np.ndarray) – The fractional populations [Nstage].

  • dfj (np.ndarray) – The derivatives of the fractional populations [Nstage].

lte_ionisation(atmos: Atmosphere) ndarray[source]

Compute the population of the species in each ionisation stage in a given atmosphere.

Parameters:

atmos (Atmosphere) – The atmosphere in which to compute the populations.

Returns:

pops – The LTE ionisation populations [Nstage x Nspace].

Return type:

np.ndarray

class lightweaver.atomic_table.KuruczPfTable(atomicAbundance: AtomicAbundance | None = None, kuruczPfPath: str | None = None)[source]

Bases: object

Container for all of the Kurucz partition function data, allowing different paths and AtomicAbundances to be used. Serves to construct the KuruczPf objects if used.

Parameters:
  • atomicAbundance (AtomicAbundance, optional) – The abundance data to use, if non-standard.

  • kuruczPfPath (str) – The path to the Kurucz parition function data in RH’s XDR format, if non-standard.

class lightweaver.atomic_table.PeriodicTableData[source]

Bases: object

Container and accessor for the periodic table data. Not intended to be instantiated by users, instead use the pre-instantiated PeriodicTable instance.

property elements

Return a sorted list of Elements by proton number Z.

classmethod get_isotopes(e: Element) List[Isotope][source]

Get all isotopes associated with a certain Element.

property isotopes

Return a sorted list of Isotopes by proton number Z.

property nuclides

Return a list of all nuclides (Elements and Isotopes).

lightweaver.atomic_table.load_periodic_table_data()[source]

Internal use function to load data from the AtomicMassesNames.pickle data file.

lightweaver.atomic_table.normalise_atom_name(n: str) str[source]

Normalises Element names to be two characters long with an uppercase first letter, and lower case second, or a space in the case of single character names.

Parameters:

n (str) – The name to normalise.

Returns:

result – The normalised name.

Return type:

str

lightweaver.barklem module

class lightweaver.barklem.Barklem[source]

Bases: object

Storage for all three Barklem cross-section cases and application via the get_active_cross_section function.

classmethod get_active_cross_section(atom: AtomicModel, line: AtomicLine, vals: Sequence[float]) Sequence[float][source]

Returns the cross section data for use in the Van der Waals collisional broadening routines. See:

  • Anstee & O’Mara 1995, MNRAS 276, 859-866

  • Barklem & O’Mara 1998, MNRAS 300, 863-871

  • Unsold:

    • Traving 1960, “Uber die Theorie der Druckverbreiterung von Spektrallinien”, p 91-97

    • Mihalas 1978, p. 282ff, and Table 9-1/

Returns:

result – Barklem cross-section, Barklem alpha, Helium contribution following Unsold (always 1.0)

Return type:

list of 3 float

exception lightweaver.barklem.BarklemCrossSectionError[source]

Bases: Exception

Raised if the Barklem cross-section cannot be applied to the atom in question.

class lightweaver.barklem.BarklemTable(path: str, neff0: Tuple[float, float])[source]

Bases: object

Storage for each table of Barklem data for Van der Waals approximation.

lightweaver.benchmark module

lightweaver.benchmark.benchmark(Niter=50, Nrep=3, verbose=True, writeConfig=True, warmUp=True)[source]

Benchmark the various SIMD implementations for Lightweaver’s formal solver and iteration functions.

Parameters:
  • Niter (int, optional) – The number of iterations to use for each scheme. (Default: 50)

  • Nrep (int, optional) – The number of repetitions to average for each scheme. (Default: 3)

  • verbose (bool, optional) – Whether to print information as the function runs. (Default: True)

  • writeConfig (bool, optional) – Whether to writ the optimal method to the user’s config file. (Default: True)

  • warmUp (bool, optional) – Whether to run a Context first (discarded) to ensure that all numba jit code is jitted and warm. (Default: True)

lightweaver.broadening module

class lightweaver.broadening.HydrogenLinearStarkBroadening[source]

Bases: StandardLineBroadener

Linear Stark broadening for the case of Hydrogen from Sutton 1978 (like RH).

broaden(atmos: Atmosphere, eqPops: SpeciesStateTable) ndarray[source]

The function that is called by LineBroadening.

Parameters:
  • atmos (Atmosphere) – The atmosphere in which to compute the broadening.

  • eqPops (SpeciesStateTable) – The populations to use for computing the broadening.

Returns:

broad – An array detailing the broadening at each location in the atmosphere [Nspace].

Return type:

np.ndarray

class lightweaver.broadening.LineBroadener[source]

Bases: object

Base class for broadening implementations. To be used if your broadener does something special and can’t just return an array.

class lightweaver.broadening.LineBroadening(natural: List[StandardLineBroadener], elastic: List[StandardLineBroadener], other: List[LineBroadener] | None = None)[source]

Bases: object

Standard component of AtomicLine to compute the broadening parameters in a flexible way.

For most Voigt-like situations, this should be usable without modifications, but this class can be inherited from to make substantial modifications.

Parameters:
  • natural (list of StandardLineBroadener) – List of broadening terms that are not elastic collisions (separated for PRD).

  • elastic (list of StandardLineBroadener) – List of elastic broadening terms.

  • other (list of LineBroadener, optional) – List of other broadening terms, not used by the VoigtLine by default, but existing to provide _options_ (default: None)

broaden(atmos: Atmosphere, eqPops: SpeciesStateTable) LineBroadeningResult[source]

Computes the broadening, this function is called by the AtomicLine object.

static compute_other_broadening(broadeners: List[LineBroadener] | None, atmos: Atmosphere, eqPops: SpeciesStateTable) List | None[source]

Returns a list of the computed broadening terms.

static sum_broadening_list(broadeners: List[StandardLineBroadener], atmos: Atmosphere, eqPops: SpeciesStateTable) ndarray | None[source]

Sums a list of StandardLineBroadeners.

class lightweaver.broadening.LineBroadeningResult(natural: ndarray, Qelast: ndarray, other: List | None = None)[source]

Bases: object

Result expected from instances of LineBroadening.broaden.

class lightweaver.broadening.MultiplicativeStarkBroadening(coeff: float)[source]

Bases: StandardLineBroadener

Simple expression for multiplicative Stark broadening, assumes that this can be expresed as a constant * ne.

broaden(atmos: Atmosphere, eqPops: SpeciesStateTable) ndarray[source]

The function that is called by LineBroadening.

Parameters:
  • atmos (Atmosphere) – The atmosphere in which to compute the broadening.

  • eqPops (SpeciesStateTable) – The populations to use for computing the broadening.

Returns:

broad – An array detailing the broadening at each location in the atmosphere [Nspace].

Return type:

np.ndarray

class lightweaver.broadening.QuadraticStarkBroadening(coeff: float)[source]

Bases: StandardLineBroadener

Lindholm theory result for Quadratic Stark broadening by electrons and singly ionised particles. Follows HM2014 pp. 238-239, uses C4 from Traving 1960 via RH.

broaden(atmos: Atmosphere, eqPops: SpeciesStateTable) ndarray[source]

The function that is called by LineBroadening.

Parameters:
  • atmos (Atmosphere) – The atmosphere in which to compute the broadening.

  • eqPops (SpeciesStateTable) – The populations to use for computing the broadening.

Returns:

broad – An array detailing the broadening at each location in the atmosphere [Nspace].

Return type:

np.ndarray

class lightweaver.broadening.RadiativeBroadening(gamma: float)[source]

Bases: StandardLineBroadener

Simple constant radiative broadening with coefficient gamma.

broaden(atmos: Atmosphere, eqPops: SpeciesStateTable) ndarray[source]

The function that is called by LineBroadening.

Parameters:
  • atmos (Atmosphere) – The atmosphere in which to compute the broadening.

  • eqPops (SpeciesStateTable) – The populations to use for computing the broadening.

Returns:

broad – An array detailing the broadening at each location in the atmosphere [Nspace].

Return type:

np.ndarray

class lightweaver.broadening.ScaledExponentBroadening(scaling: float, temperatureExp: float, hydrogenExp: float, electronExp: float)[source]

Bases: StandardLineBroadener

Broadening implementation following the CRTAF ScaledExponents recipe.

scaling * T**a * n_H(0)**b * n_e**c

scalingfloat

Scalar multiplication term

temperatureExpfloat

Temperature exponent

hydrogenExpfloat

Neutral (ground-state) hydrogen density exponent

electronExpfloat

Electron density exponent

broaden(atmos: Atmosphere, eqPops: SpeciesStateTable) ndarray[source]

The function that is called by LineBroadening.

Parameters:
  • atmos (Atmosphere) – The atmosphere in which to compute the broadening.

  • eqPops (SpeciesStateTable) – The populations to use for computing the broadening.

Returns:

broad – An array detailing the broadening at each location in the atmosphere [Nspace].

Return type:

np.ndarray

class lightweaver.broadening.StandardLineBroadener[source]

Bases: LineBroadener

Standard base class for broadening implementations. Unless you need to do something weird, inherit from this one.

class lightweaver.broadening.VdwApprox(vals: Sequence[float])[source]

Bases: StandardLineBroadener

Base class for van der Waals approximation using a list of coefficients (vals).

class lightweaver.broadening.VdwBarklem(vals: Sequence[float])[source]

Bases: VdwApprox

Implementation of the Barklem method for van der Waals broadening.

broaden(atmos: Atmosphere, eqPops: SpeciesStateTable) ndarray[source]

The function that is called by LineBroadening.

Parameters:
  • atmos (Atmosphere) – The atmosphere in which to compute the broadening.

  • eqPops (SpeciesStateTable) – The populations to use for computing the broadening.

Returns:

broad – An array detailing the broadening at each location in the atmosphere [Nspace].

Return type:

np.ndarray

class lightweaver.broadening.VdwUnsold(vals: Sequence[float])[source]

Bases: VdwApprox

Implementation of the Unsold method for van der Waals broadening.

broaden(atmos: Atmosphere, eqPops: SpeciesStateTable) ndarray[source]

The function that is called by LineBroadening.

Parameters:
  • atmos (Atmosphere) – The atmosphere in which to compute the broadening.

  • eqPops (SpeciesStateTable) – The populations to use for computing the broadening.

Returns:

broad – An array detailing the broadening at each location in the atmosphere [Nspace].

Return type:

np.ndarray

lightweaver.collisional_rates module

class lightweaver.collisional_rates.Ar85Cdi(j: int, i: int, cdi: Sequence[Sequence[float]])[source]

Bases: CollisionalRates

Collisional ionisation rates based on Arnaud & Rothenflug (1985, ApJS 60). Units remain in CGS as per paper.

class lightweaver.collisional_rates.Burgess(j: int, i: int, fudge: float = 1.0)[source]

Bases: CollisionalRates

Collisional ionisation from excited states from Burgess & Chidichimo (1983, MNRAS 203, 1269). Fudge parameter is dimensionless.

class lightweaver.collisional_rates.CE(j: int, i: int, temperature: Sequence[float], rates: Sequence[float])[source]

Bases: TemperatureInterpolationRates

Collisional (de-)excitation of neutrals by electrons. Units: s^-1 K^-1/2 m^3 Rate scales as sqrt(T) exp(DeltaE)

class lightweaver.collisional_rates.CH(j: int, i: int, temperature: Sequence[float], rates: Sequence[float])[source]

Bases: TemperatureInterpolationRates

Collisions with neutral hydrogen. Units: s^-1 m^3

class lightweaver.collisional_rates.CI(j: int, i: int, temperature: Sequence[float], rates: Sequence[float])[source]

Bases: TemperatureInterpolationRates

Collisional ionisation by electrons. Units: s^-1 K^-1/2 m^3 Rate scales as sqrt(T) exp(DeltaE)

class lightweaver.collisional_rates.CP(j: int, i: int, temperature: Sequence[float], rates: Sequence[float])[source]

Bases: TemperatureInterpolationRates

Collisional (de-)excitation by protons. Units: s^-1 m^3

class lightweaver.collisional_rates.ChargeExchangeNeutralH(j: int, i: int, temperature: Sequence[float], rates: Sequence[float])[source]

Bases: TemperatureInterpolationRates

Charge exchange with neutral hydrogen. Units: s^-1 m^3 Note: downward rate only.

class lightweaver.collisional_rates.ChargeExchangeProton(j: int, i: int, temperature: Sequence[float], rates: Sequence[float])[source]

Bases: TemperatureInterpolationRates

Charge exchange with protons. Units: s^-1 m^3 Note: upward rate only.

class lightweaver.collisional_rates.CollisionalRates(j: int, i: int)[source]

Bases: object

Base class for all CollisionalRates.

class lightweaver.collisional_rates.Omega(j: int, i: int, temperature: Sequence[float], rates: Sequence[float])[source]

Bases: TemperatureInterpolationRates

Collisional (de-)excitation of ions by electrons (dimensionless). Omega as in Seaton’s collision strength. Rate scales as 1/(sqrt(T)) exp(DeltaE).

class lightweaver.collisional_rates.TemperatureInterpolationRates(j: int, i: int, temperature: Sequence[float], rates: Sequence[float])[source]

Bases: CollisionalRates

Base class for rates defined by interpolating a coefficient on a temperature grid.

lightweaver.config module

lightweaver.config.get_config_path() str | None[source]

Returns the path to the lightweaverrc configuration file, or None if one cannot be found.

lightweaver.config.get_home_config_path() str[source]

Return the location where the user’s configuration data should be stored, whether it is currently present or not.

lightweaver.config.set_most_advanced_simd_impl()[source]

Picks the most advanced SIMD extensions (as detected by NumPy), that can be used on this system. This does not guarantee fastest (see the note in update_config_dict).

lightweaver.config.update_config_dict(configPath: str | None)[source]

Updates the configuration dict (lightweaver.ConfigDict), from the config file. If there is no config file, the defaults are used, and the most advanced instruction set is chosen for the SimdImpl. If the SimdImpl in the config file is too advanced for the current CPU, the maximum available is chosen.

Parameters:

configPath (str, optional) – The path to the config file, or None.

lightweaver.config.update_config_file(configPath: str)[source]

Updates the config file to the current values of the config dict.

Parameters:

configPath (str) – The path to the config file.

lightweaver.fal module

lightweaver.iterate_ctx module

class lightweaver.iterate_ctx.ConvergenceCriteria(ctx: Context, JTol: float, popsTol: float, rhoTol: float | None)[source]

Bases: object

Abstract base class for determining convergence inside iterate_ctx_se. A derived variant of this class will be instantiated by iterate_ctx_se dependent upon its arguments. The default implementation is DefaultConvergenceCriteria.

Parameters:
  • ctx (Context) – The context being iterated.

  • JTol (float) – The value of JTol passed to iterate_ctx_se.

  • popsTol (float) – The value of popsTol passed to iterate_ctx_se.

  • rhoTol (float or None) – The value of rhoTol passed to iterate_ctx_se.

is_converged(JUpdate: IterationUpdate, popsUpdate: IterationUpdate, prdUpdate: IterationUpdate | None) bool[source]

This function takes the IterationUpdate objects from ctx.formal_sol_gamma_matrices and ctx.stat_equil and optionally from ctx.prd_redistribute (or None). Should return a bool indicated whether the Context is sufficiently converged.

class lightweaver.iterate_ctx.DefaultConvergenceCriteria(ctx: Context, JTol: float, popsTol: float, rhoTol: float | None)[source]

Bases: ConvergenceCriteria

Default ConvergenceCriteria implementation. Usually sufficient for statistical equilibrium problems, but you may occasionally need to override this.

Parameters:
  • ctx (Context) – The context being iterated.

  • JTol (float) – The value of JTol passed to iterate_ctx_se.

  • popsTol (float) – The value of popsTol passed to iterate_ctx_se.

  • rhoTol (float or None) – The value of rhoTol passed to iterate_ctx_se.

is_converged(JUpdate: IterationUpdate, popsUpdate: IterationUpdate, prdUpdate: IterationUpdate | None) bool[source]

Returns whether the context is converged.

lightweaver.iterate_ctx.iterate_ctx_se(ctx: Context, Nscatter: int = 3, NmaxIter: int = 2000, prd: bool = False, JTol: float = 0.005, popsTol: float = 0.001, rhoTol: float | None = None, prdIterTol: float = 0.01, maxPrdSubIter: int = 3, printInterval: float = 0.2, quiet: bool = False, convergence: Type[ConvergenceCriteria] | None = None, returnFinalConvergence: bool = False)[source]

Iterate a configured Context towards statistical equilibrium solution.

Parameters:
  • ctx (Context) – The context to iterate.

  • Nscatter (int, optional) – The number of lambda iterations to perform for an initial estimate of J (default: 3).

  • NmaxIter (int, optional) – The maximum number of iterations (including Nscatter) to take (default: 2000).

  • prd (bool, optional) – Whether to perform PRD subiterations to estimate rho for PRD lines (default: False).

  • JTol (float, optional) – The maximum relative change in J from one iteration to the next (default: 5e-3).

  • popsTol (float, optional) – The maximum relative change in an atomic population from one iteration to the next (default: 1e-3).

  • rhoTol (float, optional) – The maximum relative change in rho for a PRD line on the final subiteration from one iteration to the next. If None, the change in rho will not be considered in judging convergence (default: None).

  • prdIterTol (float, optional) – The maximum relative change in rho for a PRD line below which PRD subiterations will cease for this iteration (default: 1e-2).

  • maxPrdSubIter (int, optional) – The maximum number of PRD subiterations to make, whether or not rho has reached the tolerance of prdIterTol (which isn’t necessary every iteration). (Default: 3)

  • printInterval (float, optional) – The interval between printing the update size information in seconds. A value of 0.0 will print every iteration (default: 0.2).

  • quiet (bool, optional) – Overrides any other print arguments and iterates silently if True. (Default: False).

  • convergence (derived ConvergenceCriteria class, optional) – The ConvergenceCriteria version to be used in determining convergence. Will be instantiated by this function, and the is_converged method will then be used. (Default: DefaultConvergenceCriteria).

  • returnFinalConvergence (bool, optional) – Whether to return the IterationUpdate objects used in the final convergence decision, if True, these will be returned in a list as the second return value. (Default: False).

Returns:

  • it (int) – The number of iterations taken.

  • finalIterationUpdates (List[IterationUpdate], optional) – The final IterationUpdates computed, if requested by returnFinalConvergence.

lightweaver.iteration_update module

class lightweaver.iteration_update.IterationUpdate(ctx: Context, crsw: float = 1.0, updatedJ: bool = False, dJMax: float = 0.0, dJMaxIdx: int = 0, updatedPops: bool = False, dPops: ~typing.List[float] = <factory>, dPopsMaxIdx: ~typing.List[int] = <factory>, ngAccelerated: bool = False, updatedNe: bool = False, dNeMax: float = 0.0, dNeMaxIdx: int = 0, updatedRho: bool = False, NprdSubIter: int = 0, dRho: ~typing.List[float] = <factory>, dRhoMaxIdx: ~typing.List[int] = <factory>, updatedJPrd: bool = False, dJPrdMax: ~typing.List[float] = <factory>, dJPrdMaxIdx: ~typing.List[int] = <factory>)[source]

Bases: object

Stores the results of an iteration of one of the backend functions, and determines how to format this for printing. All changes refer to relative change.

ctx

The context with which this update is associated.

Type:

Context

crsw

The current value of the collisional radiative switching parameter.

Type:

float

updatedJ

Whether the iteration affected the global J grid.

Type:

bool

dJMax

The maximum change in J.

Type:

float

dJMaxIdx

The index of the maximum change of J in a flattened array of J.

Type:

int

updatedPops

Whether the active atomic populations were modified by the iteration.

Type:

bool

dPops

The maximum change in each active population.

Type:

List[float]

dPopsMaxIdx

The location of the maximum change in each population in the flattened population array.

Type:

List[int]

ngAccelerated

Whether the atomic populations were modified by Ng Acceleration.

Type:

bool

updatedNe

Whether the electron density in the atmosphere was affected by the iteration.

Type:

bool

dNeMax

The maximum change in the electron density.

Type:

float

dNeMaxIdx

The location of the maximum change in the electron density array.

Type:

int

updatedRho

Whether the iteration affected the value of rhoPrd on PRD lines.

Type:

bool

NprdSubIter

The number of PRD sub-iterations taken (if multiple),

Type:

int

dRho

The maximum change in rho for each spectral line treated with PRD, in the order of the lines on each activeAtom. These values are repeated for each sub-iteration < NprdSubIter.

Type:

List[float]

dRhoMaxIdx

The location of the maximum change in rho for each PRD line.

Type:

List[int]

updatedJPrd

Whether the PRD iteration affected J.

Type:

bool

dJPrdMax

The maximum change in J during each PRD sub-iteration.

Type:

float

dJPrdMaxIdx

The location of the maximum change in J for each PRD sub-iteration.

Type:

int

dPopsMax

The maximum population change (including ne) over the iteration (read-only property).

Type:

float

dRhoMax

The maximum change in the PRD rho value for any line in the final subiteration (read-only property).

Type:

float

compact_representation()[source]

Produce a compact string representation of the object (similar to Lightweaver < v0.8).

lightweaver.molecule module

class lightweaver.molecule.MolecularTable(paths: List[str] | None = None)[source]

Bases: object

Stores a set of molecular models, can be indexed by model name to return the associated model (as a string).

class lightweaver.molecule.Molecule(filePath: str)[source]

Bases: object

Simple class for working with RH molecule definitions.

Parameters:

filePath (str) – Path from which to load molecular data. Use get_default_molecule_path for the path to the default RH distribution of molecule files (C2, CH, CN, CO, CaH, H2+, H2, H2O, HF, LiH, MgH, N2, NH, NO, O2, OH, SiO, TiO) all of which have the ‘.molecule’ extension.

lightweaver.multi module

class lightweaver.multi.MultiMetadata(name: str, logG: float)[source]

Bases: object

Metadata that is stored in a MULTI atmosphere, but doesn’t really belong in a Lightweaver atmosphere.

lightweaver.multi.read_multi_atmos(filename: str) Tuple[MultiMetadata, Atmosphere][source]

Load a MULTI atmosphere definition from a file for use in Lightweaver.

Parameters:

filename (str) – The path to load the file from.

Returns:

  • meta (MultiMetadata) – Additional metadata from the MULTI metadata that doesn’t appear in a Lightweaver atmosphere (name, gravitational acceleration used).

  • atmos (Atmosphere) – The atmosphere loaded into the standard Lightweaver format.

Raises:

ValueError – if file isn’t found, or cannot be parsed correctly.

lightweaver.nr_update module

lightweaver.nr_update.nr_post_update(self, fdCollisionRates=True, hOnly=False, timeDependentData=None, chunkSize=5, ngUpdate=None, printUpdate=None, extraParams=None)[source]

Compute the Newton-Raphson terms for updating the electron density through charge conservation. Is attached to the Context object. :param fdCollisionRates: Whether to use a finite difference approximation to the collisional

rates to find the population gradient WRT ne (default: True i.e. use finite-difference, if False, collisional rates are ignored for this process.)

Parameters:
  • hOnly (bool, optional) – Ignore atoms other than Hydrogen (the primary electron contributor) (default: False)

  • timeDependentData (dict, optional) – The presence of this argument indicates that the time-dependent formalism should be used. Should contain the keys ‘dt’ with a floating point timestep, and ‘nPrev’ with a list of population vectors from the start of the time integration step, in the order of the active atoms. This latter term can be obtained from the previous state provied by Context.time_dep_update.

  • chunkSize (int, optional) – Not currently used.

  • ngUpdate (bool, optional) – Whether to apply Ng Acceleration (default: None, to apply automatic behaviour), will only accelerate if the counter on the Ng accelerator has seen enough steps since the previous acceleration (set in Context initialisation).

  • printUpdate (bool, optional) – Whether to print information on the size of the update (default: None, to apply automatic behaviour).

  • extraParams (dict, optional) – Dict of extra parameters to be converted through the dict2ExtraParams function and passed onto the C++ core.

Returns:

dPops – The maximum relative change of any of the NLTE populations in the atmosphere.

Return type:

float

lightweaver.rh_atoms module

lightweaver.simd_management module

lightweaver.simd_management.filter_usable_simd_impls(implLibs: List[str]) List[str][source]

Filter a list of SimdImpl library names, returning those that are usable on the current machine.

lightweaver.simd_management.get_available_simd_suffixes() List[str][source]

Verifies the necessary flags against the features NumPy indicates are available, and returns a list of available LightweaverSimdImpls

lightweaver.utils module

exception lightweaver.utils.ConvergenceError[source]

Bases: Exception

Raised by some iteration schemes, can also be used in user code.

class lightweaver.utils.CrswIterator(initVal=1000.0)[source]

Bases: object

Basic iterator to be used for controlling the scale of the collisional radiative switching (of Hummer & Voels) multiplicative paramter. Can be inherited to provide different behaviour. By default starts from a factor of 1e3 and scales this factor by 0.1**(1.0/value) each iteration, as is the default behaviour in RH.

exception lightweaver.utils.ExplodingMatrixError[source]

Bases: Exception

Raised by the linear system matrix solver in the case of unsolvable systems.

class lightweaver.utils.InitialSolution(value)[source]

Bases: Enum

Initial solutions to use for atomic populations, either LTE, Zero radiation (not yet supported), or second order escape probability.

class lightweaver.utils.NgOptions(Norder: int = 0, Nperiod: int = 0, Ndelay: int = 0)[source]

Bases: object

Container for the options related to Ng acceleration. .. attribute:: Norder

The order of the extrapolation to use (default: 0, i.e. none).

type:

int, optional

Nperiod

The number of iterations to run between extrapolations.

Type:

int, optional

Ndelay

The number of iterations to run before starting acceleration.

Type:

int, optional

class lightweaver.utils.UnityCrswIterator[source]

Bases: CrswIterator

A specific case representing no collisional radiative switching (i.e. parameter always 1).

lightweaver.utils.air_to_vac(wavelength: ndarray) ndarray[source]

Convert air wavelength to vacuum.

Parameters:

wavelength (float or array-like or astropy.Quantity) – If no units then the wavelength is assumed to be in [nm], otherwise the provided units are used.

Returns:

result – The converted wavelength in [nm].

Return type:

float or array-like or astropy.Quantity

lightweaver.utils.check_shape_exception(a: ndarray, shape: int | Tuple[int], ndim: int | None = 1, name: str | None = 'array')[source]

Ensure that an array matches the expected number of dimensions and shape. Raise a ValueError if not, quoting the array’s name (if provided)

Parameters:
  • a (np.ndarray) – The array to verify.

  • shape (int or Tuple[int]) – The length (for a 1D array), or shape for a multi-dimensional array.

  • ndim (int, optional) – The expected number of dimensions (default: 1)

  • name (str, optional) – The name to in any exception (default: array)

lightweaver.utils.compute_contribution_fn(ctx, mu: int = -1, outgoing: bool = True) ndarray[source]

Computes the contribution function for all wavelengths in the simulation, for a chosen angular index.

Parameters:
  • ctx (Context) – A context with the full depth-dependent data (i.e. ctx.depthData.fill = True set before the most recent formal solution).

  • mu (Optional[int]) – The angular index to use (corresponding to the order of the angular quadratures in atmosphere), default: -1.

  • outgoing (Optional[bool]) – Whether to compute the contribution for outgoing or incoming radiation (wrt to the atmosphere). Default: outgoing==True, i.e. to observer.

Returns:

cfn – The contribution function in terms of depth and wavelength.

Return type:

np.ndarray

lightweaver.utils.compute_height_edges(ctx) ndarray[source]

Compute the edges of the height bins associated with the stratified altitude array used in a simulation, typically used in conjunction with a plot using pcolormesh.

Parameters:

ctx (Context) – The context from which to construct the height edges.

Returns:

heightEdges – The edges of the height bins.

Return type:

np.ndarray

lightweaver.utils.compute_radiative_losses(ctx) ndarray[source]

Compute the radiative gains and losses for each wavelength in the grid used by the context. Units of J/s/m3/Hz. Includes background/contributions from overlapping lines. Convention of positive => radiative gain, negative => radiative loss.

Parameters:

ctx (Context) – A context with full depth-dependent data (i.e. ctx.depthData.fill = True set before the most recent formal solution).

Returns:

loss – The radiative gains losses for each depth and wavelength in the simulation.

Return type:

np.ndarray

lightweaver.utils.compute_wavelength_edges(ctx) ndarray[source]

Compute the edges of the wavelength bins associated with the wavelength array used in a simulation, typically used in conjunction with a plot using pcolormesh.

Parameters:

ctx (Context) – The context from which to construct the wavelength edges.

Returns:

wlEdges – The edges of the wavelength bins.

Return type:

np.ndarray

lightweaver.utils.convert_specific_intensity(wavelength: ndarray, specInt: ndarray, outUnits) Quantity[source]

Convert a specific intensity between different units.

Parameters:
  • wavelength (np.ndarray or astropy.Quantity) – If no units are provided then this is assumed to be in nm.

  • specInt (np.ndarray or astropy.Quantity) – If no units are provided then this is assumed to be in J/s/m2/sr/Hz, the default for Lightweaver.

  • outUnits (str or astropy.Unit) – The units to convert specInt to e.g. ‘erg/s/cm2/sr/Angstrom’

Returns:

result – specInt converted to the desired units.

Return type:

astropy.Quantity

lightweaver.utils.filter_fs_iter_libs(libs: Sequence[str], exts: Sequence[str]) Sequence[str][source]

Filter a list of libraries (e.g. SimdImpl_{SimdType}.{pep3149}.so) with a valid collection of extensions. (As .so is a valid extension, we can’t just check the end of the file name).

lightweaver.utils.gaunt_bf(wvl, nEff, charge) float[source]

Gaunt factor for bound-free transitions, from Seaton (1960), Rep. Prog. Phys. 23, 313, as used in RH.

Parameters:
  • wvl (float or array-like) – The wavelength at which to compute the Gaunt factor [nm].

  • nEff (float) – Principal quantum number.

  • charge (float) – Charge of free state.

Returns:

result – Gaunt factor for bound-free transitions.

Return type:

float or array-like

lightweaver.utils.get_code_location()[source]

Returns the directory containing the Lightweaver Python source.

lightweaver.utils.get_data_path()[source]

Returns the location of the Lightweaver support data.

lightweaver.utils.get_default_molecule_path()[source]

Returns the location of the default molecules taken from RH.

lightweaver.utils.get_fs_iter_libs() Sequence[str][source]

Returns the paths of the default FsIterationScheme libraries usable on the current machine (due to available SIMD optimisations – these are detected by NumPy).

lightweaver.utils.integrate_line_losses(ctx, loss: ndarray, lines: AtomicLine | Sequence[AtomicLine], extendGridNm: float = 0.0) Sequence[ndarray] | ndarray[source]

Integrate the radiative gains and losses over the band associated with a line or list of lines. Units of J/s/m3. Includes background/contributions from overlapping lines. Convention of positive => radiative gain, negative => radiative loss.

Parameters:
  • ctx (Context) – A context with the full depth-dependent data (i.e. ctx.depthData.fill = True set before the most recent formal solution).

  • loss (np.ndarray) – The radiative gains/losses for each wavelength and depth computed by compute_radiative_losses.

  • lines (AtomicLine or list of AtomicLine) – The lines for which to compute losses.

  • extendGridNm (float, optional) – Set this to a positive value to add an additional point at each end of the integration range to include a wider continuum/far-wing contribution. Units: nm, default: 0.0.

Returns:

linesLosses – The radiative gain/losses per line at each depth.

Return type:

array or list of array

lightweaver.utils.planck(temp, wav)[source]

Planck black-body function B_nu(T) from wavelength.

Parameters:
  • temp (float or array-like) – Temperature [K]

  • wav (float or array-like) – The wavelength at which to compute B_nu [nm].

Returns:

result – B_nu(T)

Return type:

float or array-like

lightweaver.utils.sequence_repr(x: Sequence) str[source]

Uniform representation of arrays and lists as lists for use in round-tripping AtomicModels.

lightweaver.utils.vac_to_air(wavelength: ndarray) ndarray[source]

Convert vacuum wavelength to air.

Parameters:

wavelength (float or array-like or astropy.Quantity) – If no units then the wavelength is assumed to be in [nm], otherwise the provided units are used.

Returns:

result – The converted wavelength in [nm].

Return type:

float or array-like or astropy.Quantity

lightweaver.utils.view_flatten(x: ndarray) ndarray[source]

Return a flattened view over an array, will raise an Exception if it cannot be represented as a flat array without copy.

lightweaver.utils.voigt_H(a, v)[source]

Scalar Voigt profile.

Parameters:
  • a (float or array-like) – The a damping parameter to be used in the Voigt profile.

  • v (float or array-like) – The position in the line profile in Doppler units.

lightweaver.witt module

lightweaver.zeeman module

class lightweaver.zeeman.ZeemanComponents(alpha: ndarray, strength: ndarray, shift: ndarray)[source]

Bases: object

Storage for communicating the Zeeman components between functions, also shared with the backend, giving a slightly tighter contract than usual: all arrays must be contiguous and alpha must be of dtype np.int32.

lightweaver.zeeman.compute_zeeman_components(line: AtomicLine) ZeemanComponents | None[source]

Computes, if possible, the set of Zeeman components for an atomic line.

If gLandeEff is specified on the line, then basic three-component Zeeman splitting will be computed directly. Otherwise, if both the lower and upper levels of the line support LS-coupling (i.e. J, L, and S all specified, and J <= L + S), then the LS-coupling formalism is applied to compute the components of “anomalous” Zeeman splitting. If neither of these cases are fulfilled, then None is returned.

Parameters:

line (AtomicLine) – The line to attempt to compute the Zeeman components from.

Returns:

components – The Zeeman splitting components, if possible.

Return type:

ZeemanComponents or None

lightweaver.zeeman.effective_lande(line: AtomicLine)[source]

Computes the effective Lande g-factor for an atomic line.

lightweaver.zeeman.fraction_range(start: Fraction, stop: Fraction, step: Fraction = Fraction(1, 1)) Iterator[Fraction][source]

Works like range, but with Fractions. Does no checking, so best to make sure the range you’re asking for is sane and divides down properly.

lightweaver.zeeman.lande_factor(J: Fraction, L: int, S: Fraction) float[source]

Computes the Lande g-factor for an atomic level from the J, L, and S quantum numbers.

lightweaver.zeeman.zeeman_strength(Ju: Fraction, Mu: Fraction, Jl: Fraction, Ml: Fraction) float[source]

Computes the strength of a Zeeman component, following del Toro Iniesta (p. 137) albeit larger by a factor of 2 which is corrected by normalisation. Takes J upper and lower (u and l respectively), and M upper and lower.

lightweaver.LwCompiled module

lightweaver.LwCompiled.BC_to_enum(bc)

Returns the C++ enum associated with the type of python BoundaryCondition object.

class lightweaver.LwCompiled.BackgroundProvider

Bases: object

Base class for implementing background packages. Inherit from this to implement a new background scheme.

Parameters:
  • eqPops (SpeciesStateTable) – The populations of all species present in the simulation.

  • radSet (RadiativeSet) – The atomic models and configuration data.

  • wavelength (np.ndarray) – The array of wavelengths at which to compute the background.

compute_background(atmos, chi, eta, sca)

The function called by the backend to compute the background.

Parameters:
  • atmos (LwAtmosphere) – The atmospheric model.

  • chi (np.ndarray) – Array in which to store the background opacity [Nlambda, Nspace].

  • eta (np.ndarray) – Array in which to store the background emissivity [Nlambda, Nspace].

  • sca (np.ndarray) – Array in which to store the background scattering [Nlambda, Nspace].

class lightweaver.LwCompiled.BasicBackground

Bases: BackgroundProvider

Basic background implementation used by default in Lightweaver; equivalent to RH’s treatment i.e. H- opacity, CH, OH, H2 continuum opacities if present, continua from all passive atoms in the RadiativeSet, Thomson and Rayleigh scattering (from H and He).

bf_opacities(atmos, chi, eta)
compute_background(atmos, chiIn, etaIn, scaIn)

The function called by the backend to compute the background.

Parameters:
  • atmos (LwAtmosphere) – The atmospheric model.

  • chi (np.ndarray) – Array in which to store the background opacity [Nlambda, Nspace].

  • eta (np.ndarray) – Array in which to store the background emissivity [Nlambda, Nspace].

  • sca (np.ndarray) – Array in which to store the background scattering [Nlambda, Nspace].

rayleigh_scattering(atmos, sca)
class lightweaver.LwCompiled.FastBackground

Bases: BackgroundProvider

A faster implementation (due to C++ implementations) of BasicBackground supporting multiple threads.

compute_background(atmos, chiIn, etaIn, scaIn)

The function called by the backend to compute the background.

Parameters:
  • atmos (LwAtmosphere) – The atmospheric model.

  • chi (np.ndarray) – Array in which to store the background opacity [Nlambda, Nspace].

  • eta (np.ndarray) – Array in which to store the background emissivity [Nlambda, Nspace].

  • sca (np.ndarray) – Array in which to store the background scattering [Nlambda, Nspace].

class lightweaver.LwCompiled.LwAtmosphere

Bases: object

Storage for the C++ class, ensuring all of the arrays remained pinned from python. Usually constructed by the Context.

Parameters:
  • atmos (Atmosphere) – The python atmosphere object.

  • Nwavelengths (int) – The number of wavelengths used in the wavelength grid.

B

The magnetic field structure for the atmosphereic model (flat array).

Ndim

The dimensionality of the atmosphere.

Nrays

The number of rays in the angular quadrature.

Nspace

The number of points in the atmosphere.

Nx

The number of points along the x dimension.

Ny

The number of points along the y dimension.

Nz

The number of points along the z dimension.

chiB

Magnetic field azimuth

compute_bcs(spect)
configure_bcs(atmos)
cos2chi

cosine of 2*chi

cosGamma

cosine of gammaB

gammaB

Magnetic field co-altitude.

height

The z (altitude) grid.

mux

Cosine of angle with x-axis for each ray.

muy

Cosine of angle with y-axis for each ray.

muz

Cosine of angle with z-axis for each ray.

nHTot

Total hydrogen number density strucutre.

ne

The electron density structure of the atmospheric model (flat array).

pyAtmos
sin2chi

sine of 2*chi

temperature

The temperature structure of the atmospheric model (flat array).

update_projections()

Update all arrays of projected terms in the atmospheric model.

vlos

The z-velocity structure of the atmospheric model for 1D atmospheres (flat array).

vlosMu

The projected line of sight veloctity for each ray in the atmosphere.

vturb

Microturbelent velocity structure of the atmospheric model.

vx

The x-velocity structure of the atmospheric model (flat array).

vy

The y-velocity structure of the atmospheric model (flat array).

vz

The z-velocity structure of the atmospheric model (flat array).

wmu

Integration weights for angular quadrature.

x

The x grid.

y

The y grid.

z

The z grid.

class lightweaver.LwCompiled.LwAtom

Bases: object

Storage and access to computational atomic data used by backend. Sets up the computations transitions (LwTransition) present on the model. Instantiated by Context.

atomicModel

The atomic model object associated with this computational atom.

Type:

AtomicModel

modelPops

The population data for this species, in a python accessible form.

Type:

AtomicState

Parameters:
  • atom (AtomicModel) – The atomic model object associated with this computational atom.

  • atmos (LwAtmosphere) – The computational atmosphere to be used in the simulation.

  • eqPops (SpeciesStateTable) – The population of species present in the simulation.

  • spect (SpectrumConfiguration) – The configuration of the spectral grids.

  • background (LwBackground) – The background opacity terms, currently only used in the case of escape probability initial solution.

  • detailed (bool, optional) – Whether the atom is in detailed static or fully active mode (default: False).

  • initSol (InitialSolution, optional) – The initial solution to use for the atomic populations (default: LTE).

  • ngOptions (NgOptions, optional) – The Ng acceleration options (default: None)

  • conserveCharge (bool, optional) – Whether to conserve charge whilst setting populations from escape probability (ignored otherwise) (default: False).

  • fsIterSchemeProperties (dict, optional) – The properties of the FsIterScheme used as a dict, can be obtained from the FsIterSchemeManager. Only necessary keys are boolean defaultWlaGijStorage and defaultPerAtomStorage to determine allocation of wla, gij, eta, U, and chi on the underlying object. If not supplied, both of these default to True.

C

The collisional rates matrix [Nlevel, Nlevel, Nspace]. This is filled s.t. C_{ji} is C[i, j] to facilitate addition to Gamma.

Gamma

The Gamma iteration matrix [Nlevel, Nlevel, Nspace].

Nlevel

The number of levels in the atomic model.

Ntrans

The number of transitions in the atomic model.

atomicModel
compute_collisions(fillDiagonal=False)
compute_profiles(polarised=False)

Compute the line profiles for the spectral lines on the model.

Parameters:

polarised (bool, optional) – If True, and the lines are polarised, then the full Stokes line profiles will be computed, otherwise the scalar case will be computed (default: False).

element

The element identifier for this atomic model.

load_pops_rates_prd_from_state(prevState, popsOnly=False, preserveProfiles=False)
modelPops
n

The atomic populations (NLTE if in use) [Nlevel, Nspace].

nStar

The LTE populations for this species in this atmosphere [Nlevel, Nspace].

nTotal

The total number density of the model throughout the atmosphere.

set_pops_escape_probability(a, bg, conserveCharge=False, Niter=100)
setup_wavelength(la)

Initialise the wavelength dependent arrays for the wavelength at index la.

stages

The ionisation stage of each level of this model.

trans

List of computational transitions (LwTransition).

vBroad

The broadening velocity associated with this atomic model in this atmosphere.

class lightweaver.LwCompiled.LwBackground

Bases: object

Storage and driver for the background computations in Lightweaver. The storage is allocated and managed by this class, before being passed to C++ when necessary. This class is also responsible for calling the BackgroundProvider instance used (by default FastBackground with one thread).

chi

The background opacity [Nlambda, Nspace].

eta

The background eta [Nlambda, Nspace].

sca

The background scattering [Nlambda, Nspace].

update_background(atmos)

Recompute the background opacities, perhaps in the case where, for example, the atmospheric parameters have been updated.

Parameters:

atmos (LwAtmosphere) – The atmosphere in which to compute the background opacities and emissivities.

class lightweaver.LwCompiled.LwContext

Bases: object

Context that configures and drives the backend. Whilst the class is named LwContext (to avoid cython collisions) it is exposed from the lightweaver package as Context.

kwargs

A dictionary of all inputs provided to the context, stored as the under the argument names to __init__.

Type:

dict

eqPops

The populations of each species in the simulation.

Type:

SpeciesStateTable

conserveCharge

Whether charge is being conserved in the calculations.

Type:

bool

nrHOnly

Whether H is the only element included in charge conservation

Type:

bool

detailedAtomPrd

Whether PRD emission rho is computed for PRD lines with detailed static populations.

Type:

bool

crswCallback

The object controlling the value of the Collisional Radiative Switching term.

Type:

CrswIterator

crswDone

Indicates whether CRSW is done (i.e. the parameter has reached 1).

Type:

bool

Parameters:
  • atmos (Atmosphere) – The atmospheric structure object.

  • spect (SpectrumConfiguration) – The configuration of wavelength grids and active atoms/transitions.

  • eqPops (SpeciesStateTable) – The initial populations and storage for these populations during the simulation.

  • ngOptions (NgOptions, optional) – The parameters for Ng acceleration in the simulation (default: No acceleration).

  • initSol (InitialSolution, optional) – The starting solution for the population of all active species (default: LTE).

  • conserveCharge (bool, optional) – Whether to conserve charge in the simulation (default: False).

  • nrHOnly (bool, optional) – Only include hydrogen in charge conservation calculations (default: False).

  • hprd (bool, optional) – Whether to use the Hybrid PRD method to account for velocity shifts in the atmosphere (if PRD is used otherwise, then it is angle-averaged).

  • detailedAtomPrd (bool, optional) – Whether to compute the PRD emission coefficient rho for PRD lines on atoms with detailed static populations (default, True).

  • crswCallback (CrswIterator, optional) – An instance of CrswIterator (or derived thereof) to control collisional radiative swtiching (default: None for UnityCrswIterator i.e. no CRSW).

  • Nthreads (int, optional) – Number of threads to use in the computation of the formal solution, default 1.

  • backgroundProvider (BackgroundProvider, optional) – Implementation for the background, if non-standard. Must follow the BackgroundProvider interface.

  • formalSolver (str, optional) – Name of formalSolver registered with the FormalSolvers object.

  • interpFn (str, optional) – Name of interpolation function to use in the multi-dimensional formal solver. Must be registered with InterpFns.

Nthreads

The number of threads used by the formal solver. A new value can be assigned to this, and the necessary support structures will be automatically allocated.

activeAtoms

All active computational atomic models (LwAtom).

atmos

The atmospheric model storage object (LwAtmosphere).

background

The background storage object (LwBackground).

clear_ng()

Resets Ng acceleration objects on all active atoms.

compute_profiles(polarised=False)

Compute the line profiles for the spectral lines on all active and detailed atoms.

Parameters:

polarised (bool, optional) – If True, and the lines are polarised, then the full Stokes line profiles will be computed, otherwise the scalar case will be computed (default: False).

compute_rays(wavelengths=None, mus=None, stokes=False, updateBcs=None, upOnly=True, returnCtx=False, refinePrd=False, squeeze=True)

Compute the formal solution through a converged simulation for a particular ray (or set of rays). The wavelength range can be adjusted to focus on particular lines.

Parameters:
  • wavelengths (np.ndarray, optional) – The wavelengths at which to compute the solution (default: None, i.e. the original grid).

  • mus (float or sequence of float or dict) – The cosines of the angles between the rays and the z-axis to use, if a float or sequence of float then these are taken as muz. If a dict, then it is expected to be dictionary unpackable (double-splat) into atmos.rays, and can then be used for multi-dimensional atmospheres.

  • stokes (bool, optional) – Whether to compute a full Stokes solution (default: False).

  • updateBcs (Callable[[Atmosphere], None]) – Function to be applied to the Atmosphere (intended to update the boundary conditions if needed) before constructing the new Context for these rays. If a ray doesn’t intersect the boundary (i.e. x and y boundaries for muz == 1), then the boundary condition can be ignored.

  • upOnly (bool, optional) – Whether to only compute upgoing rays. Mostly affects the handling of boundary conditions for 2D atmospheres. (Default: True).

  • returnCtx (bool, optional) – Whether to return the Context used to compute the formal solution for these rays. If true, it will be returned as the second value. Default: False.

  • refinePrd (bool, optional) – Whether to update the rhoPrd term by reevaluating the scattering integral on the new wavelength grid. This can sometimes visually improve the final solution, but is quite computationally costly. (default: False i.e. not reevaluated, instead if the wavelength grid is different, rhoPrd is interpolated onto the new grid).

  • squeeze (bool, optional) – Whether to squeeze singular dimensions from the output array (default: True).

Returns:

intensity – The outgoing intensity for the chosen rays. If stokes=True then the first dimension indicates, in order, the I, Q, U, V components.

Return type:

np.ndarray

conserveCharge
static construct_from_state_dict_with(sd, atmos=None, spect=None, eqPops=None, ngOptions=None, initSol=None, conserveCharge=None, nrHOnly=None, detailedAtomPrd=None, hprd=None, preserveProfiles=False, fromScratch=False, backgroundProvider=None)

Construct a new Context informed by a state dictionary with changes provided to this function. This function is primarily aimed at making similar versions of a Context, as this can be duplicated much more easily by deepcopy or pickle.loads(pickle.dumps(ctx)). For example, wanting to replace the SpectrumConfiguration to run a different set of active atoms in the same atmospheric model. N.B. stateDict will not be deepcopied by this function – do that yourself if needed.

Parameters:
  • sd (dict) – The state dictionary, from Context.state_dict.

  • atmos (Atmosphere, optional) – Atmospheric model to use instead of the one present in stateDict.

  • spect (SpectrumConfiguration, optional) – Spectral configuration to use instead of the one present in stateDict.

  • eqPops (SpeciesStateTable, optional) – Species population object to use instead of the one present in stateDict.

  • ngOptions (NgOptions, optional) – Ng acceleration options to use.

  • initSol (InitialSolution, optional) – Initial solution to use, only matters if fromScratch is True.

  • conserveCharge (bool, optional) – Whether to conserve charge.

  • nrHOnly (bool, optional) – Whether to only consider Hydrogen in charge conservation calculations.

  • detailedAtomPrd (bool, optional) – Whether to compute the PRD emission coefficient rho for PRD lines on detailed atoms.

  • hprd (bool, optional) – Whether to use Hybrid-PRD.

  • preserveProfiles (bool, optional) – Whether to copy the current line profiles, or compute new ones (default: recompute).

  • fromScratch (bool, optional) – Whether to construct the new Context, but not make any modifications, such as copying profiles and rates.

  • backgroundProvider (BackgroundProvider, optional) – The background package to use instead of the one present in stateDict.

Returns:

ctx – The new context object for the simulation.

Return type:

Context

crswCallback
crswDone
depthData

Configuration and storage for full depth-dependent data of large parameters (LwDepthData).

detailedAtomPrd
detailedAtoms

All detailed static computational atomic models (LwAtom).

eqPops
formal_sol(upOnly=True, extraParams=None)

Compute the formal solution across all wavelengths (used by compute_rays). Only computes upgoing rays by default, which has implication on boundary conditions in 2D.

Parameters:
  • upOnly (bool, optional) – Only compute upgoing rays, (default: True)

  • extraParams (dict, optional) – Dict of extra parameters to be converted through the dict2ExtraParams function and passed onto the C++ core.

Returns:

update – An object representing the updates to the model. See IterationUpdate for details.

Return type:

IterationUpdate

formal_sol_gamma_matrices(fixCollisionalRates=False, lambdaIterate=False, printUpdate=None, extraParams=None)

Compute the formal solution across all wavelengths and fill in the Gamma matrix for each active atom, allowing the populations to then be updated using the radiative information.

Will use Nthreads for the formal solution.

Parameters:
  • fixCollisionalRates (bool, optional) – Whether to not recompute the collisional rates (default: False i.e. recompute them).

  • lambdaIterate (bool, optional) – Whether to use Lambda iteration (setting the approximate Lambda term to zero), may be useful in certain unstable situations (default: False).

  • printUpdate (bool, optional) – Whether to print the maximum relative change in J and any changes in CRSW (default: True). (Deprecated)

  • extraParams (dict, optional) – Dict of extra parameters to be converted through the dict2ExtraParams function and passed onto the C++ core.

Returns:

update – An object representing the updates to the model. See IterationUpdate for details.

Return type:

IterationUpdate

get_fs_iter_scheme_properties(fsIterScheme)
hprd

Whether PRD calculations are using the Hybrid PRD mode.

kwargs
nrHOnly
prd_redistribute(maxIter=3, tol=0.01, printUpdate=None, extraParams=None)

Update emission profile ratio rho by computing the scattering integral for each prd line. Does not affect the populations, interleave before each formal solution for a standard problem.

Parameters:
  • maxIter (int, optional) – The maximum number of iterations of updating rho to be taken (Default: 3).

  • tol (float, optional) – The default stopping tolerance for relative changes in rho. If the relative change in rho falls below this threshold then this function returns i.e. maxIter iterations do not need to be taken (Default: 1e-2).

  • printUpdate (bool, optional) – Whether to print information about the iteration process i.e. the size of the update to rho and the number of iterations taken (Default: True). Deprecated.

  • extraParams (dict, optional) – Dict of extra parameters to be converted through the dict2ExtraParams function and passed onto the C++ core.

Returns:

update – An object representing the updates to the model. See IterationUpdate for details.

Return type:

IterationUpdate

rel_diff_ng_accelerate(printUpdate=None)

Internal.

rel_diff_pops(printUpdate=None)

Internal.

set_formal_solver(formalSolver, inConstructor=False)

For internal use. Set the formal solver through the constructor.

set_fs_iter_scheme(fsIterScheme)
set_interp_fn(interpFn)

For internal use. Set the interpolation function through the constructor.

setup_stokes(recompute=False)

Configure the Context for Full Stokes radiative transfer.

Parameters:

recompute (bool, optional) – If previously called, and called again with recompute = True the line profiles will be recomputed.

single_stokes_fs(recompute=False, updateJ=False, upOnly=True, extraParams=None)

Compute a full Stokes formal solution across all wakelengths in the grid, setting up the Context first (it is rarely necessary to call setup_stokes directly).

The full Stokes formal solution is not currently multi-threaded, as it is usually only called once at the end of a simulation, however this could easily be changed.

Parameters:
  • recompute (bool, optional) – If previously called, and called again with recompute = True the line profiles will be recomputed. (Default: False)

  • updateJ (bool, optional) – Whether to update J on the Context during the calculation (Default: False)

  • upOnly (bool, optional) – Whether to compute the formal solver only for upgoing rays (used in final synthesis).

  • extraParams (dict, optional) – Dict of extra parameters to be converted through the dict2ExtraParams function and passed onto the C++ core.

Returns:

update – An object representing the updates to the model. See IterationUpdate for details.

Return type:

IterationUpdate

spect

The spectrum storage object (LwSpectrum).

stat_equil(printUpdate=None, chunkSize=20, extraParams=None)

Update the populations of active atoms using the current values of their Gamma matrices. This function solves the time-independent statistical equilibrium equations.

Parameters:
  • printUpdate (bool, optional) – Whether to print information on the size of the update (default: True). Deprecated.

  • chunkSize (int, optional) – Not currently used.

  • extraParams (dict, optional) – Dict of extra parameters to be converted through the dict2ExtraParams function and passed onto the C++ core.

Returns:

update – An object representing the updates to the model. See IterationUpdate for details.

Return type:

IterationUpdate

state_dict()

Return the state dictionary for the Context, which can be used to serialise the entire Context and/or reconstruct it.

time_dep_restore_prev_pops(prevTimePops)

Restore the populations to their state prior to the time-dependent updates for this timestep. Also resets I and J to 0. May be useful in cases where a problem was encountered.

Parameters:

prevTimePops (list of np.ndarray) – prevTimePops returned by time_dep_update.

time_dep_update(dt, prevTimePops=None, ngUpdate=None, printUpdate=None, chunkSize=20, extraParams=None)

Update the populations of active atoms using the current values of their Gamma matrices. This function solves the time-dependent kinetic equilibrium equations (ignoring advective terms). Currently uses a fully implicit (theta = 0) integrator.

Parameters:
  • dt (float) – The timestep length [s].

  • prevTimePops (list of np.ndarray or None) – The NLTE populations for each active atom at the start of the timestep (order matching that of Context.activeAtoms). This does not need to be provided the first time time_dep_update is called for a timestep, as if this parameter is None then this list will be constructed, and returned as the second return value, and can then be passed in again for additional iterations on a timestep.

  • ngUpdate (bool, optional) – Whether to apply Ng Acceleration (default: None, to apply automatic behaviour), will only accelerate if the counter on the Ng accelerator has seen enough steps since the previous acceleration (set in Context initialisation).

  • printUpdate (bool, optional) – Whether to print information on the size of the update (default: None, to apply automatic behaviour). Deprecated.

  • chunkSize (int, optional) – Not currently used.

  • extraParams (dict, optional) – Dict of extra parameters to be converted through the dict2ExtraParams function and passed onto the C++ core.

Returns:

  • update (IterationUpdate) – An object representing the updates to the model. See IterationUpdate for details.

  • prevTimePops (list of np.ndarray) – The input needed as prevTimePops if this function is to be called again for this timestep.

update_deps(temperature=True, ne=True, vturb=True, vlos=True, B=True, background=True, hprd=True, quiet=True)

Update various dependent parameters in the simulation after changes to different components. If a component has not been adjust then its associated argument can be set to False. By default, all standard dependent components are recomputed (e.g. projected velocities, line profiles, LTE populations, background terms).

Parameters:
  • temperature (bool, optional) – Whether the temperature has been modified.

  • ne (bool, optional) – Whether the electron density has been modified.

  • vturb (bool, optional) – Whether the microturbulent velocity has been modified.

  • vlos (bool, optional) – Whether the bulk velocity field has been modified.

  • B (bool, optional) – Whether the magnetic field has been modified.

  • background (bool, optional) – Whether the background needs updating.

  • hprd (bool, optional) – Whether the hybrid PRD terms need updating.

  • quiet (bool, optional) – Whether to print any update information from these functions (default: True).

update_hprd_coeffs()

Update the values of the H-PRD coefficients, this needs to be called if changes are made to the atmospheric structure to ensure that the interpolation parameters are correct.

update_projections()

Update all arrays of projected terms in the atmospheric model.

update_threads()

Internal.

class lightweaver.LwCompiled.LwDepthData

Bases: object

Simple object to lazily hold data that isn’t usually stored during the Formal Solution (full angularly dependent emissivity, opacity, and intensity at every point), due to the high memory cost. This is a part of the Context and doesn’t need to be instantiated directly.

I

Full depth dependent intensity [Nlambda, Nmu, Up/Down, Nspace].

chi

Full depth dependent opacity [Nlambda, Nmu, Up/Down, Nspace].

eta

Full depth dependent emissivity [Nlambda, Nmu, Up/Down, Nspace].

fill

Set this to True to fill the arrays, this will take care of allocating the space if not previously done.

class lightweaver.LwCompiled.LwFormalSolverManager

Bases: object

Storage and enumeration of the different formal solvers loaded for use in Lightweaver. There is no need to instantiate this class directly, instead there is a single instance of it instantiated as FormalSolvers, which should be used.

paths

The currently loaded paths.

Type:

list of str

names

The names of all available formal solvers, each of which can be passed to the Context constructor.

Type:

list of str

default_formal_solver(Ndim)

Returns the name of the default formal solver for a given dimensionality.

Parameters:

Ndim (int) – The dimensionality of the simulation.

Returns:

name – The name of the default formal solver.

Return type:

str

load_fs_from_path(path)

Attempt to load a formal solver, following the Lightweaver API, from a shared library at path.

Parameters:

path (str) – The path from which to load the formal solver.

names
paths
class lightweaver.LwCompiled.LwFsIterationManager

Bases: object

default_scheme()
default_scheme_name()
load_fns_from_path(path)
names
paths
scheme_properties(name)
class lightweaver.LwCompiled.LwInterpFnManager

Bases: object

Storage and enumeration of the different interpolation functions for multi-dimensional formal solvers loaded for use in Lightweaver. There is no need to instantiate this class directly, instead there is a single instance of it instantiated as InterpFns, which should be used.

paths

The currently loaded paths.

Type:

list of str

names

The names of all available interpolation functions, each of which can be passed to the Context constructor.

Type:

list of str

default_interp(Ndim)

Returns the name of the default interpolation function for a given dimensionality.

Parameters:

Ndim (int) – The dimensionality of the simulation.

Returns:

name – The name of the default interpolation function.

Return type:

str

load_interp_fn_from_path(path)

Attempt to load an interpolation function, following the Lightweaver API, from a shared library at path.

Parameters:

path (str) – The path from which to load the interpolation function.

names
paths
class lightweaver.LwCompiled.LwSpectrum

Bases: object

Storage and access to spectrum data used by backend. Instantiated by Context.

Parameters:
  • wavelength (np.ndarray) – The wavelength grid used in the simulation [nm].

  • Nrays (int) – The number of rays in the angular quadrature.

  • Nspace (int) – The number of points in the atmospheric model.

  • Noutgoing (int) – The number of outgoing point in the atmosphere, essentially max(Ny*Nx, Nx, 1), (when used in an array these elements will be ordered as a flattened array of [Ny, Nx]).

I

Intensity [J/s/m2/sr/Hz], shape is squeeze([Nlambda, Nmu, Noutgoing]).

J

Angle-averaged intensity [J/s/m2/sr/Hz], shape is squeeze([Nlambda, Nspace]).

Quv

Q, U and V Stokes parameters [J/s/m2/sr/Hz], shape is squeeze([3, Nlambda, Nmu, Noutgoing]).

interp_J_from_state(prevSpect)
setup_stokes()
wavelength

Wavelength grid used [nm].

class lightweaver.LwCompiled.LwTransition

Bases: object

Storage and access to transition data used by backend. Instantiated by Context.

Parameters:
  • trans (AtomicTransition) – The transition model object.

  • compAtom (LwAtom) – The computational atom to which this computational transition belongs.

  • atmos (LwAtmosphere) – The computational atmosphere in which this transition is to be used.

  • spect (SpectrumConfiguration) – The spectral configuration of the simulation.

transModel

The transition model object.

Type:

AtomicTransition

Aji

Einstein A for transition.

Bij

Einstein Bij for transition.

Bji

Einstein Bji for transition.

Nblue

Index into global wavelength grid where this transition’s local grid starts.

Qelast

The elastic collision rate for this transition in the atmosphere, needed for PRD.

Rij

Upwards radiative rates for the transition throughout the atmosphere.

Rji

Downward radiative rates for the transition throughout the atmosphere.

aDamp

The Voigt damping parameter for this transition in the atmosphere.

active

The active wavelength mask for this transition.

alpha

The wavelength-dependent cross-section for a continuum. AttributeError for lines.

atom
compute_phi()

Computes the line profile phi (phi_num in the technical report), by calling compute_phi on the line object. Provides a callback to the default Voigt implementation used in the backend. Does nothing if called on a continuum.

compute_polarised_profiles()

Compute the polarised line profiles (all of phi, phi_{Q, U, V}, and psi_{Q, U, V}) for a Voigt line, this currently doesn’t support non-standard line profile types, but could do so quite simply by following compute_phi. Does nothing if the transitions is a continuum or the line is not polarisable. By calling this and iterating the Context as usual, a field

i

Index of lower level.

iLevel

Access the lower level on the model object.

j

Index of upper level.

jLevel

Access the upper level on the model object.

lambda0

Line rest wavelength or continuum edge wavelength.

load_rates_prd_from_state(prevState, preserveProfiles=True)
phi

Numerical line profile. AttributeError for continua.

phiQ
phiU
phiV
polarisable

The polarisability of the transition, based on model data.

psiQ
psiU
psiV
recompute_gII()

Triggers lazy recalculation of gII for this line (if PRD).

rhoPrd

Ratio of emission to absorption profiles throughout the atmosphere, in the case of PRD lines.

transModel
type

The type of transition (Line or Continuum) as a str.

uv(la, mu, toObs, Uji, Vij, Vji)

Thin wrapper for computing U and V using the core. Must be called after atom.setup_wavelength(la), and Uji, Vij, Vji must tbe the proper size, as no verification is performed.

Parameters:
  • la (int) – The wavelength index at which to compute U and V.

  • mu (int) – The angle index at which to compute U and V.

  • Uji (np.ndarray) – Storage arrays for the result.

  • Vij (np.ndarray) – Storage arrays for the result.

  • Vji (np.ndarray) – Storage arrays for the result.

wavelength

The transition’s local wavelength grid.

wphi

Multiplicative inverse of integrated line profile at each location in the atmosphere, used to ensure terms based on integration across the entire line profile (e.g. in the Gamma matrix) are correctly normalised.

class lightweaver.LwCompiled.LwZeemanComponents

Bases: object

Stores the Zeeman components to be passed to the backend, only exists transiently.

class lightweaver.LwCompiled.RayleighScatterer

Bases: object

For computing Rayleigh scattering, used by BasicBackground.

scatter(wavelength, sca)