Promweaver API Documentation
Submodules
promweaver.bc_provider module
- class promweaver.bc_provider.DynamicContextPromBcProvider(ctx: Context)[source]
Bases:
PromBcProvider
Class for computing the incident radiation from a Lightweaver Context
- Parameters:
ctx (lw.Context) – The Lightweaver Context representing the chromospheric model to take into account, ready for compute_rays to be called.
- compute_I(wavelength: ndarray, mu: float | None)[source]
Return the intensity at the requested wavelengths and mu. If mu is None, return zeros.
- compute_multi_I(wavelength: ndarray, mus: List[float | None])[source]
Method to compute_I for multiple mus. It can be helpful to override this if your method has a significant overhead on compute_I that can be reduced by batching (e.g. this is very much the case for DynamicContextPromBcProvider).
- Returns:
I – The requested intensity, shape: [wavelength, mu]
- Return type:
np.ndarray
- class promweaver.bc_provider.PromBcProvider(specialised_multi_I: bool = False)[source]
Bases:
object
Abstract base class for providing the intensity to prominence BCs.
- Parameters:
specialised_multi_I (bool) – A bool that can be checked by boundary conditions to check if it is worth calling compute_multi_I (i.e. there are significant optimisations to doing so).
- compute_I(wavelength: ndarray, mu: float | None)[source]
Return the intensity at the requested wavelengths and mu. If mu is None, return zeros.
- compute_multi_I(wavelength: ndarray, mus: List[float | None])[source]
Method to compute_I for multiple mus. It can be helpful to override this if your method has a significant overhead on compute_I that can be reduced by batching (e.g. this is very much the case for DynamicContextPromBcProvider).
- Returns:
I – The requested intensity, shape: [wavelength, mu]
- Return type:
np.ndarray
- class promweaver.bc_provider.TabulatedPromBcProvider(wavelength, mu_grid, I, interp=None)[source]
Bases:
PromBcProvider
Class for tabulated incident radiation for prominence BC. Mu grids are interpolated linearly, wavelength grids using WENO4.
- Parameters:
wavelength (array-like) – The wavelengths at which the intensity is defined [nm]
mu_grid (array-like) – The mus at which the intensity is defined (typically 0.01 - 1.0)
I (array-like (2D)) – The intensity associated with the above 2 grids [wavelength, mu], SI units.
interp (interpolation function, optional) – The interpolation function to use over mu at each wavelength, if None, linear. If a function is provided then it must have the same signature as np.interp/weno4
- compute_I(wavelength: ndarray, mu: float | None)[source]
Return the intensity at the requested wavelengths and mu. If mu is None, return zeros.
- compute_multi_I(wavelength: ndarray, mus: List[float | None])[source]
Method to compute_I for multiple mus. It can be helpful to override this if your method has a significant overhead on compute_I that can be reduced by batching (e.g. this is very much the case for DynamicContextPromBcProvider).
- Returns:
I – The requested intensity, shape: [wavelength, mu]
- Return type:
np.ndarray
promweaver.compute_bc module
- promweaver.compute_bc.compute_falc_bc_ctx(active_atoms: List[str], atomic_models: List[AtomicModel] | None = None, prd: bool = False, vz: float | None = None, Nthreads: int = 1, quiet: bool = False, ctx_kwargs: dict | None = None) Context [source]
Configures and iterates a lightweaver Context with a FALC atmosphere, the selected atomic models, and active atoms. This can then be used with a DynamicContextPromBcProvider.
- promweaver.compute_bc.tabulate_bc(ctx: Context, wavelength: ndarray | None = None, mu_grid: ndarray | None = None, compute_rays_kwargs: dict | None = None, max_rays: int = 10)[source]
Computes the necessary data for a TabulatedPromBcProvider from a Context (e.g. from compute_falc_bc_ctx).
- Parameters:
ctx (lw.Context) – The Context that will be used to compute the boundary condition (should already be converged, as no iteration occurs here).
wavelength (array-like, optional) – The wavelength grid to use for the boundary condition table. Default: the one used in the Context.
mu_grid (array-like, optional) – The grid of viewing angles (in the form of mu) to be used for the boundary condition table. Default: np.linspace(0.01, 1.0, 100).
compute_rays_kwargs (dict, optional) – Any extra kwargs to be passed to compute_rays (e.g. refinePrd).
max_rays (int, optional) – The maximum number of rays to be computed simultaneously. The higher this number is the more efficient the process will be, but memory consumption scales linearly (with a big coefficient with this number), default 10. If set <= 0 the limit is taken to be infinite.
- Returns:
data – Dict with keys ‘wavelength’, ‘mu_grid’ and ‘I’ that can be splatted directly into a TabulatedPromBcProvider or pickled.
- Return type:
dict
promweaver.cone_prom_bc module
- class promweaver.cone_prom_bc.ConePromBc(projection: str, bc_provider: PromBcProvider, altitude_m: float, prom_upper_lower: str | None = None, Nphi_nodes: int = 50, Ncone_rays: int = 10)[source]
Bases:
PromBc
Complex boundary condition that supplies an averaged \(I_\nu\) distinct for each ray of the boundary condition. The primary use case for this is better handling complex velocity fields. In the filament view, this is simple as mu_z for the filament and the underlying atmosphere are aligned. For the prominence case, we instead need to average over 1/4 of a cone (in the sense of a polar angle \(phi\) in the prominence frame) and its angles of intersection with the solar surface. This slightly breaks the axisymmetric assumptions of one-dimensional models, but the averaging ensures that it’s correct from an energy standpoint. In final_synthesis mode (by requesting a single compute_rays from the model, or setting final_synthesis to true on the lw.Atmosphere used in the model), radiation is directly sampled along the ray path, as the converged populations have been determined.
This style of boundary condition is similar to the Gouttebroze 2005, and based on those described in Jenkins, Osborne & Keppens 2022 (in prep.).
- Parameters:
projection (str) – filament or prominence
bc_provider (PromBcProvider) – The provider to be used for computing the necessary intensity values.
altitude_m (float) – The altitude of the prominence above the solar surface [m]
prom_upper_lower (str, optional) – Whether this is the ‘upper’ or ‘lower’ z-boundary for a prominence (not used in filament cases).
Nphi (int, optional) – The number of trapezoidal phi rays to be used in each cone. Default: 50.
Ncone_rays (int, optional) – The number of Gauss-Legendre rays per cone to supersample the range of mus associated with this cone. This prevents sharp changes in intensity when one ray (out of few) starts to miss the Sun with increasing altitude. Default: 10.
- compute_bc(atmos, spect)[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
- class promweaver.cone_prom_bc.ConeRadParams(wavelength: ndarray, z_mu_cone: ndarray, w_mu_cone: ndarray, phi_nodes: ndarray, w_phi_nodes: ndarray, alt: float, projection: str, provider: PromBcProvider)[source]
Bases:
object
Dataclass to compact the parameters to compute_I_cone which were getting out of hand.
promweaver.iso_prom module
- class promweaver.iso_prom.IsoPromModel(projection, temperature, pressure, thickness, vturb, altitude, active_atoms: List[str], detailed_atoms: List[str] = None, detailed_pops: Dict[str, ndarray] | None = None, atomic_models=None, Nhalf_points=45, Nrays=3, Nthreads=1, prd=False, vlos: float | None = None, vrad: float | None = None, ctx_kwargs=None, BcType: Type[PromBc] | None = None, bc_kwargs=None, bc_provider=None, add_vertical_ray: bool = False, add_extra_rays: Dict[str, List[float] | Tuple[float] | ndarray] = None)[source]
Bases:
PromModel
Class for “Iso” prominence simulations. Iso implies isothermal and isobaric.
- Parameters:
projection (str) – Whether the object is to be treated as a “filament” or a “prominence”.
temperature (float) – The temperature of the prominence [K].
pressure (float) – The pressure of the prominence [Pa].
thickness (float) – The thickness of the prominence [m].
vturb (float) – The microturbulent velocity inside the prominence [m/s].
altitude (float) – The altitude of the prominence above the solar surface [m].
active_atoms (list of str) – The element names to make “active” i.e. consider in non-LTE.
detailed_atoms (list of str, optional) – The element names to make “detailed static” i.e. consider in detail radiatively. If ‘H’ is set detailed then the model will be run with a fixed electron density (set to balance pressure) and iterative charge conservation disabled. Provide a compatible bc_provider as these will not be taken account in the default one computed on the fly.
detailed_pops (dict, optional) – A mapping from names in detailed_atoms to an array of their population distribution of shape suitable for the model.
atomic_models (list of lw.AtomicModels, optional) – The atomic models to use, a default set will be chosen if none are specified.
Nhalf_points (int, optional) – The number of points in half of the slab. Default: 45.
Nrays (int, optional) – The number of Gauss-Legendre angular quadrature rays to use in the model. Default: 3. This number will need to be set higher (e.g. 10) if using the ConePromBc.
Nthreads (int, optional) – The number of CPU threads to use when solving the radiative transfer equations. Default: 1.
prd (bool, optional) – Whether to consider the effects of partial frequency redistribution. Default: False.
vlos (float, optional) – The z-projected velocity to apply to the prominence model. Default: None, i.e. 0.
vrad (float, optional) – The Doppler-dimming radial velocity to apply. Note that for filaments this is the same as vlos and that should be used instead. Not fully supported in boundary conditions yet, (i.e. you should interpolate the wavelength grid first). Default: None, i.e. 0.
ctx_kwargs (dict, optional) – Extra kwargs to be passed when constructing the Context.
BcType (Constructor for a type of PromBc, optional) – The base type to be used for constructing the boundary conditions. Default: UniformJPromBc.
bc_kwargs (dict, optional) – Extra kwargs to be passed to the construction of the boundary conditions.
bc_provider (PromBcProvider, optional) – The provider to use for computing the radiation in the boundary conditions. Default: DynamicContextPromBcProvider using an lw.Context configured to match the current model. Note that the default is note very performant, but is convenient for experimenting. When running a grid of models, consider creating a TabulatedPromBcProvider using compute_falc_bc_ctx and tabulate_bc, since the default performs quite a few extra RT calculations.
add_vertical_ray (bool, optional) – Whether to add a non-weighted vertical ray to the model. This doesn’t participate in energy balance but allows for directly extracting mu=1, without needing to do a boundary condition. Default: False. Equivalent to add_extra_mus={“muz”: [1.0], “mux”: [0.0]}.
add_extra_rays (dict, optional) – Extra rays to add to the quadrature along which the solution should be sampled. These rays will not have a weight for integration, but can serve as output for different viewing angles. Should be a dict with keys muz and mux as iterables.
- iterate_se(*args, update_model_kwargs: dict | None = None, **kwargs)[source]
Iterate the context towards statistical equilibrium solution. Slightly modified variant of the core lw.iterate_ctx_se to add and handle the update_model callback.
- Parameters:
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).
update_model (callable(PromModel, bool, **kwargs), optional) – A function to use to update model parameters based on iteration (e.g. scale populations to maintain pressure), receives the model, and whether it should print based on the current iteration progress loop.
update_model_kwargs (dict, optional) – Extra kwargs to be passed to update_model.
- Returns:
it (int) – The number of iterations taken.
finalIterationUpdates (List[IterationUpdate], optional) – The final IterationUpdates computed, if requested by returnFinalConvergence.
promweaver.stratified_prom module
- class promweaver.stratified_prom.StratifiedPromModel(projection, z, temperature, vlos, vturb, altitude, active_atoms: List[str], detailed_atoms: List[str] = None, detailed_pops: Dict[str, ndarray] | None = None, pressure=None, nh_tot=None, ne=None, initial_ionisation_fraction=1.0, conserve_charge=True, conserve_pressure=True, atomic_models=None, Nrays=3, Nthreads=1, prd=False, vrad: float | None = None, ctx_kwargs=None, BcType: Type[PromBc] | None = None, bc_kwargs=None, bc_provider=None, add_vertical_ray: bool = False, add_extra_rays: Dict[str, List[float] | Tuple[float] | ndarray] = None, temperature_threshold: float | None = None)[source]
Bases:
PromModel
Class for stratified prominence simulations. These fixed stratifications are usually produced by simulations (such as those with AMRVAC), or as intermediate steps during inversions.
- Parameters:
projection (str) – Whether the object is to be treated as a “filament” or a “prominence”.
z (array) – Height stratification (usual RT upside-down orientation) [m].
temperature (array) – The temperature stratification of the prominence [K].
vlos (array) – The z-projected velocity stratification to apply to the prominence model [m/s].
vturb (array) – The microturbulent velocity inside the prominence [m/s].
altitude (float) – The altitude of the prominence above the solar surface [m].
active_atoms (list of str) – The element names to make “active” i.e. consider in non-LTE.
detailed_atoms (list of str, optional) – The element names to make “detailed static” i.e. consider in detail radiatively. If ‘H’ is set detailed then the model will be run with a fixed electron density (set to balance pressure) and iterative charge conservation disabled. Provide a compatible bc_provider as these will not be taken account in the default one computed on the fly.
detailed_pops (dict, optional) – A mapping from names in detailed_atoms to an array of their population distribution of shape suitable for the model.
pressure (array, optional) – The pressure stratification of the prominence [Pa]. Note: pressure should be provided or at least 2 of pressure, nh_tot, and ne. If all three are provided, they are assumed to be self consistent as per P = (ne + lw.DefaultAtomicAbundance.totalAbundance * nh_tot) * k_B * T.
nh_tot (array, optional) – The hydrogen number density stratification of the prominence [m-3]. Note: pressure should be provided or at least 2 of pressure, nh_tot, and ne. If all three are provided, they are assumed to be self consistent as per P = (ne + lw.DefaultAtomicAbundance.totalAbundance * nh_tot) * k_B * T.
ne (array, optional) – Electron density stratification [m-3] Note: pressure should be provided or at least 2 of pressure, nh_tot, and ne. If all three are provided, they are assumed to be self consistent as per P = (ne + lw.DefaultAtomicAbundance.totalAbundance * nh_tot) * k_B * T.
initial_ionisation_fraction (float or array, optional) – The initial ionisation fraction to assume when computing ne from pressure if not provided. Default: 1.0
conserve_charge (bool, optional) – Whether to perform additional iterative procedure to balance charge (will affect pressure). Default: True
conserve_pressure (bool, optional) – Whether to perform additional iterations to rebalance pressure changes from variation in ionisation to match the initial stratification (pressure if provided, or (nh_tot + ne) * k_B * T otherwise). Default: True
atomic_models (list of lw.AtomicModels, optional) – The atomic models to use, a default set will be chosen if none are specified.
Nrays (int, optional) – The number of Gauss-Legendre angular quadrature rays to use in the model. Default: 3. This number will need to be set higher (e.g. 10) if using the ConePromBc.
Nthreads (int, optional) – The number of CPU threads to use when solving the radiative transfer equations. Default: 1.
prd (bool, optional) – Whether to consider the effects of partial frequency redistribution. Default: False.
vrad (float, optional) – The Doppler-dimming radial velocity to apply. Note that for filaments this is the same as vlos and that should be used instead. Not fully supported in boundary conditions yet, (i.e. you should interpolate the wavelength grid first). Default: None, i.e. 0.
ctx_kwargs (dict, optional) – Extra kwargs to be passed when constructing the Context.
BcType (Constructor for a type of PromBc, optional) – The base type to be used for constructing the boundary conditions. Default: UniformJPromBc.
bc_kwargs (dict, optional) – Extra kwargs to be passed to the construction of the boundary conditions.
bc_provider (PromBcProvider, optional) – The provider to use for computing the radiation in the boundary conditions. Default: DynamicContextPromBcProvider using an lw.Context configured to match the current model. Note that the default is note very performant, but is convenient for experimenting. When running a grid of models, consider creating a TabulatedPromBcProvider using compute_falc_bc_ctx and tabulate_bc, since the default performs quite a few extra RT calculations.
add_vertical_ray (bool, optional) – Whether to add a non-weighted vertical ray to the model. This doesn’t participate in energy balance but allows for directly extracting mu=1, without needing to do a boundary condition. Default: False. Equivalent to add_extra_mus={“muz”: [1.0], “mux”: [0.0]}.
add_extra_rays (dict, optional) – Extra rays to add to the quadrature along which the solution should be sampled. These rays will not have a weight for integration, but can serve as output for different viewing angles. Should be a dict with keys muz and mux as iterables.
temperature_threshold (float, optional) – Whether to ignore regions of the model above a particular temperature (e.g. 250 kK). This can dramatically speed up models which are primarily coronal, in addition to reducing memory pressure. Default: None, i.e. don’t
- iterate_se(*args, update_model_kwargs: dict | None = None, **kwargs)[source]
Iterate the context towards statistical equilibrium solution. Slightly modified variant of the core lw.iterate_ctx_se to add and handle the update_model callback.
- Parameters:
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).
update_model (callable(PromModel, bool, **kwargs), optional) – A function to use to update model parameters based on iteration (e.g. scale populations to maintain pressure), receives the model, and whether it should print based on the current iteration progress loop.
update_model_kwargs (dict, optional) – Extra kwargs to be passed to update_model.
- Returns:
it (int) – The number of iterations taken.
finalIterationUpdates (List[IterationUpdate], optional) – The final IterationUpdates computed, if requested by returnFinalConvergence.
promweaver.j_prom_bc module
- class promweaver.j_prom_bc.UniformJPromBc(projection: str, bc_provider: PromBcProvider, altitude_m: float, prom_upper_lower: str | None = None, Nrays: int = 50)[source]
Bases:
PromBc
Simple boundary condition that supplies \(J_\nu\) at each incoming ray and frequency for the boundary condition. This value is uniform over each ray. Putting the boundary condition into final_synthesis mode by requesting a single compute_rays from the model, or setting final_synthesis to true on the lw.Atmosphere used in the model samples the radiation directly along the ray path, as the converged populations have been determined.
- Parameters:
projection (str) – filament or prominence
bc_provider (PromBcProvider) – The provider to be used for computing the necessary intensity values.
altitude_m (float) – The altitude of the prominence above the solar surface [m]
prom_upper_lower (str, optional) – Whether this is the ‘upper’ or ‘lower’ z-boundary for a prominence (not used in filament cases).
Nrays (int, optional) – The number of Gauss-Legendre rays to be used to compute J. Default: 50.
- compute_bc(atmos, spect)[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
promweaver.limb_darkening module
- promweaver.limb_darkening.outgoing_chromo_ray_mu(mu, altitude_m)[source]
Compute the outgoing mu (relative to the surface normal at the location of the ray hit). Assumes Sun is a sphere with Radius 695,700 km. Radius from: 2008ApJ…675L..53H
- Parameters:
mu (float) – The mu of the ray of interest for the prominence.
altitude_m (float) – The altitude of the prominence [m].
- Returns:
mu_out – The requested outgoing mu, or None if the ray does not hit the Sun.
- Return type:
float, optional
promweaver.pctr_prom module
- class promweaver.pctr_prom.PctrPromModel(projection, cen_temperature, tr_temperature, cen_pressure, tr_pressure, max_cmass, vturb, altitude, gamma, active_atoms: List[str], detailed_atoms: List[str] = None, detailed_pops: Dict[str, ndarray] | None = None, atomic_models=None, Nhalf_points=100, Ncmass_decades=6.0, Nrays=3, prd=False, vlos: float | None = None, vrad: float | None = None, Nthreads=1, ctx_kwargs=None, BcType: Type[PromBc] | None = None, bc_kwargs=None, bc_provider=None, add_vertical_ray: bool = False, add_extra_rays: Dict[str, List[float] | Tuple[float] | ndarray] = None)[source]
Bases:
PromModel
Class for prominence with prominence-corona transition region (PCTR) simulations. Uses the stratifications first described by Anzer & Heinzel (1999), and later adopted by Labrosse et al. in PROM.
As a function of column mass \(m\), with maximum value along the line of sight \(M\), the pressure and temperature stratifications in the first half of the slab \(m \in [0, M/2]\) are given by
\[\begin{split}p(m) = 4 * (p_{cen} - p_{tr}) \frac{m}{M} \left(1 - \frac{m}{M}\right) + p_{tr}\\ T(m) = T_{cen} + (T_{tr} - T_{cen}) * \left(1 - 4 \frac{m}{M} \left(1 - \frac{m}{M}\right)\right)^\gamma\end{split}\]This is then reflected in the second half of the slab.
- Parameters:
projection (str) – Whether the object is to be treated as a “filament” or a “prominence”.
cen_temperature (float) – The central temperature of the prominence [K].
tr_temperature (float) – The transition-region temperature of the prominence [K].
cen_pressure (float) – The central pressure of the prominence [Pa].
tr_pressure (float) – The transition-region pressure of the prominence [Pa].
max_cmass (float) – The maximum column mass of the prominence model [kg/m2].
vturb (float) – The microturbulent velocity inside the prominence [m/s].
altitude (float) – The altitude of the prominence above the solar surface [m].
gamma (float) – The slope of the temperature gradient in the PCTR (must be >=2.0). Higher values are steeper.
active_atoms (list of str) – The element names to make “active” i.e. consider in non-LTE.
detailed_atoms (list of str, optional) – The element names to make “detailed static” i.e. consider in detail radiatively. If ‘H’ is set detailed then the model will be run with a fixed electron density (set to balance pressure) and iterative charge conservation disabled. Provide a compatible bc_provider as these will not be taken account in the default one computed on the fly.
detailed_pops (dict, optional) – A mapping from names in detailed_atoms to an array of their population distribution of shape suitable for the model.
atomic_models (list of lw.AtomicModels, optional) – The atomic models to use, a default set will be chosen if none are specified.
Ncmass_decades (float, optional) – The number of column mass decades to span in each half of the slab, i.e. the value of \(log10(M/2) - log10(m_0)\). Default: 6.
Nhalf_points (int, optional) – The number of points in half of the slab. Default: 100. Total slab size will be 2N-1 for symmetry reasons.
Nrays (int, optional) – The number of Gauss-Legendre angular quadrature rays to use in the model. Default: 3. This number will need to be set higher (e.g. 10) if using the ConePromBc.
Nthreads (int, optional) – The number of CPU threads to use when solving the radiative transfer equations. Default: 1.
prd (bool, optional) – Whether to consider the effects of partial frequency redistribution. Default: False.
vlos (float, optional) – The z-projected velocity to apply to the prominence model. Default: None, i.e. 0.
vrad (float, optional) – The Doppler-dimming radial velocity to apply. Note that for filaments this is the same as vlos and that should be used instead. Not fully supported in boundary conditions yet, (i.e. you should interpolate the wavelength grid first). Default: None, i.e. 0.
ctx_kwargs (dict, optional) – Extra kwargs to be passed when constructing the Context.
BcType (Constructor for a type of PromBc, optional) – The base type to be used for constructing the boundary conditions. Default: UniformJPromBc.
bc_kwargs (dict, optional) – Extra kwargs to be passed to the construction of the boundary conditions.
bc_provider (PromBcProvider, optional) – The provider to use for computing the radiation in the boundary conditions. Default: DynamicContextPromBcProvider using an lw.Context configured to match the current model. Note that the default is note very performant, but is convenient for experimenting. When running a grid of models, consider creating a TabulatedPromBcProvider using compute_falc_bc_ctx and tabulate_bc, since the default performs quite a few extra RT calculations.
add_vertical_ray (bool, optional) – Whether to add a non-weighted vertical ray to the model. This doesn’t participate in energy balance but allows for directly extracting mu=1, without needing to do a boundary condition. Default: False. Equivalent to add_extra_mus={“muz”: [1.0], “mux”: [0.0]}.
add_extra_rays (dict, optional) – Extra rays to add to the quadrature along which the solution should be sampled. These rays will not have a weight for integration, but can serve as output for different viewing angles. Should be a dict with keys muz and mux as iterables.
- iterate_se(*args, update_model_kwargs: dict | None = None, **kwargs)[source]
Iterate the context towards statistical equilibrium solution. Slightly modified variant of the core lw.iterate_ctx_se to add and handle the update_model callback.
- Parameters:
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).
update_model (callable(PromModel, bool, **kwargs), optional) – A function to use to update model parameters based on iteration (e.g. scale populations to maintain pressure), receives the model, and whether it should print based on the current iteration progress loop.
update_model_kwargs (dict, optional) – Extra kwargs to be passed to update_model.
- Returns:
it (int) – The number of iterations taken.
finalIterationUpdates (List[IterationUpdate], optional) – The final IterationUpdates computed, if requested by returnFinalConvergence.