Implemented Astrophysical Environments
This page summarizes the available astrophysical environments in which the photon-ALP propagation
can be calculated.
The environments can be added to the propagation by invoking the add_propagation()
function, see the Core Modules page for more information, as well as the Tutorials.
The astrophysical environments combine the magnetic field models (see Magnetic Field Models) and models for the electron density \(n_\mathrm{el}(r)\) (see Electron Density Models) and are summarized in the table below.
A minimal example for initializing mixing in the intergalactic magnetic field (IGMF) would like this:
from gammaALPs.core import Source, ALP, ModuleList
src = Source(z=0.536) # initialize a source with redshift z=0.536
alp = ALP(m=1., g=1.)
ml = ModuleList(alp, src)
ml.add_propagaion('IGMF')
px, py, pa = ml.run()
The other supported environments are listed in the table below.
Environment Name |
Environment Class |
Intended Use |
Used \(B\)-field model |
Used \(n_\mathrm{el}\) model |
---|---|---|---|---|
IGMF |
Mixing in the intergalactic magnetic field with a cell-like structure |
constant (evolves with redshift) |
||
ICMCell |
Mixing in a galaxy cluster magnetic field with a cell-like structure that decreases with growing distance from cluster center following \(n_\mathrm{el}(r)\) |
|||
ICMGaussTurb |
Mixing in a galaxy cluster magnetic field with Gaussian turbulence that decreases with growing distance from cluster center following \(n_\mathrm{el}(r)\) |
|||
ICMStructured |
Mixing in an intra-cluster cavity field. |
:py:class:`~gammaALPs.bfields.struc.structured_field |
||
Jet |
Mixing in the toroidal magnet field of an AGN jet, where the B field and \(n_\mathrm{el}(r)\) decrease as a power law with increasing distance from the central supermassive black hole |
|||
JetHelicalTangled |
Mixing in helical and tangled field of an AGN jet, where the B field and \(n_\mathrm{el}(r)\) decrease as a power law with increasing distance from the central supermassive black hole. |
|||
GMF |
Mixing in the magnetic field of the Milky Way. Currently, no model for \(n_\mathrm{el}\) is implemented |
constant |
||
File |
Environment initalized from B field and \(n_\mathrm{el}\) from a file |
— |
— |
|
Array |
Environment initalized from B field and \(n_\mathrm{el}\) from numpy arrays |
— |
— |
|
EBL |
|
No mixing, only absorption on the EBL |
— |
— |
Reference / API
- class gammaALPs.base.environs.MixFromArray(alp, Btrans, psi, nel, dL, **kwargs)[source]
Bases:
GammaALPTransfer
- __init__(alp, Btrans, psi, nel, dL, **kwargs)[source]
Initialize mixing in environment given by numpy arrays
- Parameters:
alp (
ALP
) –ALP
object with ALP parametersBtrans (array-like) – n-dim or m x n-dim array with absolute value of transversal B field, in muG if m x n-dimensional, m realizations are assumed
psi (array-like) – n-dim or m x n-dim array with angles between transversal direction and one polarization, if m x n-dimensional, m realizations are assumed
nel (array-like) – n-dim or m x n-dim array electron density, in cm^-3, if m x n-dimensional, m realizations are assumed
dL (array-like) – n-dim array with bin widths along line of sight in kpc, if m x n-dimensional, m realizations are assumed
EGeV (array-like) – Gamma-ray energies in GeV
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
chi (array-like) – n x m -dim numpy array with photon dispersion rate at energy E and distance L
- class gammaALPs.base.environs.MixFromFile(alp, filename, **kwargs)[source]
Bases:
GammaALPTransfer
- __init__(alp, filename, **kwargs)[source]
Initialize mixing in environment given by a data file. Data file has to have 5 columns: 1: distance along l.o.s. (z-axis) in kpc (bin bounds) 2: electron density in cm^-3 at bin bounds 3: Temperature in K 4: Bx component in muG at bin bounds 5: By component in muG at bin bounds 6: Bz component in muG at bin bounds
- Parameters:
alp (
ALP
) –ALP
object with ALP parametersfilename (str) – full path to file with electron density and B field
EGeV (array-like) – Gamma-ray energies in GeV
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
- class gammaALPs.base.environs.MixGMF(alp, source, **kwargs)[source]
Bases:
GammaALPTransfer
- property Bfield_model
- Bgmf_calc(l=0.0, b=0.0)[source]
compute GMF at (s,l,b) position where origin at self.d along x-axis in GC coordinates is assumed
- Parameters:
- Returns:
Btrans, Psin – (3,N)-dim arrays containing GMF for all domains in galactocentric cylindrical coordinates (rho, phi, z) and N-dim array with angles between propagation direction and line of sight.
- Return type:
tuple with
ndarray
- __init__(alp, source, **kwargs)[source]
Initialize mixing in the coherent component of the Galactic magnetic field
- Parameters:
alp (
ALP
) –ALP
object with ALP parameterssource (
Source
) –Source
object with source parametersEGeV (array-like) – Gamma-ray energies in GeV
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
int_steps (int (default = 100)) – Number of integration steps
rbounds (array-like) – bin bounds for steps along line of sight in kpc, default: lin range between end of Galactic Bfield and 0. with int_step steps
chi (~scipy.interpolate.RectBivariateSpline) – Spline function in (E [GeV], r [kpc]) giving values of photon-photon dispersion chi along the line of sight
parameters (Electron density) –
ra (float) – R.A. of the source (J2000)
dec (float) – Declination of the source (J2000)
galactic (float) – Distance of source to sun in kpc. If -1, source is extragalactic
parameters –
rho_max (float) – maximal rho of GMF in kpc default: 20 kpc
zmax (float) – maximal z of GMF in kpc default: 50 kpc
model (str) – (default = jansson) GMF model that is used. Currently the model by Jansson & Farrar (2012) (also with updates from Planck measurements) and Pshirkov et al. (2011) are implemented. Usage: model=[jansson12, jansson12b, jansson12c, pshirkov]
model_sym (str) – (default = ASS) Only applies if pshirkov model is chosen: you can choose between the axis- and bisymmetric version by setting model_sym to ASS or BSS.
parameters –
n0 (float) – Electron density in cm^-3 (default: 10). NE2001 code implementation still missing.
- property galactic
- property nel_model
- property r
- property rbounds
- class gammaALPs.base.environs.MixICMCell(alp, **kwargs)[source]
Bases:
GammaALPTransfer
- property Bfield_model
- __init__(alp, **kwargs)[source]
Initialize mixing in the intracluster magnetic field, assuming that it follows a domain-like structure.
- Parameters:
alp (
ALP
) –ALP
object with ALP parametersEGeV (array-like) – Gamma-ray energies in GeV
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
rbounds (array-like) – bin bounds for steps along line of sight in kpc, default: linear range between 0. and r_abell with L0 as step size
B0 (float) – ICM at r = 0 in muG
L0 (float) – Coherence length in kpc
nsim (int) – number of B field realizations
r_abell (float) – Abell radius of cluster (radius until which oscillation probability is computed)
kwargs (ICM) –
n0 (float) – electron density in cm**-3 (default 1e-3)
r_core (float) – core radius in kpc (default 10.)
beta (float) – exponent of density profile (default: 1.)
eta (float) – exponent for scaling of B field with electron density (default = 2./3.)
n2 (float) – if > 0., use profile with this second density component
r_core2 (float) – if > 0., use profile with this second r_core value
beta2 (float) – if > 0., use profile with this second beta value as for NGC1275
chi (~scipy.interpolate.RectBivariateSpline) – Spline function in (E [GeV], r [kpc]) giving values of photon-photon dispersion chi along the line of sight through the cluster
- property nel_model
- property r
- property rbounds
- class gammaALPs.base.environs.MixICMGaussTurb(alp, **kwargs)[source]
Bases:
GammaALPTransfer
- property Bfield_model
- __init__(alp, **kwargs)[source]
Initialize mixing in the intracluster magnetic field, assuming that it follows a Gaussian turbulence
- Parameters:
alp (
ALP
) –ALP
object with ALP parametersEGeV (array-like) – Gamma-ray energies in GeV
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
rbounds (larray-like) – bin bounds for steps along line of sight in kpc. Default: linear range between 0. and r_abell with 1/kH (min turbulence scale) as step size
thinning (int) – if > 1, thin out array of r
kwargs (ICM) –
B0 (float) – ICM at r = 0 in muG
kH (float) – upper wave number cutoff, should be at at least > 1. / osc. wavelength (default = 1 / (1 kpc))
kL (float) – lower wave number cutoff, should be of same size as the system (default = 1 / (100 kpc))
q (float) – power-law turbulence spectrum (default: q = 11/3 is Kolmogorov type spectrum)
kMin (float) – minimum wave number in 1. / kpc, default 1e-3 * kL (the k interval runs from kMin to kH)
dkType (string) – either linear, log, or random. Determine the spacing of the dk intervals
dkSteps (int) – number of dkSteps. For log spacing, number of steps per decade / number of decades ~ 10 should be chosen.
nsim (int) – number of B field realizations
kwargs –
n0 (float) – electron density in cm**-3 (default 1e-3)
r_core (float) – core radius in kpc (default 10.)
beta (float) – exponent of density profile (default: 1.)
eta (float) – exponent for scaling of B field with electron density (default = 2./3.)
n2 (float) – if > 0., use profile with this second density component
r_core2 (float) – if > 0., use profile with this second r_core value
beta2 (float) – if > 0., use profile with this second beta value as for NGC1275
chi (~scipy.interpolate.RectBivariateSpline) – Spline function in (E [GeV], r [kpc]) giving values of photon-photon dispersion chi along the line of sight through the cluster
- property nel_model
- property r
- property rbounds
- class gammaALPs.base.environs.MixICMStructured(alp, **kwargs)[source]
Bases:
GammaALPTransfer
- property Bfield_model
- __init__(alp, **kwargs)[source]
Initialize mixing in the intracluster magnetic field, assuming that it follows a structured field, see https://arxiv.org/pdf/1008.5353.pdf.
- Parameters:
alp (
ALP
) –ALP
object with ALP parametersEGeV (array-like) – Gamma-ray energies in GeV
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
B0 (float) – ICM at r = 0 in muG (default 1).
R (float) – Radius of cavity, B field vanishes at r=R, in kpc (default 100).
theta (float) – Angle of B field symmetry axis w.r.t. line of sight in degrees (default 0).
radians (bool) – If True, theta is considered to be expressed in radians (default False).
pa (float) – Position angle of symmetry axis in galactic coordinates (default 0).
pa_rad (bool) – True if pa given in radians (default False).
cell_num (int) – Number of cells B field is divided into for propagation of polarization density matrix (default 1000).
ICM kwargs:
- n0: float
electron density in cm**-3 (default 1e-3)
- r_core: float
core radius in kpc (default 10.)
- beta: float
exponent of density profile (default: 1.)
- n2: float
if > 0., use profile with this second density component
- r_core2: float
if > 0., use profile with this second r_core value
- beta2: float
if > 0., use profile with this second beta value as for NGC1275
- chi: ~scipy.interpolate.RectBivariateSpline
Spline function in (E [GeV], r [kpc]) giving values of photon-photon dispersion chi along the line of sight through the cluster.
- property nel_model
- property r
- class gammaALPs.base.environs.MixIGMFCell(alp, source, **kwargs)[source]
Bases:
GammaALPTransfer
- property Bfield_model
- __init__(alp, source, **kwargs)[source]
Initialize mixing in the intergalactic magnetic field (IGMF).
- Parameters:
alp (
ALP
) –ALP
object with ALP parameterssource (
Source
) –Source
object with source parametersEGeV (array-like) – Gamma-ray energies in GeV
dL (array-like) – Domain lengths. If not given, they will be automatically generated.
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
B0 (float) – IGMF at z = 0 in muG
L0 (float) – Coherence length at z = 0 in kpc
n0 (float) – electron density of intergalactic medium at z=0 in cm^-3
chi (~scipy.interpolate.RectBivariateSpline) – Spline function in (E [GeV], z) giving values of photon-photon dispersion chi at redshift z.
eblmodel (string) – name of the used EBL model (default: Dominguez et al. 2011 Model)
cosmo (~astropy.cosmology.core.FlatLambdaCDM) – chosen cosmology, default is H0 = 70, Om0 = 0.3
nsim (int) – number of B field realizations
- property nel_model
- property t
- class gammaALPs.base.environs.MixJet(alp, source, **kwargs)[source]
Bases:
GammaALPTransfer
- property Bfield_model
- property Rjet
- __init__(alp, source, **kwargs)[source]
Initialize mixing in the magnetic field of the jet, assumed here to be coherent
- Parameters:
alp (
ALP
) –ALP
object with ALP parameterssource (
Source
) –Source
object with source parametersEGeV (array-like) – Gamma-ray energies in GeV
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
rgam (float) – distance of gamma-ray emitting region to BH in pc (default: 0.1)
sens (float) – sens > 0 and sens < 1., sets the number of domains, for the B field in the n-th domain, it will have changed by B_n = sens * B_{n-1}
rbounds (list or ~numpy.ndarray) – bin bounds for steps along line of sight in pc, default: log range between rgam and Rjet with step size chosen such that B field changes by sens parameter in each step
kwargs (Jet) –
B0 (float) – Jet field at r = R0 in G (default: 0.1)
r0 (float) – distance from BH where B = B0 in pc (default: 0.1)
alpha (float) – exponent of toroidal mangetic field (default: -1.)
psi (float) – Angle between one photon polarization state and B field. Assumed constant over entire jet. (default: pi / 4)
helical (bool) – if True, use helical magnetic-field model from Clausen-Brown et al. (2011). In this case, the psi kwarg is treated is as the phi angle of the photon trajectory in the cylindrical jet coordinate system (default: True)
kwargs –
n0 (float) – electron density at R0 in cm**-3 (default 1e3)
beta (float) – exponent of electron density (default = -2.)
equipartition (bool) – if true, assume equipartition between electrons and the B field. This will overwrite beta = 2 * alpha and set n0 given the minimum electron lorentz factor set with gamma_min
gamma_min (float) – minimum lorentz factor of emitting electrons, only used if equipartition = True
gamma_max (float) – maximum lorentz factor of emitting electrons, only used if equipartition = True by default assumed to be gamma_min * 1e4
kwargs –
Rjet (float) – maximum jet length in pc (default: 1000.)
theta_obs (float) – Angle between l.o.s. and jet axis in degrees (default: 3.)
bulk_lorentz (float) – bulk lorentz factor of gamma-ray emitting plasma (default: 10.)
theta_jet (float) – Jet opening angle in degrees. If not given, assumed to be 1./bulk_lorentz
chi (~scipy.interpolate.RectBivariateSpline) – Spline function in (E [GeV], r [pc]) giving values of photon-photon dispersion chi down the jet
- property nel_model
- property r
- property rbounds
- class gammaALPs.base.environs.MixJetHelicalTangled(alp, source, **kwargs)[source]
Bases:
GammaALPTransfer
- property Bfield_model
- property Rjet
- __init__(alp, source, **kwargs)[source]
Initialize mixing in the magnetic field of the jet
- Parameters:
alp (
ALP
) –ALP
object with ALP parameterssource (
Source
) –Source
object with source parametersEGeV (array-like) – Gamma-ray energies in GeV
restore (str or None) – if str, identifier for files to restore environment. If None, initialize mixing with new B field
restore_path (str) – full path to environment files
ndom (Number of domains in jet model. (default 400)) –
kwargs (Jet) –
ft (float) – fraction of magnetic field energy density in tangled field
r_T (float) – radius at which helical field becomes toroidal in pc
Bt_exp (float) – exponent of the transverse component of the helical field at r<=r_T. i.e. sin(pitch angle) ~ r^Bt_exp while r<r_T and pitch angle = pi/2 at r=r_T
B0 (float) – Bfield strength in G
r0 (float) – radius where B field is equal to b0 in pc
rvhe (float) – distance of gamma-ray emission region from BH in pc
rjet (float) – jet length in pc
rem (float) – distance of emission region from BH if different from large-scale jet transition region in pc
alpha (float) – power-law index of electron energy distribution function
l_tcor (float) – tangled field coherence average length in pc
jwf (float) – jet width factor used when calculating l_tcor = jwf*jetwidth
jwf_dist (string) – type of distribution for jet width factors (jwf) when calculating l_tcor with jwf*jetwidth
tseed (float) – seed for random tangled domains
kwargs –
n0 (float) – electron density at R0 in cm**-3 (default 1e3)
beta (float) – exponent of electron density (default = -2.)
kwargs –
gmin (float) – jet lorenz factor at rjet
chi (~scipy.interpolate.RectBivariateSpline) – Spline function in (E [GeV], r [pc]) giving values of photon-photon dispersion chi down the jet
Gamma (<class 'scipy.interpolate.fitpack2.RectBivariateSpline'>) – Spline function in (E [GeV], r [pc]) giving values of absorption rate Gamma in kpc^-1 down the jet
- property gammas
- static jet_gammas_scaled_gg(rs, rvhe, rjet, gmin, gmax)[source]
Function to get jet lorentz factors. The shape of the gammas vs. r from PC Jet model, scaled to r0, gmin, gmax and rjet. Jet accelerates in the parabolic base (up to rvhe), then logarithmically decelerates in the conical jet.
- property nel_model
- property r
- property rbounds