agnpy package

Modules

agnpy.spectra module

class agnpy.spectra.BrokenPowerLaw(k=<Quantity 1.e-13 1 / cm3>, p1=2.0, p2=3.0, gamma_b=1000.0, gamma_min=10, gamma_max=10000000.0, mass=<<class 'astropy.constants.codata2018.CODATA2018'> name='Electron mass' value=9.1093837015e-31 uncertainty=2.8e-40 unit='kg' reference='CODATA 2018'>, integrator=<function trapz>)[source]

Bases: ParticleDistribution

Class describing a broken power-law particle distribution. When called, the particle density \(n(\gamma)\) in \(\mathrm{cm}^{-3}\) is returned.

\[n(\gamma') = k \left[ \left(\frac{\gamma'}{\gamma'_b}\right)^{-p_1} \, H(\gamma'; \gamma'_{\rm min}, \gamma'_b) + \left(\frac{\gamma'}{\gamma'_b}\right)^{-p_2} \, H(\gamma'; \gamma'_{b}, \gamma'_{\rm max}) \right]\]
Parameters:
  • k (Quantity) – spectral normalisation

  • p1 (float) – spectral index before the break (positive by definition)

  • p2 (float) – spectral index after the break (positive by definition)

  • gamma_b (float) – Lorentz factor at which the change in spectral index is occurring

  • gamma_min (float) – minimum Lorentz factor of the particle distribution

  • gamma_max (float) – maximum Lorentz factor of the particle distribution

  • mass (Quantity) – particle mass, default is the electron mass

  • integrator (func) – function to be used for integration, default is trapz

SSA_integrand(gamma)[source]
static evaluate(gamma, k, p1, p2, gamma_b, gamma_min, gamma_max)[source]
static evaluate_SSA_integrand(gamma, k, p1, p2, gamma_b, gamma_min, gamma_max)[source]

Analytical integrand for the synchrotron self-absorption: \(\gamma'^2 \frac{d}{d \gamma'} \left(\frac{n_e(\gamma)}{\gamma'^2}\right)\).

property parameters
class agnpy.spectra.ExpCutoffPowerLaw(k=<Quantity 1.e-13 1 / cm3>, p=2.1, gamma_c=1000.0, gamma_min=10, gamma_max=100000.0, mass=<<class 'astropy.constants.codata2018.CODATA2018'> name='Electron mass' value=9.1093837015e-31 uncertainty=2.8e-40 unit='kg' reference='CODATA 2018'>, integrator=<function trapz>)[source]

Bases: ParticleDistribution

Class describing a power-law with an exponetial cutoff particle distribution. When called, the particle density \(n_e(\gamma)\) in \(\mathrm{cm}^{-3}\) is returned.

\[n(\gamma'_c) = k \, \gamma'^{-p} exp(-\gamma'/\gamma_c) \, H(\gamma'; \gamma'{\rm min}, \gamma'{\rm max})\]
Parameters:
  • k (Quantity) – spectral normalisation

  • p (float) – spectral index, note it is positive by definition, will change sign in the function

  • gamma_c (float) – cutoff Lorentz factor of the particle distribution

  • gamma_min (float) – minimum Lorentz factor of the particle distribution

  • gamma_max (float) – maximum Lorentz factor of the particle distribution

  • mass (Quantity) – particle mass, default is the electron mass

  • integrator (func) – function to be used for integration, default is trapz

SSA_integrand(gamma)[source]
static evaluate(gamma, k, p, gamma_c, gamma_min, gamma_max)[source]
static evaluate_SSA_integrand(gamma, k, p, gamma_c, gamma_min, gamma_max)[source]

(analytical) integrand for the synchrotron self-absorption: \(\gamma'^2 \frac{d}{d \gamma'} \left(\frac{n(\gamma)}{\gamma'^2}\right)\)

property parameters
class agnpy.spectra.InterpolatedDistribution(gamma, n, norm=1, mass=<<class 'astropy.constants.codata2018.CODATA2018'> name='Electron mass' value=9.1093837015e-31 uncertainty=2.8e-40 unit='kg' reference='CODATA 2018'>, integrator=<function trapz>)[source]

Bases: ParticleDistribution

Class describing a particle distribution with an arbitrary shape. The spectrum is interpolated from an array of Lorentz factor and densities.

Parameters:
  • gamma (ndarray) – array of Lorentz factors where the density has been computed

  • n (Quantity) – array of densities to be interpolated

  • norm (float) – parameter to scale the density

  • mass (Quantity) – particle mass, default is the electron mass

  • integrator (func) – function to be used for integration, default is trapz

SSA_integrand(gamma)[source]

Integrand for the synchrotron self-absorption. It is

\[\gamma^2 \frac{d}{d \gamma} (\frac{n_e(\gamma)}{\gamma^2}) = ( \frac{dn_e(\gamma)}{d\gamma}+\frac{2n_e(\gamma)}{\gamma})\]

The derivative is:

\[\frac{dn_e(\gamma)}{d\gamma} = \frac{d 10^{f(u(\gamma))}}{d\gamma} = \frac{d10^{f(u)}}{du} \cdot \frac{du(\gamma)}{d\gamma}\]

where we have \(\frac{d 10^{f(u(\gamma))}}{d\gamma} = \frac{d10^{f(u)}}{du} \cdot \frac{du(\gamma)}{d\gamma}\), where \(u\) is the \(log_{10}(\gamma)\). This is equal to \(\frac{d 10^{f(u(\gamma))}}{d\gamma} = 10^{f(u)} \cdot \frac{df(u)}{du} \cdot \frac{1}{\gamma}\)

evaluate(gamma, norm, gamma_min, gamma_max)[source]
log10_interpolation()[source]

Returns the function interpolating in log10 the particle spectrum. TODO: make possible to pass arguments to CubicSpline.

class agnpy.spectra.LogParabola(k=<Quantity 1.e-13 1 / cm3>, p=2.0, q=0.1, gamma_0=1000.0, gamma_min=10, gamma_max=10000000.0, mass=<<class 'astropy.constants.codata2018.CODATA2018'> name='Electron mass' value=9.1093837015e-31 uncertainty=2.8e-40 unit='kg' reference='CODATA 2018'>, integrator=<function trapz>)[source]

Bases: ParticleDistribution

Class describing a log-parabolic particle distribution. When called, the particle density \(n(\gamma)\) in \(\mathrm{cm}^{-3}\) is returned.

\[n(\gamma') = k \, \left(\frac{\gamma'}{\gamma'_0}\right)^{-(p + q \log_{10}(\gamma' / \gamma'_0))}\]
Parameters:
  • k (Quantity) – spectral normalisation

  • p (float) – spectral index, note it is positive by definition, will change sign in the function

  • q (float) – spectral curvature, note it is positive by definition, will change sign in the function

  • gamma_0 (float) – reference Lorentz factor

  • gamma_min (float) – minimum Lorentz factor of the particle distribution

  • gamma_max (float) – maximum Lorentz factor of the particle distribution

  • mass (Quantity) – particle mass, default is the electron mass

  • integrator (func) – function to be used for integration, default is trapz

SSA_integrand(gamma)[source]
static evaluate(gamma, k, p, q, gamma_0, gamma_min, gamma_max)[source]
static evaluate_SSA_integrand(gamma, k, p, q, gamma_0, gamma_min, gamma_max)[source]

Analytical integrand for the synchrotron self-absorption: \(\gamma'^2 \frac{d}{d \gamma'} \left(\frac{n_e(\gamma)}{\gamma'^2}\right)\).

property parameters
class agnpy.spectra.ParticleDistribution(mass=<<class 'astropy.constants.codata2018.CODATA2018'> name='Electron mass' value=9.1093837015e-31 uncertainty=2.8e-40 unit='kg' reference='CODATA 2018'>, integrator=<function trapz>)[source]

Bases: object

Base class grouping common functionalities to be used by all particles distributions.

Parameters:
  • mass (~astropy.units.Quantity) – particle mass, default is the electron mass

  • integrator (function) – function to be used to integrate the particle distribution

classmethod from_density_at_gamma_1(n_gamma_1, mass, **kwargs)[source]

Set the normalisation of the particle distribution, \(k [{\rm cm}^{-3}]\), such that norm = \(n(\gamma=1)\).

Parameters:
  • n_gamma_1 (Quantity) – value \(n(\gamma)\) should have at \(\gamma=1\), in cm-3

  • mass (Quantity) – particle mass

classmethod from_total_density(n_tot, mass, **kwargs)[source]

Set the normalisation of the particle distribution, \(k [{\rm cm}^{-3}]\), from the total particle density \(n_{\rm tot} [{\rm cm}^{-3}]\).

Parameters:
  • n_tot (Quantity) – total particle density (integral of \(n(\gamma)\)), in cm-3

  • mass (Quantity) – particle mass

classmethod from_total_energy(W, V, mass, **kwargs)[source]

Set the normalisation of the particle distribution, \(k [{\rm cm}^{-3}]\), based on the total energy in particles \(W = m c^2 \, \int {\rm d}\gamma \, \gamma \, n(\gamma)\).

Parameters:
  • W (Quantity) – total energy in particles, in erg

  • V (Quantity) – volume of the emission region, in cm^3

  • mass (~astropy.units.Quantity) – particle mass

classmethod from_total_energy_density(u_tot, mass, **kwargs)[source]

Set the normalisation of the particle distribution, \(k [{\rm cm}^{-3}]\), from the total energy density \(u_{\rm tot} [{\rm erg}{\rm cm}^{-3}]\), Eq. 6.64 in [DermerMenon2009].

Parameters:
  • u_tot (Quantity) – total energy density (integral of \(\gamma\,n(\gamma)\)), in erg cm-3

  • mass (Quantity) – particle mass

static integral(self, gamma_low, gamma_up, gamma_power=0, integrator=<function trapz>, **kwargs)[source]

Integral of the particle distribution over the range gamma_low, gamma_up for a general set of parameters.

Parameters:
  • gamma_low (float) – lower integration limit

  • gamma_up (float) – higher integration limit

  • gamma_power (int) – power of gamma to raise the particle distribution before integration

  • integrator (func) – function to be used for integration, default is trapz

  • kwargs (dict) – parameters of the particle distribution

integrate(gamma_low, gamma_up, gamma_power=0)[source]

Integral of this particular particle distribution over the range gamma_low, gamma_up.

Parameters:
  • gamma_low (float) – lower integration limit

  • gamma_up (float) – higher integration limit

plot(gamma=None, gamma_power=0, ax=None, **kwargs)[source]

Plot the particle energy distribution.

Parameters:
  • gamma (ndarray) – array of Lorentz factors over which to plot the SED

  • gamma_power (float) – power of gamma to raise the electron distribution

  • ax (Axes, optional) – Axis

class agnpy.spectra.PowerLaw(k=<Quantity 1.e-13 1 / cm3>, p=2.1, gamma_min=10, gamma_max=100000.0, mass=<<class 'astropy.constants.codata2018.CODATA2018'> name='Electron mass' value=9.1093837015e-31 uncertainty=2.8e-40 unit='kg' reference='CODATA 2018'>, integrator=<function trapz>)[source]

Bases: ParticleDistribution

Class describing a power-law particle distribution. When called, the particle density \(n(\gamma)\) in \(\mathrm{cm}^{-3}\) is returned.

\[n(\gamma') = k \, \gamma'^{-p} \, H(\gamma'; \gamma'_{\rm min}, \gamma'_{\rm max})\]
Parameters:
  • k (Quantity) – spectral normalisation

  • p (float) – spectral index, note it is positive by definition, will change sign in the function

  • gamma_min (float) – minimum Lorentz factor of the particle distribution

  • gamma_max (float) – maximum Lorentz factor of the particle distribution

  • mass (Quantity) – particle mass, default is the electron mass

  • integrator (func) – function to be used for integration, default is trapz

SSA_integrand(gamma)[source]
static evaluate(gamma, k, p, gamma_min, gamma_max)[source]
static evaluate_SSA_integrand(gamma, k, p, gamma_min, gamma_max)[source]

Analytical integrand for the synchrotron self-absorption: \(\gamma'^2 \frac{d}{d \gamma'} \left(\frac{n(\gamma)}{\gamma'^2}\right)\).

property parameters

agnpy.emission_regions module

class agnpy.emission_regions.Blob(R_b=<Quantity 1.e+16 cm>, z=0.069, delta_D=10, Gamma=10, B=<Quantity 1. G>, n_e=<agnpy.spectra.spectra.PowerLaw object>, n_p=None, xi=1.0, gamma_e_size=200, gamma_p_size=200)[source]

Bases: object

Simple spherical emission region.

Note: all these quantities are defined in the comoving frame so they are actually primed quantities, when referring the notation in [DermerMenon2009].

All the quantities returned from the base attributes (e.g. volume of the emission region and variability time scale are derived from the blob radius) will be returned as properties, such that their value will be updated if the parameter from which they are derived is updated.

Parameters:
  • R_b (Quantity) – radius of the blob

  • z (float) – redshift of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • Gamma (float) – Lorentz factor of the relativistic outflow

  • B (Quantity) – magnetic field in the blob (Gauss)

  • n_e (ParticleDistribution) – electron distribution contained in the blob

  • n_p (ParticleDistribution) – proton distribution contained in the blob

  • xi (float) – acceleration coefficient \(\xi\) for first-order Fermi acceleration \((\mathrm{d}E/\mathrm{d}t \propto v \approx c)\) used to compute limits on the maximum Lorentz factor via \((\mathrm{d}E/\mathrm{d}t)_{\mathrm{acc}} = \xi c E / R_L\)

  • gamma_e_size (int) – size of the array of electrons Lorentz factors

  • gamma_p_size (int) – size of the array of protons Lorentz factors

property B_cgs

Magnetic field decomposed in Gaussian-cgs units.

property Beta

Bulk Lorentz factor of the Blob.

N_e(gamma)[source]

Number of electrons as a function of the Lorentz factor, \(N_{\rm e}(\gamma') = V_{\rm b}\,n_{\rm e}(\gamma')\).

Parameters:

gamma (ndarray) – array of Lorentz factor over which to evaluate the number of electrons

property N_e_tot

Total number of electrons

\[N_{\rm e,\,tot} = \int^{\gamma'_{\rm max}}_{\gamma'_{\rm min}} {\rm d}\gamma' N_{\rm e}(\gamma').\]
N_p(gamma)[source]

Number of protons as a function of the Lorentz factor, \(N_{\rm p}(\gamma') = V_{\rm b}\,n_{\rm p}(\gamma')\).

Parameters:

gamma (ndarray) – array of Lorentz factor over which to evaluate the number of electrons

property N_p_tot

total number of electrons

\[N_{\rm p,\,tot} = \int^{\gamma'_{\rm max}}_{\gamma'_{\rm min}} {\rm d}\gamma' N_{\rm p}(\gamma').\]
property P_jet_B

Jet power in magnetic field

\[P_{\mathrm{jet},\,B} = 2 \pi R_{\rm b}^2 \beta \Gamma^2 c \frac{B^2}{8\pi}.\]
property P_jet_ke

Total jet power in kinetic energy of the particles

\[P_{{\rm jet},\,{\rm ke}} = 2 \pi R_{\rm b}^2 \beta \Gamma^2 c (u_{\rm e} + u_{\rm p}).\]
property U_B

Energy density of magnetic field

\[U_B = B^2 / (8 \pi)\]
property V_b

Volume of the blob.

property W_e

Total energy in electrons

\[W_{\rm e} = m_{\rm e} c^2\,\int^{\gamma'_{\rm max}}_{\gamma'_{\rm min}} {\rm d}\gamma' \gamma' N_{\rm e}(\gamma').\]
property W_p

Total energy in protons

\[W_{\rm p} = m_{\rm p} c^2\,\int^{\gamma'_{\rm max}}_{\gamma'_{\rm min}} {\rm d}\gamma' \gamma' N_{\rm p}(\gamma').\]
property d_L

Luminosity distance.

property gamma_e

Array of electrons Lorentz factors, to be used for integration in the reference frame comoving with the emission region.

property gamma_e_external_frame

Array of electrons Lorentz factors, to be used for integration in the reference frame external to the emission region.

property gamma_p

Array of protons Lorentz factors, to be used for integration in the reference frame comoving with the emission region.

property k_eq

ratio between totoal particle energy density and magnetic field energy density, Eq. 7.75 of [DermerMenon2009]

Type:

Equipartition parameter

property mu_s

Cosine of the viewing angle from the jet axis to the observer.

property n_e

Electron distribution.

property n_e_tot

Total density of electrons

\[n_{\rm e,\,tot} = \int^{\gamma'_{\rm max}}_{\gamma'_{\rm min}} {\rm d}\gamma' n_{\rm e}(\gamma').\]
property n_p

Proton distribution.

property n_p_tot

Total density of protons

\[n_{\rm p,\,tot} = \int^{\gamma'_{\rm max}}_{\gamma'_{\rm min}} {\rm d}\gamma' n_{\rm p}(\gamma').\]
set_delta_D(Gamma, theta_s)[source]

Set the Doppler factor by specifying the bulk Lorentz factor of the outflow and the viewing angle.

Parameters:
  • Gamma (float) – Lorentz factor of the relativistic outflow

  • theta_s (Quantity) – viewing angle of the jet

set_gamma_e(gamma_size, gamma_min=1, gamma_max=100000000.0)[source]

Set the array of Lorentz factors for the electrons.

set_gamma_p(gamma_size, gamma_min=1, gamma_max=100000000.0)[source]

Set the array of Lorentz factors for the protons.

property t_var

Variability time scale, defined as \(t_{\rm var} = \frac{(1 + z) R_{\rm b}}{c \delta_{\rm D}}\).

property theta_s

Viewing angle from the jet axis to the observer.

property u_e

Total energy density of electrons

\[u_{\rm e} = m_{\rm e} c^2\,\int^{\gamma'_{\rm max}}_{\gamma'_{\rm min}} {\rm d}\gamma' \gamma' n_{\rm e}(\gamma').\]
property u_p

Total energy density of protons

\[u_{\rm p} = m_{\rm p} c^2\,\int^{\gamma'_{\rm max}}_{\gamma'_{\rm min}} {\rm d}\gamma' \gamma' n_{\rm p}(\gamma').\]
property u_ph_synch

Energy density of the synchrotron photons energy losses are:

\[(\mathrm{d}E/\mathrm{d}t)_{\mathrm{synch}} = 4 / 3 \sigma_T c U_B \gamma^2\]

the radiation stays an average time of \((3/4) (R_b/c)\) (the factor of 3/4 cames from averaging over a sphere), so an e- with Lorentz factor \(\gamma\) produces:

\[0.75\,(\mathrm{d}E/\mathrm{d}t)_{\mathrm{synch}}\,(R_b/c)\,/\,V_b\]

of radiation. We need to integrate over the electron spectrum (and multiply back by V_b)

\[0.75\,\int n_e(\gamma) (\mathrm{d}E/\mathrm{d}t)_{\mathrm{synch}} R_b \mathrm{d}\gamma\]

so

\[u_{\mathrm{synch}} = \sigma_T U_B R_b \int n_e(\gamma) \, \gamma^2 \mathrm{d}\gamma\]

WARNING: this does not take into account SSA!

agnpy.targets module

class agnpy.targets.CMB(z)[source]

Bases: object

Cosmic Microwave Background radiation, approximated as an isotropic monochromatic target.

Parameters:

z (float) – redshift at which the CMB is considered

u(blob=None)[source]

integral energy density of the CMB

Parameters:

blob (Blob) – if provided, the energy density is computed in a reference frame comvoing with the blob

class agnpy.targets.PointSourceBehindJet(L_0, epsilon_0)[source]

Bases: object

Monochromatic point source behind the jet.

Parameters:
  • L_0 (Quantity) – luminosity of the source

  • epsilon_0 (float) – dimensionless monochromatic energy of the source

u(r, blob=None)[source]

integral energy density of the point source at distance r along the jet axis

Parameters:
  • r (Quantity) – array of distances along the jet axis

  • blob (Blob) – if provided, the energy density is computed in a reference frame comvoing with the blob

class agnpy.targets.RingDustTorus(L_disk, xi_dt, T_dt, R_dt=None)[source]

Bases: object

Dust Torus as infinitesimally thin annulus, from [Finke2016]. For the Compton scattering monochromatic emission at the peak energy of the Black Body spectrum is considered.

Parameters:
  • L_disk (Quantity) – Luminosity of the disk whose radiation is being reprocessed by the Torus

  • xi_dt (float) – fraction of the disk radiation reprocessed

  • T_dt (Quantity) – peak temperature of the black body emission of the Torus

  • R_dt (Quantity) – radius of the Torus, if not specified the saturation radius of Eq. 96 in [Finke2016] will be used

static evaluate_bb_norm_sed(nu, z, L_dt, T_dt, R_dt, d_L)[source]

evaluate the torus black-body SED such that its integral luminosity is equal to the torus luminosity (xi_dt * L_disk)

static evaluate_bb_sed(nu, z, T_dt, R_dt, d_L)[source]

evaluate the black body SED corresponding to the torus temperature

sed_flux(nu, z)[source]

evaluate the black-body SED for the Dust Torus

u(r, blob=None)[source]

Density of radiation produced by the Torus at the distance r along the jet axis. Integral over the solid angle of Eq. 85 in [Finke2016]

Parameters:
  • r (Quantity) – array of distances along the jet axis

  • blob (Blob) – if provided, the energy density is computed in a reference frame comvoing with the blob

class agnpy.targets.SSDisk(M_BH, L_disk, eta, R_in, R_out, R_g_units=False)[source]

Bases: object

[Shakura1973] accretion disk as modeled by [Dermer2002]

Parameters:
  • M_BH (Quantity) – Black Hole mass

  • L_disk (Quantity) – luminosity of the disk

  • eta (float) – accretion efficiency

  • R_in (Quantity / float) – inner disk radius

  • R_out (Quantity / float) – outer disk radius

  • R_g_units (bool) – whether or not input radiuses are specified in units of the gravitational radius

T(R)[source]

temperature of the disk at radius \(R\). Eq. 64 in [Dermer2009].

epsilon(R_tilde)[source]

dimensionless energy emitted at the disk radius \(\tilde{R}\)

epsilon_mu(mu, r_tilde)[source]
static evaluate_T(R, M_BH, m_dot, R_in)[source]

black body temperature (K) at the radial coordinate R

static evaluate_epsilon(L_disk, M_BH, eta, R_tilde)[source]

evaluate the dimensionless energy emitted at the radius R_tilde Eq. 65 [Dermer2009]

static evaluate_epsilon_mu(L_disk, M_BH, eta, mu, r_tilde)[source]

same as evaluate_epsilon() but considering the cosine of the subtended zenith mu and the height above the disk r instead of the radius R_tilde

static evaluate_mu_from_r_tilde(R_in_tilde, R_out_tilde, r_tilde, size=100)[source]

array of cosine angles, spanning from \(R_{\mathrm{in}}\) to \(R_{\mathrm{out}}\), viewed from a given height \(\tilde{r}\) above the disk, Eq. 72 and 73 in [Finke2016].

static evaluate_multi_T_bb_norm_sed(nu, z, L_disk, M_BH, m_dot, R_in, R_out, d_L, mu_s=1)[source]

Evaluate a normalised, multi-temperature black body SED as in evaluate_multi_T_bb_sed(), but the integral luminosity is set to be equal to L_disk.

static evaluate_multi_T_bb_sed(nu, z, M_BH, m_dot, R_in, R_out, d_L, mu_s=1)[source]

Evaluate a multi-temperature black body SED in the case of the SS Disk. The SED is calculated for an observer far away from the disk (\(d_L \gg R\)) with the following:

\[\nu F_{\nu} \approx \mu_s \, \nu \frac{2 \pi}{d_L^2} \int_{R_{\rm in}}^{R_{\rm out}}{\rm d}R \, R \, I_{\nu}(T(R)),\]

where \(I_{\nu}\) is Planck’s law, \(R\) the radial coordinate along the disk, and \(d_L\) the luminosity distance. \(\mu_s\) is the cosine of the angle between the disk axis and the observer’s line of sight.

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • M_BH (Quantity) – Black Hole mass

  • m_dot (float) – mass accretion rate

  • R_in (Quantity) – inner disk radius

  • R_out (Quantity) – outer disk radius

  • d_L (Quantity) – luminosity of the source

  • mu_s (float) – cosine of the angle between the observer line of sight and the disk axis

static evaluate_phi_disk_mu(mu, R_in_tilde, r_tilde)[source]

dependency of the radiant surface-energy flux from the disk radius, here obtained from the cosine of the zenith mu and the height above the disk r_tilde (in graviational radius units), Eq. 63 [Dermer2009]

phi_disk(R_tilde)[source]
phi_disk_mu(mu, r_tilde)[source]
sed_flux(nu, z, mu_s=1)[source]

evaluate the multi-temperature black body SED for this SS Disk, refer to evaluate_multi_T_bb_sed() and to evaluate_multi_T_bb_norm_sed()

u(r, blob=None)[source]

integral energy density of radiation produced by the Disk at the distance r along the jet axis. Integral over the solid angle of Eq. 69 in [Dermer2009].

Parameters:
  • r (Quantity) – array of distances along the jet axis

  • blob (Blob) – if provided, the energy density is computed in a reference frame comvoing with the blob

class agnpy.targets.SphericalShellBLR(L_disk, xi_line, line, R_line)[source]

Bases: object

Spherical Shell Broad Line Region, from [Finke2016]. Each line is emitted from an infinitesimally thin spherical shell.

Parameters:
  • L_disk (Quantity) – Luminosity of the disk whose radiation is being reprocessed by the BLR

  • xi_line (float) – fraction of the disk radiation reprocessed by the BLR

  • line (string) – type of line emitted

  • R_line (Quantity) – radius of the BLR spherical shell

print_lines_list()[source]

Print the list of the available spectral lines. The dictionary with the possible emission lines is taken from Table 5 in [Finke2016] and contains the value of the line wavelength and the ratio of its radius to the radius of the \(H_{\beta}\) shell, not used at the moment.

u(r, blob=None)[source]

Density of radiation produced by the BLR at the distance r along the jet axis. Integral over the solid angle of Eq. 80 in [Finke2016].

Parameters:
  • r (Quantity) – array of distances along the jet axis

  • blob (Blob) – if provided, the energy density is computed in a reference frame comvoing with the blob

agnpy.synchrotron module

class agnpy.synchrotron.ProtonSynchrotron(blob, integrator=<function trapz>)[source]

Bases: object

Class for synchrotron radiation computation

Parameters:
  • blob (Blob) – emitting region and proton distribution

  • ssa (bool) – whether or not to consider synchrotron self absorption (SSA). The absorption factor will be taken into account in com_sed_emissivity(), in order to be propagated to sed_luminosity() and sed_flux().

  • integrator (func) – function to be used for integration (default = np.trapz)

static evaluate_sed_flux(nu, z, d_L, delta_D, B, R_b, n_p, *args, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]))[source]

Evaluates the flux SED (\(\nu F_{\nu}\)) due to synchrotron radiation, for a general set of model parameters. As for electrons, we implement Eq. 21 in [Finke2008] with just a change in the mass value (we are using the proton mass now). For a reference on proton synchrotron and other hadronic processes see [Cerruti2015].

Note parameters after *args need to be passed with a keyword

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • d_L (Quantity) – luminosity distance of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • B (Quantity) – magnetic field in the blob

  • R_b (Quantity) – size of the emitting region (spherical blob assumed)

  • n_p (ProtonDistribution) – proton energy distribution

  • *args – parameters of the proton energy distribution (k_e, p, …)

  • ssa (#) –

  • self-absorption (# whether to consider or not the) –

  • false (default) –

  • integrator (func) – which function to use for integration, default numpy.trapz

  • gamma (ndarray) – array of Lorentz factor over which to integrate the proton distribution

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

sed_flux(nu)[source]

Evaluates the synchrotron flux SED for a Synchrotron object built from a Blob.

sed_luminosity(nu)[source]

Evaluates the synchrotron luminosity SED \(\nu L_{\nu} \, [\mathrm{erg}\,\mathrm{s}^{-1}]\) for a a Synchrotron object built from a blob.

sed_peak_flux(nu)[source]

provided a grid of frequencies nu, returns the peak flux of the SED

sed_peak_nu(nu)[source]

provided a grid of frequencies nu, returns the frequency at which the SED peaks

class agnpy.synchrotron.Synchrotron(blob, ssa=False, integrator=<function trapz>)[source]

Bases: object

Class for synchrotron radiation computation

Parameters:
  • blob (Blob) – emitting region and electron distribution

  • ssa (bool) – whether or not to consider synchrotron self absorption (SSA). The absorption factor will be taken into account in com_sed_emissivity(), in order to be propagated to sed_luminosity() and sed_flux().

  • integrator (func) – function to be used for integration (default = np.trapz)

static evaluate_sed_flux(nu, z, d_L, delta_D, B, R_b, n_e, *args, ssa=False, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]))[source]

Evaluates the flux SED (\(\nu F_{\nu}\)) due to synchrotron radiation, for a general set of model parameters. Eq. 21 in [Finke2008].

Note parameters after *args need to be passed with a keyword

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • d_L (Quantity) – luminosity distance of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • B (Quantity) – magnetic field in the blob

  • R_b (Quantity) – size of the emitting region (spherical blob assumed)

  • n_e (ElectronDistribution) – electron energy distribution

  • *args – parameters of the electron energy distribution (k_e, p, …)

  • ssa (bool) – whether to consider or not the self-absorption, default false

  • integrator (func) – which function to use for integration, default numpy.trapz

  • gamma (ndarray) – array of Lorentz factor over which to integrate the electron distribution

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

static evaluate_sed_flux_delta_approx(nu, z, d_L, delta_D, B, R_b, n_e, *args)[source]

Synchrotron flux SED using the delta approximation for the synchrotron radiation Eq. 7.70 [DermerMenon2009].

static evaluate_tau_ssa(nu, z, d_L, delta_D, B, R_b, n_e, *args, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]))[source]

Computes the syncrotron self-absorption opacity for a general set of model parameters, see evaluate_sed_flux() for parameters defintion. Eq. before 7.122 in [DermerMenon2009].

sed_flux(nu)[source]

Evaluates the synchrotron flux SED for a Synchrotron object built from a Blob.

sed_flux_delta_approx(nu)[source]

Evaluates the synchrotron flux SED using the delta approximation for a Synchrotron object built from a blob.

sed_luminosity(nu)[source]

Evaluates the synchrotron luminosity SED \(\nu L_{\nu} \, [\mathrm{erg}\,\mathrm{s}^{-1}]\) for a a Synchrotron object built from a blob.

sed_peak_flux(nu)[source]

provided a grid of frequencies nu, returns the peak flux of the SED

sed_peak_nu(nu)[source]

provided a grid of frequencies nu, returns the frequency at which the SED peaks

agnpy.compton module

class agnpy.compton.ExternalCompton(blob, target, r=None, integrator=<function trapz>)[source]

Bases: object

class for External Compton radiation computation

Parameters:
  • blob (Blob) – emission region and electron distribution hitting the photon target

  • target (targets) – class describing the target photon field

  • r (Quantity) – distance of the blob from the Black Hole (i.e. from the target photons)

  • integrator (func) – function to be used for integration (default = np.trapz)

static evaluate_sed_flux_blr(nu, z, d_L, delta_D, mu_s, R_b, L_disk, xi_line, epsilon_line, R_line, r, n_e, *args, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]), mu=array([-1., -0.97979798, -0.95959596, -0.93939394, -0.91919192, -0.8989899, -0.87878788, -0.85858586, -0.83838384, -0.81818182, -0.7979798, -0.77777778, -0.75757576, -0.73737374, -0.71717172, -0.6969697, -0.67676768, -0.65656566, -0.63636364, -0.61616162, -0.5959596, -0.57575758, -0.55555556, -0.53535354, -0.51515152, -0.49494949, -0.47474747, -0.45454545, -0.43434343, -0.41414141, -0.39393939, -0.37373737, -0.35353535, -0.33333333, -0.31313131, -0.29292929, -0.27272727, -0.25252525, -0.23232323, -0.21212121, -0.19191919, -0.17171717, -0.15151515, -0.13131313, -0.11111111, -0.09090909, -0.07070707, -0.05050505, -0.03030303, -0.01010101, 0.01010101, 0.03030303, 0.05050505, 0.07070707, 0.09090909, 0.11111111, 0.13131313, 0.15151515, 0.17171717, 0.19191919, 0.21212121, 0.23232323, 0.25252525, 0.27272727, 0.29292929, 0.31313131, 0.33333333, 0.35353535, 0.37373737, 0.39393939, 0.41414141, 0.43434343, 0.45454545, 0.47474747, 0.49494949, 0.51515152, 0.53535354, 0.55555556, 0.57575758, 0.5959596, 0.61616162, 0.63636364, 0.65656566, 0.67676768, 0.6969697, 0.71717172, 0.73737374, 0.75757576, 0.77777778, 0.7979798, 0.81818182, 0.83838384, 0.85858586, 0.87878788, 0.8989899, 0.91919192, 0.93939394, 0.95959596, 0.97979798, 1.        ]), phi=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the flux SED (\(\nu F_{\nu}\)) for External Compton on the photon field of a spherical shell BLR, for a general set of model parameters

Note parameters after *args need to be passed with a keyword

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • d_L (Quantity) – luminosity distance of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • L_disk (Quantity) – Luminosity of the disk whose radiation is being reprocessed by the BLR

  • xi_line (float) – fraction of the disk radiation reprocessed by the BLR

  • epsilon_line (string) – dimensionless energy of the emitted line

  • R_line (Quantity) – radius of the BLR spherical shell

  • r (Quantity) – distance between the Broad Line Region and the blob

  • n_e (ElectronDistribution) – electron energy distribution

  • *args – parameters of the electron energy distribution (k_e, p, …)

  • integrator (func) – which function to use for integration, default numpy.trapz

  • gamma (ndarray) – array of Lorentz factor over which to integrate the electron distribution

  • mu (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

  • phi (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

static evaluate_sed_flux_dt(nu, z, d_L, delta_D, mu_s, R_b, L_disk, xi_dt, epsilon_dt, R_dt, r, n_e, *args, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]), phi=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the flux SED (\(\nu F_{\nu}\)) for External Compton on the photon field of a ring dust torus, for a general set of model parameters

Note parameters after *args need to be passed with a keyword

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • d_L (Quantity) – luminosity distance of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • L_disk (Quantity) – Luminosity of the disk whose radiation is being reprocessed by the BLR

  • xi_dt (float) – fraction of the disk radiation reprocessed by the disk

  • epsilon_dt (string) – peak (dimensionless) energy of the black body radiated by the torus

  • R_dt (Quantity) – radius of the ting-like torus

  • r (Quantity) – distance between the Broad Line Region and the blob

  • n_e (ElectronDistribution) – electron energy distribution

  • *args – parameters of the electron energy distribution (k_e, p, …)

  • integrator (func) – which function to use for integration, default numpy.trapz

  • gamma (ndarray) – array of Lorentz factor over which to integrate the electron distribution

  • mu (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

  • phi (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

static evaluate_sed_flux_iso_mono(nu, z, d_L, delta_D, mu_s, R_b, epsilon_0, u_0, n_e, *args, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]), mu=array([-1., -0.97979798, -0.95959596, -0.93939394, -0.91919192, -0.8989899, -0.87878788, -0.85858586, -0.83838384, -0.81818182, -0.7979798, -0.77777778, -0.75757576, -0.73737374, -0.71717172, -0.6969697, -0.67676768, -0.65656566, -0.63636364, -0.61616162, -0.5959596, -0.57575758, -0.55555556, -0.53535354, -0.51515152, -0.49494949, -0.47474747, -0.45454545, -0.43434343, -0.41414141, -0.39393939, -0.37373737, -0.35353535, -0.33333333, -0.31313131, -0.29292929, -0.27272727, -0.25252525, -0.23232323, -0.21212121, -0.19191919, -0.17171717, -0.15151515, -0.13131313, -0.11111111, -0.09090909, -0.07070707, -0.05050505, -0.03030303, -0.01010101, 0.01010101, 0.03030303, 0.05050505, 0.07070707, 0.09090909, 0.11111111, 0.13131313, 0.15151515, 0.17171717, 0.19191919, 0.21212121, 0.23232323, 0.25252525, 0.27272727, 0.29292929, 0.31313131, 0.33333333, 0.35353535, 0.37373737, 0.39393939, 0.41414141, 0.43434343, 0.45454545, 0.47474747, 0.49494949, 0.51515152, 0.53535354, 0.55555556, 0.57575758, 0.5959596, 0.61616162, 0.63636364, 0.65656566, 0.67676768, 0.6969697, 0.71717172, 0.73737374, 0.75757576, 0.77777778, 0.7979798, 0.81818182, 0.83838384, 0.85858586, 0.87878788, 0.8989899, 0.91919192, 0.93939394, 0.95959596, 0.97979798, 1.        ]), phi=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the flux SED (\(\nu F_{\nu}\)) for external Compton on a monochromatic isotropic target photon field, for a general set of model parameters

Note parameters after *args need to be passed with a keyword

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • d_L (Quantity) – luminosity distance of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • R_b (Quantity) – size of the emitting region (spherical blob assumed)

  • epsilon_0 (float) – dimensionless energy (in electron rest mass energy units) of the target photon field

  • u_0 (Quantity) – energy density [erg cm-3] of the target photon field

  • n_e (ElectronDistribution) – electron energy distribution

  • *args – parameters of the electron energy distribution (k_e, p, …)

  • integrator (func) – which function to use for integration, default numpy.trapz

  • gamma (ndarray) – array of Lorentz factor over which to integrate the electron distribution

  • mu (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

  • phi (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

static evaluate_sed_flux_ps_behind_jet(nu, z, d_L, delta_D, mu_s, R_b, epsilon_0, L_0, r, n_e, *args, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]))[source]

Evaluates the flux SED (\(\nu F_{\nu}\)) for external Compton on a point source of photons behind the jet, for a general set of model parameters

Note parameters after *args need to be passed with a keyword

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • d_L (Quantity) – luminosity distance of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • R_b (Quantity) – size of the emitting region (spherical blob assumed)

  • epsilon_0 (float) – dimensionless energy (in electron rest mass energy units) of the target photon field

  • L_0 (Quantity) – luminosity [erg cm-3] of the point source behind the jet

  • r (Quantity) – distance between the point source and the blob

  • n_e (ElectronDistribution) – electron energy distribution

  • *args – parameters of the electron energy distribution (k_e, p, …)

  • integrator (func) – which function to use for integration, default numpy.trapz

  • gamma (ndarray) – array of Lorentz factor over which to integrate the electron distribution

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

static evaluate_sed_flux_ss_disk(nu, z, d_L, delta_D, mu_s, R_b, M_BH, L_disk, eta, R_in, R_out, r, n_e, *args, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]), mu_size=100, phi=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the flux SED (\(\nu F_{\nu}\)) for external Compton on the photon field of a Shakura Sunyaev disk, for a general set of model parameters

Note parameters after *args need to be passed with a keyword

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • d_L (Quantity) – luminosity distance of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • R_b (Quantity) – size of the emitting region (spherical blob assumed)

  • M_BH (Quantity) – Black Hole mass

  • L_disk (Quantity) – luminosity of the disk

  • eta (float) – accretion efficiency

  • R_in (Quantity) – inner disk radius

  • R_out (Quantity) – inner disk radius

  • r (Quantity) – distance between the disk and the blob

  • n_e (ElectronDistribution) – electron energy distribution

  • *args – parameters of the electron energy distribution (k_e, p, …)

  • integrator (func) – which function to use for integration, default numpy.trapz

  • gamma (ndarray) – array of Lorentz factor over which to integrate the electron distribution

  • mu_size (int) – size of the array of zenith angles to integrate over

  • phi (ndarray) – arrays of azimuth angles to integrate over

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

sed_flux(nu)[source]

SEDs for external Compton

sed_flux_blr(nu)[source]

evaluates the flux SED for External Compton on a spherical BLR

sed_flux_cmb(nu)[source]

evaluates the flux SED for External Compton on the CMB

sed_flux_dt(nu)[source]

evaluates the flux SED for External Compton on a ring dust torus

sed_flux_ps_behind_jet(nu)[source]

evaluates the flux SED for External Compton on a point source behind the jet

sed_flux_ss_disk(nu)[source]

evaluates the flux SED for External Compton on a [Shakura1973] disk

sed_luminosity(nu)[source]

Evaluates the external Compton luminosity SED \(\nu L_{\nu} \, [\mathrm{erg}\,\mathrm{s}^{-1}]\)

set_mu(mu_size=100)[source]
set_phi(phi_size=50)[source]
class agnpy.compton.SynchrotronSelfCompton(blob, ssa=False, integrator=<function trapz>)[source]

Bases: object

class for Synchrotron Self Compton radiation computation

Parameters:
  • blob (Blob) – emission region and electron distribution hitting the photon target

  • synchrotron (Synchrotron) – class describing the synchrotron photons target

  • integrator (func) – function to be used for integration (default = np.trapz)

static evaluate_sed_flux(nu, z, d_L, delta_D, B, R_b, n_e, *args, ssa=False, integrator=<function trapz>, gamma=array([1.00000000e+01, 1.09698580e+01, 1.20337784e+01, 1.32008840e+01, 1.44811823e+01, 1.58856513e+01, 1.74263339e+01, 1.91164408e+01, 2.09704640e+01, 2.30043012e+01, 2.52353917e+01, 2.76828663e+01, 3.03677112e+01, 3.33129479e+01, 3.65438307e+01, 4.00880633e+01, 4.39760361e+01, 4.82410870e+01, 5.29197874e+01, 5.80522552e+01, 6.36824994e+01, 6.98587975e+01, 7.66341087e+01, 8.40665289e+01, 9.22197882e+01, 1.01163798e+02, 1.10975250e+02, 1.21738273e+02, 1.33545156e+02, 1.46497140e+02, 1.60705282e+02, 1.76291412e+02, 1.93389175e+02, 2.12145178e+02, 2.32720248e+02, 2.55290807e+02, 2.80050389e+02, 3.07211300e+02, 3.37006433e+02, 3.69691271e+02, 4.05546074e+02, 4.44878283e+02, 4.88025158e+02, 5.35356668e+02, 5.87278661e+02, 6.44236351e+02, 7.06718127e+02, 7.75259749e+02, 8.50448934e+02, 9.32930403e+02, 1.02341140e+03, 1.12266777e+03, 1.23155060e+03, 1.35099352e+03, 1.48202071e+03, 1.62575567e+03, 1.78343088e+03, 1.95639834e+03, 2.14614120e+03, 2.35428641e+03, 2.58261876e+03, 2.83309610e+03, 3.10786619e+03, 3.40928507e+03, 3.73993730e+03, 4.10265811e+03, 4.50055768e+03, 4.93704785e+03, 5.41587138e+03, 5.94113398e+03, 6.51733960e+03, 7.14942899e+03, 7.84282206e+03, 8.60346442e+03, 9.43787828e+03, 1.03532184e+04, 1.13573336e+04, 1.24588336e+04, 1.36671636e+04, 1.49926843e+04, 1.64467618e+04, 1.80418641e+04, 1.97916687e+04, 2.17111795e+04, 2.38168555e+04, 2.61267523e+04, 2.86606762e+04, 3.14403547e+04, 3.44896226e+04, 3.78346262e+04, 4.15040476e+04, 4.55293507e+04, 4.99450512e+04, 5.47890118e+04, 6.01027678e+04, 6.59318827e+04, 7.23263390e+04, 7.93409667e+04, 8.70359136e+04, 9.54771611e+04, 1.04737090e+05, 1.14895100e+05, 1.26038293e+05, 1.38262217e+05, 1.51671689e+05, 1.66381689e+05, 1.82518349e+05, 2.00220037e+05, 2.19638537e+05, 2.40940356e+05, 2.64308149e+05, 2.89942285e+05, 3.18062569e+05, 3.48910121e+05, 3.82749448e+05, 4.19870708e+05, 4.60592204e+05, 5.05263107e+05, 5.54266452e+05, 6.08022426e+05, 6.66991966e+05, 7.31680714e+05, 8.02643352e+05, 8.80488358e+05, 9.65883224e+05, 1.05956018e+06, 1.16232247e+06, 1.27505124e+06, 1.39871310e+06, 1.53436841e+06, 1.68318035e+06, 1.84642494e+06, 2.02550194e+06, 2.22194686e+06, 2.43744415e+06, 2.67384162e+06, 2.93316628e+06, 3.21764175e+06, 3.52970730e+06, 3.87203878e+06, 4.24757155e+06, 4.65952567e+06, 5.11143348e+06, 5.60716994e+06, 6.15098579e+06, 6.74754405e+06, 7.40196000e+06, 8.11984499e+06, 8.90735464e+06, 9.77124154e+06, 1.07189132e+07, 1.17584955e+07, 1.28989026e+07, 1.41499130e+07, 1.55222536e+07, 1.70276917e+07, 1.86791360e+07, 2.04907469e+07, 2.24780583e+07, 2.46581108e+07, 2.70495973e+07, 2.96730241e+07, 3.25508860e+07, 3.57078596e+07, 3.91710149e+07, 4.29700470e+07, 4.71375313e+07, 5.17092024e+07, 5.67242607e+07, 6.22257084e+07, 6.82607183e+07, 7.48810386e+07, 8.21434358e+07, 9.01101825e+07, 9.88495905e+07, 1.08436597e+08, 1.18953407e+08, 1.30490198e+08, 1.43145894e+08, 1.57029012e+08, 1.72258597e+08, 1.88965234e+08, 2.07292178e+08, 2.27396575e+08, 2.49450814e+08, 2.73644000e+08, 3.00183581e+08, 3.29297126e+08, 3.61234270e+08, 3.96268864e+08, 4.34701316e+08, 4.76861170e+08, 5.23109931e+08, 5.73844165e+08, 6.29498899e+08, 6.90551352e+08, 7.57525026e+08, 8.30994195e+08, 9.11588830e+08, 1.00000000e+09]))[source]

Evaluates the flux SED (\(\nu F_{\nu}\)) for synchrotron self-Compton, for a general set of model parameters. Eq. 9 in [Finke2008].

Note parameters after *args need to be passed with a keyword

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • d_L (Quantity) – luminosity distance of the source

  • delta_D (float) – Doppler factor of the relativistic outflow

  • B (Quantity) – magnetic field in the blob

  • R_b (Quantity) – size of the emitting region (spherical blob assumed)

  • n_e (ElectronDistribution) – electron energy distribution

  • *args – parameters of the electron energy distribution (k_e, p, …)

  • ssa (bool) – whether to consider or not the self-absorption, default false

  • integrator (func) – which function to use for integration, default numpy.trapz

  • gamma (ndarray) – array of Lorentz factor over which to integrate the electron distribution

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

sed_flux(nu)[source]

Evaluates the SSC flux SED for a SynchrotronSelfComtpon object built from a Blob.

sed_luminosity(nu)[source]

Evaluates the SSC luminosity SED \(\nu L_{\nu} \, [\mathrm{erg}\,\mathrm{s}^{-1}]\) for a a SynchrotronSelfCompton object built from a blob.

sed_peak_flux(nu)[source]

provided a grid of frequencies nu, returns the peak flux of the SED

sed_peak_nu(nu)[source]

provided a grid of frequencies nu, returns the frequency at which the SED peaks

agnpy.absorption module

class agnpy.absorption.Absorption(target, r=None, z=0, mu_s=1)[source]

Bases: object

class to compute the absorption due to gamma-gamma pair production

Parameters:
  • blob (Blob) – emission region and electron distribution hitting the photon target

  • target (targets or class:~agnpy.emission_regions.Blob) – class describing the target photon field

  • r (Quantity) – distance of the blob from the Black Hole (i.e. from the target photons) the distance is irrelevant in the case of absorption

absorption(nu)[source]

This function returns the attenuation of the emission assuming that the optical depth tau is computed from the production place to the observer.

absorption_homogeneous(nu)[source]

This function returns the attenuation of the emission assuming that the emission is produced homogenously inside absorbing material. The calculations is only accurate for a slab of absorbing material with the total optical depth tau, but the same formula is often used also e.g. in the context of absorption of gamma-ray emission by synchrotron radiation in blobs See e.g. section 2.5.1. of Finke et al. 2008.

static evaluate_tau_blr(nu, z, mu_s, L_disk, xi_line, epsilon_line, R_line, r, l_size=50, mu=array([-1., -0.97979798, -0.95959596, -0.93939394, -0.91919192, -0.8989899, -0.87878788, -0.85858586, -0.83838384, -0.81818182, -0.7979798, -0.77777778, -0.75757576, -0.73737374, -0.71717172, -0.6969697, -0.67676768, -0.65656566, -0.63636364, -0.61616162, -0.5959596, -0.57575758, -0.55555556, -0.53535354, -0.51515152, -0.49494949, -0.47474747, -0.45454545, -0.43434343, -0.41414141, -0.39393939, -0.37373737, -0.35353535, -0.33333333, -0.31313131, -0.29292929, -0.27272727, -0.25252525, -0.23232323, -0.21212121, -0.19191919, -0.17171717, -0.15151515, -0.13131313, -0.11111111, -0.09090909, -0.07070707, -0.05050505, -0.03030303, -0.01010101, 0.01010101, 0.03030303, 0.05050505, 0.07070707, 0.09090909, 0.11111111, 0.13131313, 0.15151515, 0.17171717, 0.19191919, 0.21212121, 0.23232323, 0.25252525, 0.27272727, 0.29292929, 0.31313131, 0.33333333, 0.35353535, 0.37373737, 0.39393939, 0.41414141, 0.43434343, 0.45454545, 0.47474747, 0.49494949, 0.51515152, 0.53535354, 0.55555556, 0.57575758, 0.5959596, 0.61616162, 0.63636364, 0.65656566, 0.67676768, 0.6969697, 0.71717172, 0.73737374, 0.75757576, 0.77777778, 0.7979798, 0.81818182, 0.83838384, 0.85858586, 0.87878788, 0.8989899, 0.91919192, 0.93939394, 0.95959596, 0.97979798, 1.]), phi=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the gamma-gamma absorption produced by a spherical shell BLR for a general set of model parameters

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the tau note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • L_disk (Quantity) – Luminosity of the disk whose radiation is being reprocessed by the BLR

  • xi_line (float) – fraction of the disk radiation reprocessed by the BLR

  • epsilon_line (string) – dimensionless energy of the emitted line

  • R_line (Quantity) – radius of the BLR spherical shell

  • r (Quantity) – distance between the Broad Line Region and the blob

  • l_size (int) – size of the array of distances from the BH to integrate over

  • mu (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

  • phi (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

Returns:

array of the tau values corresponding to each frequency

Return type:

Quantity

static evaluate_tau_blr_mu_s(nu, z, mu_s, L_disk, xi_line, epsilon_line, R_line, r, u_size=100, mu=array([-1., -0.97979798, -0.95959596, -0.93939394, -0.91919192, -0.8989899, -0.87878788, -0.85858586, -0.83838384, -0.81818182, -0.7979798, -0.77777778, -0.75757576, -0.73737374, -0.71717172, -0.6969697, -0.67676768, -0.65656566, -0.63636364, -0.61616162, -0.5959596, -0.57575758, -0.55555556, -0.53535354, -0.51515152, -0.49494949, -0.47474747, -0.45454545, -0.43434343, -0.41414141, -0.39393939, -0.37373737, -0.35353535, -0.33333333, -0.31313131, -0.29292929, -0.27272727, -0.25252525, -0.23232323, -0.21212121, -0.19191919, -0.17171717, -0.15151515, -0.13131313, -0.11111111, -0.09090909, -0.07070707, -0.05050505, -0.03030303, -0.01010101, 0.01010101, 0.03030303, 0.05050505, 0.07070707, 0.09090909, 0.11111111, 0.13131313, 0.15151515, 0.17171717, 0.19191919, 0.21212121, 0.23232323, 0.25252525, 0.27272727, 0.29292929, 0.31313131, 0.33333333, 0.35353535, 0.37373737, 0.39393939, 0.41414141, 0.43434343, 0.45454545, 0.47474747, 0.49494949, 0.51515152, 0.53535354, 0.55555556, 0.57575758, 0.5959596, 0.61616162, 0.63636364, 0.65656566, 0.67676768, 0.6969697, 0.71717172, 0.73737374, 0.75757576, 0.77777778, 0.7979798, 0.81818182, 0.83838384, 0.85858586, 0.87878788, 0.8989899, 0.91919192, 0.93939394, 0.95959596, 0.97979798, 1.]), phi=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the gamma-gamma absorption produced by a spherical shell BLR for a general set of model parameters and arbitrary mu_s

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the tau note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • L_disk (Quantity) – Luminosity of the disk whose radiation is being reprocessed by the BLR

  • xi_line (float) – fraction of the disk radiation reprocessed by the BLR

  • epsilon_line (string) – dimensionless energy of the emitted line

  • R_line (Quantity) – radius of the BLR spherical shell

  • r (Quantity) – distance between the Broad Line Region and the blob

  • l_size (int) – size of the array of distances from the BH to integrate over

  • mu (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

  • phi (ndarray) – arrays of cosine of zenith and azimuth angles to integrate over

Returns:

array of the tau values corresponding to each frequency

Return type:

Quantity

static evaluate_tau_dt(nu, z, mu_s, L_disk, xi_dt, epsilon_dt, R_dt, r, l_size=50, phi=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the gamma-gamma absorption produced by a ring dust torus

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • L_disk (Quantity) – Luminosity of the disk whose radiation is being reprocessed by the BLR

  • xi_dt (float) – fraction of the disk radiation reprocessed by the disk

  • epsilon_dt (string) – peak (dimensionless) energy of the black body radiated by the torus

  • R_dt (Quantity) – radius of the ting-like torus

  • r (Quantity) – distance between the dust torus and the blob

  • l_size (int) – size of the array of distances from the BH to integrate over

  • phi (ndarray) – arrays of azimuth angles to integrate over

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

static evaluate_tau_dt_mu_s(nu, z, mu_s, L_disk, xi_dt, epsilon_dt, R_dt, r, u_size=100, phi_re=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the gamma-gamma absorption produced by a ring dust torus for the case of photon moving at an angle to the jet

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the sed note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • L_disk (Quantity) – Luminosity of the disk whose radiation is being reprocessed by the BLR

  • xi_dt (float) – fraction of the disk radiation reprocessed by the disk

  • epsilon_dt (string) – peak (dimensionless) energy of the black body radiated by the torus

  • R_dt (Quantity) – radius of the ting-like torus

  • r (Quantity) – distance between the dust torus and the blob

  • u_size (int) – size of the array of distances from the photon origin to integrate over

  • phi_re (ndarray) – arrays of azimuth angles of the dust torus to integrate over

Returns:

array of the SED values corresponding to each frequency

Return type:

Quantity

static evaluate_tau_ps_behind_blob(nu, z, mu_s, epsilon_0, L_0, r)[source]

Evaluates the absorption produced by the photon field of a point source of photons behind the blob, for a general set of model parameters

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the opacity note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • epsilon_0 (float) – dimensionless energy (in electron rest mass energy units) of the target photon field

  • L_0 (Quantity) – luminosity [erg cm-3] of the point source behind the jet

  • r (Quantity) – distance between the point source and the blob

Returns:

array of the tau values corresponding to each frequency

Return type:

Quantity

static evaluate_tau_ps_behind_blob_mu_s(nu, z, mu_s, epsilon_0, L_0, r, u_size=100)[source]

Evaluates the absorption produced by the photon field of a point source of photons behind the blob, for a general set of model parameters

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the opacity note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • epsilon_0 (float) – dimensionless energy (in electron rest mass energy units) of the target photon field

  • L_0 (Quantity) – luminosity [erg cm-3] of the point source behind the jet

  • r (Quantity) – distance between the point source and the blob

  • u_size (int) – size of the array of distances from the photon origin to integrate over

Returns:

array of the tau values corresponding to each frequency

Return type:

Quantity

static evaluate_tau_ss_disk(nu, z, mu_s, M_BH, L_disk, eta, R_in, R_out, r, R_tilde_size=100, l_tilde_size=50, phi=array([0., 0.12822827, 0.25645654, 0.38468481, 0.51291309, 0.64114136, 0.76936963, 0.8975979, 1.02582617, 1.15405444, 1.28228272, 1.41051099, 1.53873926, 1.66696753, 1.7951958, 1.92342407, 2.05165235, 2.17988062, 2.30810889, 2.43633716, 2.56456543, 2.6927937, 2.82102197, 2.94925025, 3.07747852, 3.20570679, 3.33393506, 3.46216333, 3.5903916, 3.71861988, 3.84684815, 3.97507642, 4.10330469, 4.23153296, 4.35976123, 4.48798951, 4.61621778, 4.74444605, 4.87267432, 5.00090259, 5.12913086, 5.25735913, 5.38558741, 5.51381568, 5.64204395, 5.77027222, 5.89850049, 6.02672876, 6.15495704, 6.28318531]))[source]

Evaluates the gamma-gamma absorption produced by the photon field of a Shakura-Sunyaev accretion disk

Parameters:
  • nu (Quantity) – array of frequencies, in Hz, to compute the opacity note these are observed frequencies (observer frame)

  • z (float) – redshift of the source

  • mu_s (float) – cosine of the angle between the blob motion and the jet axis

  • M_BH (Quantity) – Black Hole mass

  • L_disk (Quantity) – luminosity of the disk

  • eta (float) – accretion efficiency

  • R_in (Quantity) – inner disk radius

  • R_out (Quantity) – inner disk radius

  • R_tilde_size (int) – size of the array of disk coordinates to integrate over

  • r (Quantity) – distance between the point source and the blob

  • l_tilde_size (int) – size of the array of distances from the BH to integrate over

  • phi (ndarray) – array of azimuth angles to integrate over

Returns:

array of the tau values corresponding to each frequency

Return type:

Quantity

set_l(l_size=50)[source]
set_mu(mu_size=100)[source]
set_phi(phi_size=50)[source]

Set array of azimuth angles to integrate over

tau(nu)[source]

optical depth

\[\tau_{\gamma \gamma}(\nu)\]
Parameters:

nu (~astropy.units.Quantity) – array of frequencies, in Hz, to compute the opacity, note these are observed frequencies (observer frame).

tau_blr(nu)[source]

Evaluates the gamma-gamma absorption produced by a spherical shell BLR for a general set of model parameters

tau_blr_mu_s(nu)[source]

Evaluates the gamma-gamma absorption produced by a spherical shell BLR for a general set of model parameters and arbitrary mu_s

tau_dt(nu)[source]

evaluates the gamma-gamma absorption produced by a ring dust torus

tau_dt_mu_s(nu)[source]

evaluates the gamma-gamma absorption produced by a ring dust torus

tau_on_synchrotron(blob, nu, nu_s_size=200, delta_margin_low=0.01)[source]

Optical depth for absorption of gamma rays in synchrotron radiation of the blob. It assumes the same radiation field as the SSC class.

Parameters:
  • blob (Blob) – emission region and electron distribution hitting the photon target

  • nu (Quantity) – array of frequencies, in Hz, to compute the opacity note these are observed frequencies (observer frame)

  • nu_s_size (int) – size of the array over the synchrotron frequencies

  • delta_margin_low (float) – extension of the integration range of the synchrotron radiation beyond the delta approximation, default = 0.01, but lower value might be needed if the calculations are performed up to very high energies

tau_ps_behind_blob(nu)[source]

Evaluates the absorption produced by the photon field of a point source of photons behind the blob

tau_ps_behind_blob_mu_s(nu)[source]

Evaluates the absorption produced by the photon field of a point source of photons behind the blob

tau_ss_disk(nu)[source]

Evaluates the gamma-gamma absorption produced by the photon field of a Shakura-Sunyaev accretion disk

agnpy.constraints module

class agnpy.constraints.SpectralConstraints(blob)[source]

Bases: object

Class to describe the self-consistency constraints on the electron energy distribution

Parameters:

blob

gamma_break_EC_DT(dt, r=<Quantity 0. cm>)[source]

Simple estimation of the cooling break of electrons comparing EC time scale (see B&G 1970) with the ballistic time scale, like in gamma_break_SSC WARNING: assumes Thomson regime

\[\gamma_b = 3 m_e c^2 / 4 \sigma_T U'_{\mathrm{ext}} R_b\]
property gamma_break_SSC

Simple estimation of the cooling break of electrons comparing SSC time scale (see B&G 1970) with the ballistic time scale: WARNING: only applicable in Thomson regime

\[\begin{split}T_{\mathrm{SSC}} &= E\,/\,(\mathrm{d}E/\mathrm{d}t)_{\mathrm{SSC}} = 3 m_e c^2 / (4 \sigma_T U_{\mathrm{SSC}} \gamma) \\ T_{\mathrm{bal}} &= R_b / c \\ T_{\mathrm{SSC}} &= T_{\mathrm{bal}} \Rightarrow \gamma_b = 3 m_e c^2 / 4 \sigma_T U_{\mathrm{SSC}} R_b\end{split}\]
property gamma_break_synch

Simple estimation of the cooling break of electrons comparing synchrotron cooling time scale with the ballistic time scale:

\[\begin{split}T_{\mathrm{synch}} &= E\,/\,(\mathrm{d}E/\mathrm{d}t)_{\mathrm{synch}} = 3 m_e c^2 / (4 \sigma_T U_B \gamma) \\ T_{\mathrm{bal}} &= R_b / c \\ T_{\mathrm{synch}} &= T_{\mathrm{bal}} \Rightarrow \gamma_b = 6 \pi m_e c^2 / \sigma_T B^2 R_b\end{split}\]
gamma_max_EC_DT(dt, r=<Quantity 0. cm>)[source]

Simple estimation of maximum Lorentz factor of electrons comparing the acceleration time scale with the EC energy loss (in Thomson range, see B&G 1970), like in gamma_max_SSC WARNING: assumes Thomson regime

\[\gamma_{\mathrm{max}} = \sqrt{\frac{3 \xi e B }{ \sigma_T U'_\mathrm{ext}}}\]
property gamma_max_SSC

Simple estimation of maximum Lorentz factor of electrons comparing the acceleration time scale with the SSC energy loss (in Thomson range) WARNING: the highest energy electrons will most often scatter in Klein-Nishina range instead

\[\begin{split}(\mathrm{d}E/\mathrm{d}t)_{\mathrm{acc}} &= \xi c E / R_L \\ (\mathrm{d}E/\mathrm{d}t)_{\mathrm{SSC}} &= 4 / 3 \sigma_T c U_{\mathrm{synch}} \gamma^2 \\ (\mathrm{d}E/\mathrm{d}t)_{\mathrm{acc}} &= (\mathrm{d}E/\mathrm{d}t)_{\mathrm{SSC}} \Rightarrow \gamma_{\mathrm{max}} < \sqrt{\frac{3 \xi e B }{\sigma_T U_SSC}}\end{split}\]
property gamma_max_ballistic

Naive estimation of maximum Lorentz factor of electrons comparing acceleration time scale with ballistic time scale. For the latter we assume that the particles crosses the blob radius.

\[\begin{split}(\mathrm{d}E/\mathrm{d}t)_{\mathrm{acc}} &= \xi c E / R_L \\ T_{\mathrm{acc}} &= E \,/\,(\mathrm{d}E/\mathrm{d}t)_{\mathrm{acc}} = R_L / (\xi c) \\ T_{\mathrm{bal}} &= R_b / c \\ T_{\mathrm{acc}} &< T_{\mathrm{bal}} \Rightarrow \gamma_{\mathrm{max}} < \frac{\xi R_b e B}{m_e c^2}\end{split}\]
property gamma_max_larmor

maximum Lorentz factor of electrons that have their Larmour radius smaller than the blob radius: \(R_L < R_b\). The Larmor frequency and radius in Gaussian units read

\[\begin{split}\omega_L &= \frac{eB}{\gamma m_e c} \\ R_L &= \frac{v}{\omega_L} = \frac{\gamma m_e v c}{e B} \approx \frac{\gamma m_e c^2}{e B}\end{split}\]

therefore

\[R_L < R_b \Rightarrow \gamma_{\mathrm{max}} < \frac{R_b e B}{m_e c^2}\]
property gamma_max_synch

Simple estimation of maximum Lorentz factor of electrons comparing the acceleration time scale with the synchrotron energy loss

\[\begin{split}(\mathrm{d}E/\mathrm{d}t)_{\mathrm{acc}} &= \xi c E / R_L \\ (\mathrm{d}E/\mathrm{d}t)_{\mathrm{synch}} &= 4 / 3 \sigma_T c U_B \gamma^2 \\ (\mathrm{d}E/\mathrm{d}t)_{\mathrm{acc}} &= (\mathrm{d}E/\mathrm{d}t)_{\mathrm{synch}} \Rightarrow \gamma_{\mathrm{max}} < \sqrt{\frac{6 \pi \xi e}{\sigma_T B}}\end{split}\]

agnpy.fit module