Physical setup

import jetset
print('tested on jetset',jetset.__version__)
tested on jetset 1.2.2

In this section we describe how to build a model of jet able to reproduce SSC/EC emission processes, using the Jet class from the jet_model module. This class, through a flexible and intuitive interface, allows to access the C numerical code that provides an accurate and fast computation of the synchrotron and inverse Compton processes.

Basic setup

A jet instance can be built using the the Jet class, instantiating the object in the following way:

from jetset.jet_model import Jet
my_jet=Jet(name='test',electron_distribution='lppl',)
This instruction will create:
  • a Jet object with name test,

  • using as electron distribution the lppl model, that is a log-parabola with a low-energy power-law branch.

(a working directory jet_wd will be created, this directory can be deleted when all the processes of your script are done)

For a list of possible distribution you can run the command

Jet.available_electron_distributions()
lp: log-parabola
pl: powerlaw
lppl: log-parabola with low-energy powerlaw branch
lpep: log-parabola defined by peak energy
plc: powerlaw with cut-off
bkn: broken powerlaw
superexp: powerlaw with super-exp cut-off

to view all the paramters:

my_jet.show_pars()
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse

custom electron distributions can be created by the user as described in this section of the tutorial Custom emitters distribution

Each parameter has a default value. All the parameters listed are handled by ModelParameterArray, and each parameter is an instance of the the JetParameter. class. These parameters can be visualized by the command

my_jet.parameters
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
None

and the corresponding astropy table with units can be accessed by: my_jet.parameters.par_table

This means that you can easily convert the values in the table using the units module of astropy.

Warning

Please note, that the table is built on the fly from the ModelParameterArray and each modification you do to this table will not be reflected on the actual parameters values

To get a full description of the model you can use the instruction

my_jet.show_model()
--------------------------------------------------------------------------------
model description:
--------------------------------------------------------------------------------
type: Jet
name: test

electrons distribution:
 type: lppl
 gamma energy grid size:  201
 gmin grid : 2.000000e+00
 gmax grid : 1.000000e+06
 normalization:  True
 log-values:  False
 ratio of cold protons to relativistic electrons: 1.000000e-01

radiative fields:
 seed photons grid size:  100
 IC emission grid size:  100
 source emissivity lower bound :  1.000000e-120
 spectral components:
   name:Sum, state: on
   name:Sync, state: self-abs
   name:SSC, state: on
external fields transformation method: blob

SED info:
 nu grid size jetkernel: 1000
 nu size: 500
 nu mix (Hz): 1.000000e+06
 nu max (Hz): 1.000000e+30

flux plot lower bound   :  1.000000e-30

--------------------------------------------------------------------------------
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
--------------------------------------------------------------------------------

as you can notice, you can now access further information regarding the model, such as numerical configuration of the grid. These numerical parameters will be discussed in the :ref:`jet_numerical_guide’ section

If you want to use a cosmology model different from the default one please read the Choosing a cosmology model section.

Warning

Starting from version 1.1.0, the R parameter as default is linear and not logarithmic, please update your old scripts setting R with linear values.

Setting the parameters

Assume you want to change some of the parameters in your model, you can use two methods:

  1. using the Jet.set_par() method

my_jet.set_par('B',val=0.2)
my_jet.set_par('gamma0_log_parab',val=5E3)
my_jet.set_par('gmin',val=1E2)
my_jet.set_par('gmax',val=1E8)
my_jet.set_par('R',val=1E15)
my_jet.set_par('N',val=1E3)
  1. accessing directly the parameter

my_jet.parameters.B.val=0.2
my_jet.parameters.r.val=0.4

Investigating the electron distribution

for setting custom electron distributions can be created by the user as described in this section of the tutorial Custom emitters distribution

my_jet.show_electron_distribution()
--------------------------------------------------------------------------------
electrons distribution:
 type: lppl
 gamma energy grid size:  201
 gmin grid : 2.000000e+00
 gmax grid : 1.000000e+06
 normalization  True
 log-values  False
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testBmagnetic_fieldgauss2.000000e-011.000000e-101.000000e+10FalseFalse
testNemitters_density1 / cm31.000000e+030.000000e+00--FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testRregion_sizecm1.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*5.000000e+031.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+081.000000e+001.000000e+15FalseFalse
testgminlow-energy-cut-offlorentz-factor*1.000000e+021.000000e+001.000000e+09FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
p=my_jet.electron_distribution.plot3p()
../../../_images/Jet_example_phys_SSC_29_0.png
p=my_jet.electron_distribution.plot3p(energy_unit='eV')
../../../_images/Jet_example_phys_SSC_30_0.png
p=my_jet.electron_distribution.plot2p(energy_unit='erg')
../../../_images/Jet_example_phys_SSC_31_0.png

to obtain a loglog plot, pass loglog=True to the plot method

p=my_jet.electron_distribution.plot(energy_unit='erg',loglog=True)
../../../_images/Jet_example_phys_SSC_33_0.png
import numpy as np
p=None
for r in np.linspace(0.3,1,10):
    my_jet.parameters.r.val=r
    _l='r=%2.2f'%r
    if p is None:
        p=my_jet.electron_distribution.plot3p(label=_l)
    else:
        p=my_jet.electron_distribution.plot3p(p,label=_l)
../../../_images/Jet_example_phys_SSC_34_0.png

Using log values for electron distribution parameters

my_jet=Jet(name='test',electron_distribution='lppl',electron_distribution_log_values=True)
my_jet.show_model()
--------------------------------------------------------------------------------
model description:
--------------------------------------------------------------------------------
type: Jet
name: test

electrons distribution:
 type: lppl
 gamma energy grid size:  201
 gmin grid : 2.000000e+00
 gmax grid : 1.000000e+06
 normalization:  True
 log-values:  True
 ratio of cold protons to relativistic electrons: 1.000000e-01

radiative fields:
 seed photons grid size:  100
 IC emission grid size:  100
 source emissivity lower bound :  1.000000e-120
 spectral components:
   name:Sum, state: on
   name:Sync, state: self-abs
   name:SSC, state: on
external fields transformation method: blob

SED info:
 nu grid size jetkernel: 1000
 nu size: 500
 nu mix (Hz): 1.000000e+06
 nu max (Hz): 1.000000e+30

flux plot lower bound   :  1.000000e-30

--------------------------------------------------------------------------------
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*3.010300e-010.000000e+009.000000e+00TrueFalse
testgmaxhigh-energy-cut-offlorentz-factor*6.000000e+000.000000e+001.500000e+01TrueFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*4.000000e+000.000000e+009.000000e+00TrueFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
--------------------------------------------------------------------------------

Evaluate and plot the model

At this point we can evaluate the emission for this jet model using the instruction

my_jet.eval()
my_jet.show_pars()
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*3.010300e-010.000000e+009.000000e+00TrueFalse
testgmaxhigh-energy-cut-offlorentz-factor*6.000000e+000.000000e+001.500000e+01TrueFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*4.000000e+000.000000e+009.000000e+00TrueFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse

and plot the corresponding SED:

Warning

Starting from version 1.2.0 The rescale method as been replaced by the setlim methd. Please notice, that now jetset uses as defualt logarthmic axis rather than loglog plots, so, the correct way to use it is rescale(x_min=8)->setlim(x_min=1E8)

from jetset.plot_sedfit import PlotSED
my_plot=PlotSED()
my_jet.plot_model(plot_obj=my_plot)
my_plot.setlim(y_min=10**-17.5)
../../../_images/Jet_example_phys_SSC_43_0.png

alternatively, you can call the plot_model method without passing a Plot object

my_plot=my_jet.plot_model()
my_plot.setlim(y_min=10**-17.5)
../../../_images/Jet_example_phys_SSC_45_0.png

If you want to have more points on the IC spectrum you can set the numerical parameters for radiative fields(see :ref:`jet_numerical_guide’ section for more details):

my_jet.set_IC_nu_size(100)
my_jet.eval()
my_plot=my_jet.plot_model()
my_plot.setlim(y_min=10**-17.5)
../../../_images/Jet_example_phys_SSC_48_0.png

you can access the same plot, but in the rest frame of the black hole, or accretion disk, hence plotting the isotropic luminosity, by simply passing the frame kw to src

my_plot=my_jet.plot_model(frame='src')
my_plot.setlim(y_max=1E42,y_min=1E38)
../../../_images/Jet_example_phys_SSC_50_0.png

the my_plot object returned will be built on the fly by the plot_model method

Starting from version 1.2.0 you can also plot in the Fnu or Lnu representation adding the density=True keyword to the plot_model command

my_plot=my_jet.plot_model(frame='src',density=True)
my_plot.setlim(y_max=1E29,y_min=1E11,x_min=1E8,x_max=1E28)
../../../_images/Jet_example_phys_SSC_53_0.png

Changing the nu grid

The SED info header displayed by the .show_model() methods reports information for the SED nu_min, nu_max, nu_size and nu_grid_size. - The nu_grid_size is the internal interpolation grid used by jetkernel C code, and it should not be changed - The nu_size is the python interpolation grid used by python wrapper on top of the jetkernel one, and is used only for the SEDs production and plotting. - nu_min and nu_max, are used for the boundaries of the model, and can be changed if the custom value does not cover your expected range. Notice that if the model is below the source emissivity or flux plot lower bound, then your changes on nu_min/nu_max will have no effect.

my_jet.nu_min=1E10
my_jet.nu_size=400
my_jet.nu_max=1E30
my_jet.show_model()
--------------------------------------------------------------------------------
model description:
--------------------------------------------------------------------------------
type: Jet
name: test

electrons distribution:
 type: lppl
 gamma energy grid size:  201
 gmin grid : 2.000000e+00
 gmax grid : 1.000000e+06
 normalization:  True
 log-values:  True
 ratio of cold protons to relativistic electrons: 1.000000e-01

radiative fields:
 seed photons grid size:  100
 IC emission grid size:  100
 source emissivity lower bound :  1.000000e-120
 spectral components:
   name:Sum, state: on
   name:Sync, state: self-abs
   name:SSC, state: on
external fields transformation method: blob

SED info:
 nu grid size jetkernel: 1000
 nu size: 400
 nu mix (Hz): 1.000000e+10
 nu max (Hz): 1.000000e+30

flux plot lower bound   :  1.000000e-30

--------------------------------------------------------------------------------
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*3.010300e-010.000000e+009.000000e+00TrueFalse
testgmaxhigh-energy-cut-offlorentz-factor*6.000000e+000.000000e+001.500000e+01TrueFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*4.000000e+000.000000e+009.000000e+00TrueFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
--------------------------------------------------------------------------------
my_jet.eval()
p=my_jet.plot_model()
p.setlim(x_min=1E7,y_min=10**-17.5)
../../../_images/Jet_example_phys_SSC_57_0.png

if you want to to have interacitve plot:

  1. in a jupyter notebook use:

%matplotlib notebook
  1. in jupyter lab:

%matplotlib widget
(visit this url to setup and install: https://github.com/matplotlib/ipympl)
  1. in an ipython terminal

from matplotlib import pylab as plt
plt.ion()

Comparing models on the same plot

to compare the same model after changing a parameter

my_jet=Jet(name='test',electron_distribution='lppl',)
my_jet.set_par('B',val=0.2)
my_jet.set_par('gamma0_log_parab',val=5E3)
my_jet.set_par('gmin',val=1E2)
my_jet.set_par('gmax',val=1E8)
my_jet.set_par('R',val=10**14.5)
my_jet.set_par('N',val=1E3)

my_jet.parameters.gamma0_log_parab.val=1E4
my_jet.eval()
my_plot=my_jet.plot_model(label='gamma0_log_parab=1E4',comp='Sum')
my_jet.set_par('gamma0_log_parab',val=1.0E5)
my_jet.eval()
my_plot=my_jet.plot_model(my_plot,label='gamma0_log_parab=1E5',comp='Sum')
my_plot.setlim(y_max=1E-13,y_min=2E-17,x_min=1E8)
../../../_images/Jet_example_phys_SSC_61_0.png

Saving a plot

to save the plot

my_plot.save('jet1.png')

Saving and loading a model

Warning

starting from version 1.1.0 the saved model format has changed, if you have models saved with version<1.1.0, please update them the new models by loading the old models with the Jet.load_old_model() and then saving them again.

my_jet.save_model('test_model.pkl')
my_jet_new=Jet.load_model('test_model.pkl')
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testgminlow-energy-cut-offlorentz-factor*1.000000e+021.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+081.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+030.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+051.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
testRregion_sizecm3.162278e+141.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss2.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse

Switching on/off the particle distribution normalization

As default the electron distributions are normalized, i.e. are multiplied by a constant N_0, in such a way that :

\(\int_{\gamma_{min}}^{\gamma_{max}} n(\gamma) d\gamma =1\),

it means the the value N, refers to the actual density of emitters. If you want to chance this behavior, you can start looking at the sate of Norm_distr flag with the following command

my_jet.Norm_distr
True

and then you can switch off the normalization withe command

my_jet.switch_Norm_distr_OFF()

OR

my_jet.Norm_distr=0
my_jet.switch_Norm_distr_ON()

OR

my_jet.Norm_distr=1

Setting the particle density from observed Fluxes or Luminosities

It is possible to set the density of emitting particles starting from some observed luminosity or flux (see the method Jet.set_N_from_nuFnu(), and Jet.set_N_from_nuLnu())

my_jet=Jet(name='test',electron_distribution='lppl')

this is the initial value of N

my_jet.parameters.N.val
100.0

we now want to set the value of N in order that the observed synchrotron flux at a given frequency matches a target value. For example, assume that we wish to set N in order that the synchrotron flux at \(10^{15}\) Hz is exactly matching the desired value of \(10^{-14}\) ergs cm-2 s-1. We can accomplish this by using the method Jet.set_N_from_nuFnu() as follows:

my_jet.set_N_from_nuFnu(nuFnu_obs=1E-14,nu_obs=1E15)

This is the updated value of N, obtained in order to match the given flux at the given frequency

my_jet.get_par_by_name('N').val
272.37555111028803

OR

my_jet.parameters.N.val
272.37555111028803
my_jet.parameters.show_pars()
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm32.723756e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
my_jet.eval()
my_plot=my_jet.plot_model(label='set N from F=1E-14')
my_plot.setlim(y_max=1E-13,y_min=2E-17,x_min=1E8)
../../../_images/Jet_example_phys_SSC_92_0.png

as you can see, the synchrotron flux at \(10^{15}\) Hz, now exactly matches the desired value of \(10^{-14}\) ergs cm-2 s-1. Alternatively, the value of N can be obtained using the rest-frame luminosity and frequency, using the method Jet.set_N_from_nuLnu()

my_jet.set_N_from_nuLnu(nuLnu_src=1E43,nu_src=1E15)

where nuLnu_src is the source rest-frame isotropic luminosity in erg/s at the rest-frame frequency nu_src in Hz.

Setting the beaming factor and expression

Important

Starting from version 1.2.0, when using delta expression, the value of delta used to compute jet luminosities will be set to beam_obj. In previous version a reference value of 10 was used. In any case, if you are interseted in evaluating jet luminosities you should use the beaming_expr method

It is possible to set the beaming factor according to the relativistic BulkFactor and viewing angle, this can be done by setting the beaming_expr kw in the Jet constructor, possible choices are

  • delta to provide directly the beaming factor (default)

  • bulk_theta to provide the BulkFactor and the jet viewing angle

my_jet=Jet(name='test',electron_distribution='lppl',beaming_expr='bulk_theta')
my_jet.parameters.show_pars()
Table length=13
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testthetajet-viewing-angledeg1.000000e-010.000000e+00--FalseFalse
testBulkFactorjet-bulk-factorlorentz-factor*1.000000e+011.000000e+001.000000e+05FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse

the actual value of the beaming factor can be obtained using the Jet.get_beaming()

my_jet.get_beaming()
19.943844732554165

We can change the value of theta and get the updated value of the beaming factor

my_jet.set_par('theta',val=10.)
my_jet.get_beaming()
4.968041140891955

of course setting beaming_expr=delta we get the same beaming expression as in the default case

my_jet=Jet(name='test',electron_distribution='lppl',beaming_expr='delta')
my_jet.parameters.show_pars()
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse

Switch ON/OFF Synchrotron sefl-absorption and IC emission

my_jet.show_model()
--------------------------------------------------------------------------------
model description:
--------------------------------------------------------------------------------
type: Jet
name: test

electrons distribution:
 type: lppl
 gamma energy grid size:  201
 gmin grid : 2.000000e+00
 gmax grid : 1.000000e+06
 normalization:  True
 log-values:  False
 ratio of cold protons to relativistic electrons: 1.000000e-01

radiative fields:
 seed photons grid size:  100
 IC emission grid size:  100
 source emissivity lower bound :  1.000000e-120
 spectral components:
   name:Sum, state: on
   name:Sync, state: self-abs
   name:SSC, state: on
external fields transformation method: blob

SED info:
 nu grid size jetkernel: 1000
 nu size: 500
 nu mix (Hz): 1.000000e+06
 nu max (Hz): 1.000000e+30

flux plot lower bound   :  1.000000e-30

--------------------------------------------------------------------------------
Table length=12
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testbeam_objbeaminglorentz-factor*1.000000e+011.000000e-041.000000e+04FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
--------------------------------------------------------------------------------

as you see the state of Sync emission is self-abs, we can check accessing the specific spectral component state, and get the allowed states value

my_jet.spectral_components.Sync.show()
name                : Sync
var name            : do_Sync
state               : self-abs
allowed states : ['on', 'off', 'self-abs']
my_jet.spectral_components.Sync.state='on'

now the sate is ‘on’ with no ‘self-abs’

my_jet.eval()
p=my_jet.plot_model()
p.setlim(y_max=1E-13,y_min=5E-18,x_min=1E8)
../../../_images/Jet_example_phys_SSC_116_0.png

to re-enable

my_jet.spectral_components.Sync.state='self-abs'
my_jet.eval()
p=my_jet.plot_model()
p.setlim(y_max=1E-13,y_min=5E-18,x_min=1E8)
../../../_images/Jet_example_phys_SSC_118_0.png
my_jet.spectral_components.SSC.show()
name                : SSC
var name            : do_SSC
state               : on
allowed states : ['on', 'off']
my_jet.spectral_components.SSC.state='off'
my_jet.eval()
p=my_jet.plot_model()
p.setlim(y_max=1E-13,y_min=8E-18,x_min=1E8)
../../../_images/Jet_example_phys_SSC_120_0.png

to re-enable

my_jet.spectral_components.SSC.state='on'
my_jet.eval()
p=my_jet.plot_model()
p.setlim(y_max=1E-13,y_min=8E-18,x_min=1E8)
../../../_images/Jet_example_phys_SSC_122_0.png

Accessing individual spectral components

It is possible to access specific spectral components of our model

my_jet=Jet(name='test',electron_distribution='lppl',beaming_expr='bulk_theta')
my_jet.eval()

We can obtain this information anytime using the Jet.list_spectral_components() method

my_jet.list_spectral_components()
Sum
Sync
SSC

the on-screen message is telling us which components have been evaluated.

and we cann access a specific component using the Jet.get_spectral_component_by_name() method

Sync=my_jet.get_spectral_component_by_name('Sync')

OR

Sync=my_jet.spectral_components.Sync

and from the SED object we can extract both the nu and nuFnu array

nu_sync=Sync.SED.nu
nuFnu_sync=Sync.SED.nuFnu
print (nuFnu_sync[::10])
[0.00000000e+00 0.00000000e+00 0.00000000e+00 6.04250670e-26
 2.16351829e-24 9.84432972e-23 4.74613296e-21 2.28931297e-19
 1.09662087e-17 1.83733916e-16 4.11135769e-16 7.21745036e-16
 1.25581697e-15 2.18363181e-15 3.79383567e-15 6.57833387e-15
 1.13501032e-14 1.93585563e-14 3.21429895e-14 5.06938061e-14
 7.36908738e-14 9.77112603e-14 1.17645633e-13 1.28621805e-13
 1.26850509e-13 1.10646286e-13 7.82537850e-14 3.17631756e-14
 2.39710785e-15 8.88519981e-19 7.47780581e-29 0.00000000e+00
 0.00000000e+00 0.00000000e+00 0.00000000e+00 0.00000000e+00
 0.00000000e+00 0.00000000e+00 0.00000000e+00 0.00000000e+00
 0.00000000e+00 0.00000000e+00 0.00000000e+00 0.00000000e+00
 0.00000000e+00 0.00000000e+00 0.00000000e+00 0.00000000e+00
 0.00000000e+00 0.00000000e+00] erg / (cm2 s)

or for the src rest frame (isotropic luminosity)

nu_sync_src=Sync.SED.nu_src
nuLnu_sync_src=Sync.SED.nuLnu_src
print (nuLnu_sync_src[::10])
[0.00000000e+00 0.00000000e+00 0.00000000e+00 1.63219228e+30
 5.84406112e+31 2.65913465e+33 1.28201787e+35 6.18385569e+36
 2.96217481e+38 4.96299126e+39 1.11055338e+40 1.94956618e+40
 3.39219277e+40 5.89839143e+40 1.02478484e+41 1.77692906e+41
 3.06587177e+41 5.22910236e+41 8.68241307e+41 1.36933301e+42
 1.99052613e+42 2.63936099e+42 3.17782509e+42 3.47431170e+42
 3.42646573e+42 2.98875984e+42 2.11377876e+42 8.57981835e+41
 6.47502951e+40 2.40005601e+37 2.01989298e+27 0.00000000e+00
 0.00000000e+00 0.00000000e+00 0.00000000e+00 0.00000000e+00
 0.00000000e+00 0.00000000e+00 0.00000000e+00 0.00000000e+00
 0.00000000e+00 0.00000000e+00 0.00000000e+00 0.00000000e+00
 0.00000000e+00 0.00000000e+00 0.00000000e+00 0.00000000e+00
 0.00000000e+00 0.00000000e+00] erg / s

Moreover, you can access the corresponding astropy table

my_jet.spectral_components.build_table(restframe='obs')
t_obs=my_jet.spectral_components.table
t_obs[::10]
Table length=50
nuSumSyncSSC
Hzerg / (cm2 s)erg / (cm2 s)erg / (cm2 s)
float64float64float64float64
1000000.00.00.00.0
3026648.0593956890.00.00.0
9160598.475443710.00.00.0
27725907.598604816.042506698961876e-266.042506698961876e-260.0
83916564.428301622.1635183262103024e-242.1635182921864927e-243.4023433519560716e-32
253985906.878072929.844329885880581e-239.844329720868366e-231.6498676086277722e-30
768725952.16637214.7461330192800106e-214.746132957910235e-216.136928342027515e-29
2326662911.3314582.289313003537184e-192.289312967845199e-193.569167611890846e-27
7041989785.4492961.096620888827693e-171.0966208685756238e-172.0247793809032263e-25
............
1.7317171337233599e+252.727651989047602e-150.02.727651989047602e-15
5.2412983022060615e+252.2657229497733764e-150.02.2657229497733764e-15
1.5863565335085865e+261.376097972965168e-150.01.376097972965168e-15
4.801342923653465e+265.034239483778965e-160.05.034239483778965e-16
1.4531975242368953e+273.2430360747827925e-960.03.2430360747827925e-96
4.3983174666502106e+270.00.00.0
1.3312159025043105e+280.00.00.0
4.029122027951344e+280.00.00.0
1.2194734366967333e+290.00.00.0
3.690916910662782e+290.00.00.0

and also in the src restframe

my_jet.spectral_components.build_table(restframe='src')
t_src=my_jet.spectral_components.table
t_src[::10]
Table length=50
nuSumSyncSSC
Hzerg / serg / serg / s
float64float64float64float64
1100000.00.00.00.0
3329312.8653352580.00.00.0
10076658.3229880820.00.00.0
30498498.358465291.6321922754264707e+301.6321922754264707e+300.0
92308220.87113185.844061207893652e+315.84406111598908e+319.190355615766844e+23
279384497.565880242.659134693097609e+332.659134648524772e+334.456596079142524e+25
845598547.38300931.2820178839927871e+351.2820178674156813e+351.6576972991066742e+27
2559329202.4646046.183855784634422e+366.183855688223801e+369.64097864687752e+28
7746188763.9942262.9621748605150595e+382.9621748058104933e+385.469301780866571e+30
............
1.904888847095696e+257.367890063473198e+400.07.367890063473198e+40
5.765428132426668e+256.120134707524491e+400.06.120134707524491e+40
1.744992186859445e+263.71709390423953e+400.03.71709390423953e+40
5.2814772160188116e+261.3598407428299007e+400.01.3598407428299007e+40
1.598517276660585e+278.760037338641687e-410.08.760037338641687e-41
4.838149213315232e+270.00.00.0
1.4643374927547416e+280.00.00.0
4.432034230746478e+280.00.00.0
1.3414207803664067e+290.00.00.0
4.060008601729061e+290.00.00.0

Of cousrse, since these colums have units, you can easily convert the units of the Synchrotron luminostity form erg/s to GeV/s

t_src['Sync'][::10].to('GeV/s')
\[[0,~0,~0,~1.0187343 \times 10^{33},~3.647576 \times 10^{34},~1.6597013 \times 10^{36},~8.0017262 \times 10^{37},~3.8596591 \times 10^{39},~1.8488441 \times 10^{41},~3.0976555 \times 10^{42},~6.931529 \times 10^{42},~1.2168235 \times 10^{43},~2.1172402 \times 10^{43},~3.6814864 \times 10^{43},~6.3962039 \times 10^{43},~1.1090719 \times 10^{44},~1.9135666 \times 10^{44},~3.263749 \times 10^{44},~5.419136 \times 10^{44},~8.5467044 \times 10^{44},~1.2423887 \times 10^{45},~1.6473596 \times 10^{45},~1.9834424 \times 10^{45},~2.1684948 \times 10^{45},~2.1386317 \times 10^{45},~1.8654372 \times 10^{45},~1.3193169 \times 10^{45},~5.3551014 \times 10^{44},~4.0413955 \times 10^{43},~1.4979971 \times 10^{40},~1.260718 \times 10^{30},~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0,~0] \; \mathrm{\frac{GeV}{s}}\]

the table can be easily saved as an ascii file

t_src.write('test_SED.txt',format='ascii.ecsv',overwrite=True)

or in fits format

t_src.write('test_SED.fits',format='fits',overwrite=True)
WARNING: VerifyWarning: Keyword name 'restframe' is greater than 8 characters or contains characters not allowed by the FITS standard; a HIERARCH card will be created. [astropy.io.fits.card]

Energetic report

It is possible to get an energetic report of the jet model (updated each time that you evaluate the model). This report gives energy densities (U_) (both in the blob end disk restframe), the luminosities of the emitted components in the blob restframe (L_), and the luminosity carried by the jet (jet_L) for the radiative components, the electrons, the magnetic fields, and for the cold protons in the jet.

case of beaming expression ‘bulk_theta’

Important

In this case, the BulkFactor used for the energetic report will be the one taken from the jet parameters, in this case we set BulkFactor=15

my_jet=Jet(name='test',electron_distribution='lppl',beaming_expr='bulk_theta')
my_jet.parameters.BulkFactor.val=15
my_jet
--------------------------------------------------------------------------------
model description:
--------------------------------------------------------------------------------
type: Jet
name: test

electrons distribution:
 type: lppl
 gamma energy grid size:  201
 gmin grid : 2.000000e+00
 gmax grid : 1.000000e+06
 normalization:  True
 log-values:  False
 ratio of cold protons to relativistic electrons: 1.000000e-01

radiative fields:
 seed photons grid size:  100
 IC emission grid size:  100
 source emissivity lower bound :  1.000000e-120
 spectral components:
   name:Sum, state: on
   name:Sync, state: self-abs
   name:SSC, state: on
external fields transformation method: blob

SED info:
 nu grid size jetkernel: 1000
 nu size: 500
 nu mix (Hz): 1.000000e+06
 nu max (Hz): 1.000000e+30

flux plot lower bound   :  1.000000e-30

--------------------------------------------------------------------------------
Table length=13
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testthetajet-viewing-angledeg1.000000e-010.000000e+00--FalseFalse
testBulkFactorjet-bulk-factorlorentz-factor*1.500000e+011.000000e+001.000000e+05FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
--------------------------------------------------------------------------------
None
my_jet.eval()
my_jet.energetic_report()
Table length=34
nametypeunitsval
BulkLorentzFactorjet-bulk-factor1.500000e+01
U_eEnergy dens. blob rest. frameerg / cm31.736635e-03
U_p_coldEnergy dens. blob rest. frameerg / cm31.503276e-02
U_BEnergy dens. blob rest. frameerg / cm33.978874e-04
U_SynchEnergy dens. blob rest. frameerg / cm35.494838e-05
U_Synch_DRFEnergy dens. disk rest. frameerg / cm34.418962e+01
U_DiskEnergy dens. blob rest. frameerg / cm30.000000e+00
U_BLREnergy dens. blob rest. frameerg / cm30.000000e+00
U_DTEnergy dens. blob rest. frameerg / cm30.000000e+00
U_CMBEnergy dens. blob rest. frameerg / cm30.000000e+00
U_Disk_DRFEnergy dens. disk rest. frameerg / cm30.000000e+00
U_BLR_DRFEnergy dens. disk rest. frameerg / cm30.000000e+00
U_DT_DRFEnergy dens. disk rest. frameerg / cm30.000000e+00
U_CMB_DRFEnergy dens. disk rest. frameerg / cm30.000000e+00
L_Sync_rfLum. blob rest. frame.erg / s1.725018e+38
L_SSC_rfLum. blob rest. frame.erg / s3.362216e+36
L_EC_Disk_rfLum. blob rest. frame.erg / s0.000000e+00
L_EC_BLR_rfLum. blob rest. frame.erg / s0.000000e+00
L_EC_DT_rfLum. blob rest. frame.erg / s0.000000e+00
L_EC_CMB_rfLum. blob rest. frame.erg / s0.000000e+00
jet_L_Syncjet Lum.erg / s9.703225e+39
jet_L_SSCjet Lum.erg / s1.891247e+38
jet_L_EC_Diskjet Lum.erg / s0.000000e+00
jet_L_EC_BLRjet Lum.erg / s0.000000e+00
jet_L_EC_DTjet Lum.erg / s0.000000e+00
jet_L_EC_CMBjet Lum.erg / s0.000000e+00
jet_L_pp_gammajet Lum.erg / s0.000000e+00
jet_L_radjet Lum.erg / s9.892349e+39
jet_L_kinjet Lum.erg / s8.864278e+42
jet_L_totjet Lum.erg / s9.084493e+42
jet_L_ejet Lum.erg / s9.179824e+41
jet_L_Bjet Lum.erg / s2.103226e+41
jet_L_p_coldjet Lum.erg / s7.946295e+42
NH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-01

Important

NH_cold_to_rel_e represents the ration of cold protons to relativistic electrons used to compute the cold protons energetic

If you want to evaluate the energetic report in non verbose mode:

my_jet.energetic_report(verbose=False)
my_jet.energetic_dict
{'BulkLorentzFactor': 15.0,
 'U_e': 0.001736634756190472,
 'U_p_cold': 0.015032764261,
 'U_B': 0.00039788735772973844,
 'U_Synch': 5.4948380795748906e-05,
 'U_Synch_DRF': 44.18961605963792,
 'U_Disk': 0.0,
 'U_BLR': 0.0,
 'U_DT': 0.0,
 'U_CMB': 0.0,
 'U_Disk_DRF': 0.0,
 'U_BLR_DRF': 0.0,
 'U_DT_DRF': 0.0,
 'U_CMB_DRF': 0.0,
 'L_Sync_rf': 1.725017708419289e+38,
 'L_SSC_rf': 3.362216361563802e+36,
 'L_EC_Disk_rf': 0.0,
 'L_EC_BLR_rf': 0.0,
 'L_EC_DT_rf': 0.0,
 'L_EC_CMB_rf': 0.0,
 'jet_L_Sync': 9.7032246098585e+39,
 'jet_L_SSC': 1.8912467033796387e+38,
 'jet_L_EC_Disk': 0.0,
 'jet_L_EC_BLR': 0.0,
 'jet_L_EC_DT': 0.0,
 'jet_L_EC_CMB': 0.0,
 'jet_L_pp_gamma': 0.0,
 'jet_L_rad': 9.892349280196463e+39,
 'jet_L_kin': 8.864277658311462e+42,
 'jet_L_tot': 9.084492632274437e+42,
 'jet_L_e': 9.17982370994084e+41,
 'jet_L_B': 2.1032262468277818e+41,
 'jet_L_p_cold': 7.946295287317377e+42,
 'NH_cold_to_rel_e': 0.1}
my_jet.energetic_report_table
Table length=34
nametypeunitsval
str17str29objectfloat64
BulkLorentzFactorjet-bulk-factor1.500000e+01
U_eEnergy dens. blob rest. frameerg / cm31.736635e-03
U_p_coldEnergy dens. blob rest. frameerg / cm31.503276e-02
U_BEnergy dens. blob rest. frameerg / cm33.978874e-04
U_SynchEnergy dens. blob rest. frameerg / cm35.494838e-05
U_Synch_DRFEnergy dens. disk rest. frameerg / cm34.418962e+01
U_DiskEnergy dens. blob rest. frameerg / cm30.000000e+00
U_BLREnergy dens. blob rest. frameerg / cm30.000000e+00
U_DTEnergy dens. blob rest. frameerg / cm30.000000e+00
U_CMBEnergy dens. blob rest. frameerg / cm30.000000e+00
............
jet_L_EC_DTjet Lum.erg / s0.000000e+00
jet_L_EC_CMBjet Lum.erg / s0.000000e+00
jet_L_pp_gammajet Lum.erg / s0.000000e+00
jet_L_radjet Lum.erg / s9.892349e+39
jet_L_kinjet Lum.erg / s8.864278e+42
jet_L_totjet Lum.erg / s9.084493e+42
jet_L_ejet Lum.erg / s9.179824e+41
jet_L_Bjet Lum.erg / s2.103226e+41
jet_L_p_coldjet Lum.erg / s7.946295e+42
NH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-01
my_jet.show_model()
--------------------------------------------------------------------------------
model description:
--------------------------------------------------------------------------------
type: Jet
name: test

electrons distribution:
 type: lppl
 gamma energy grid size:  201
 gmin grid : 2.000000e+00
 gmax grid : 1.000000e+06
 normalization:  True
 log-values:  False
 ratio of cold protons to relativistic electrons: 1.000000e-01

radiative fields:
 seed photons grid size:  100
 IC emission grid size:  100
 source emissivity lower bound :  1.000000e-120
 spectral components:
   name:Sum, state: on
   name:Sync, state: self-abs
   name:SSC, state: on
external fields transformation method: blob

SED info:
 nu grid size jetkernel: 1000
 nu size: 500
 nu mix (Hz): 1.000000e+06
 nu max (Hz): 1.000000e+30

flux plot lower bound   :  1.000000e-30

--------------------------------------------------------------------------------
Table length=13
model namenamepar typeunitsvalphys. bound. minphys. bound. maxlogfrozen
testRregion_sizecm5.000000e+151.000000e+031.000000e+30FalseFalse
testR_Hregion_positioncm1.000000e+170.000000e+00--FalseTrue
testBmagnetic_fieldgauss1.000000e-011.000000e-101.000000e+10FalseFalse
testNH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-010.000000e+00--FalseTrue
testthetajet-viewing-angledeg1.000000e-010.000000e+00--FalseFalse
testBulkFactorjet-bulk-factorlorentz-factor*1.500000e+011.000000e+001.000000e+05FalseFalse
testz_cosmredshift1.000000e-010.000000e+00--FalseFalse
testgminlow-energy-cut-offlorentz-factor*2.000000e+001.000000e+001.000000e+09FalseFalse
testgmaxhigh-energy-cut-offlorentz-factor*1.000000e+061.000000e+001.000000e+15FalseFalse
testNemitters_density1 / cm31.000000e+020.000000e+00--FalseFalse
testgamma0_log_parabturn-over-energylorentz-factor*1.000000e+041.000000e+001.000000e+09FalseFalse
testsLE_spectral_slope2.000000e+00-1.000000e+011.000000e+01FalseFalse
testrspectral_curvature4.000000e-01-1.500000e+011.500000e+01FalseFalse
--------------------------------------------------------------------------------

case of beaming expression ‘delta’

Important

In this case, the BulkFactor used for the energetic report will be the beaming factor parameter beam_obj

my_jet=Jet(name='test',electron_distribution='lppl',beaming_expr='delta')
my_jet.parameters.beam_obj.val=25
my_jet.eval()
my_jet.energetic_report()
Table length=34
nametypeunitsval
BulkLorentzFactorjet-bulk-factor2.500000e+01
U_eEnergy dens. blob rest. frameerg / cm31.736635e-03
U_p_coldEnergy dens. blob rest. frameerg / cm31.503276e-02
U_BEnergy dens. blob rest. frameerg / cm33.978874e-04
U_SynchEnergy dens. blob rest. frameerg / cm35.494838e-05
U_Synch_DRFEnergy dens. disk rest. frameerg / cm32.146421e+01
U_DiskEnergy dens. blob rest. frameerg / cm30.000000e+00
U_BLREnergy dens. blob rest. frameerg / cm30.000000e+00
U_DTEnergy dens. blob rest. frameerg / cm30.000000e+00
U_CMBEnergy dens. blob rest. frameerg / cm30.000000e+00
U_Disk_DRFEnergy dens. disk rest. frameerg / cm30.000000e+00
U_BLR_DRFEnergy dens. disk rest. frameerg / cm30.000000e+00
U_DT_DRFEnergy dens. disk rest. frameerg / cm30.000000e+00
U_CMB_DRFEnergy dens. disk rest. frameerg / cm30.000000e+00
L_Sync_rfLum. blob rest. frame.erg / s1.725018e+38
L_SSC_rfLum. blob rest. frame.erg / s3.362216e+36
L_EC_Disk_rfLum. blob rest. frame.erg / s0.000000e+00
L_EC_BLR_rfLum. blob rest. frame.erg / s0.000000e+00
L_EC_DT_rfLum. blob rest. frame.erg / s0.000000e+00
L_EC_CMB_rfLum. blob rest. frame.erg / s0.000000e+00
jet_L_Syncjet Lum.erg / s2.695340e+40
jet_L_SSCjet Lum.erg / s5.253463e+38
jet_L_EC_Diskjet Lum.erg / s0.000000e+00
jet_L_EC_BLRjet Lum.erg / s0.000000e+00
jet_L_EC_DTjet Lum.erg / s0.000000e+00
jet_L_EC_CMBjet Lum.erg / s0.000000e+00
jet_L_pp_gammajet Lum.erg / s0.000000e+00
jet_L_radjet Lum.erg / s2.747875e+40
jet_L_kinjet Lum.erg / s2.465814e+43
jet_L_totjet Lum.erg / s2.527069e+43
jet_L_ejet Lum.erg / s2.553591e+42
jet_L_Bjet Lum.erg / s5.850635e+41
jet_L_p_coldjet Lum.erg / s2.210455e+43
NH_cold_to_rel_ecold_p_to_rel_e_ratio1.000000e-01