54#define G_CONST "GModelSpectralTable(GEbounds&, GModelSpectralTablePars&, "\
56#define G_FLUX "GModelSpectralTable::flux(GEnergy&, GEnergy&)"
57#define G_EFLUX "GModelSpectralTable::eflux(GEnergy&, GEnergy&)"
58#define G_MC "GModelSpectralTable::mc(GEnergy&, GEnergy&, GTime&, GRan&)"
59#define G_READ "GModelSpectralTable::read(GXmlElement&)"
60#define G_WRITE "GModelSpectralTable::write(GXmlElement&)"
61#define G_LOAD "GModelSpectralTable::load(GFilename&)"
62#define G_TABLE_PAR "GModelSpectralTable::table_par(int&)"
63#define G_LOAD_PAR "GModelSpectralTable::load_par(GFits&)"
64#define G_PAR_INDEX "GModelSpectralTable::par_index(std::string&)"
65#define G_UPDATE "GModelSpectralTable::update()"
105 const double&
norm) :
145 if (spectra.
dim() != pars.
size()+1) {
147 " dimensions but an array with "+
149 "expected. Please specify a spectra array with the "
150 "correct dimension.";
155 for (
int i = 0; i < pars.
size(); ++i) {
156 int npars = pars[i]->
size();
157 int nspec = spectra.
shape()[i];
158 if (npars != nspec) {
159 std::string msg =
"Parameter \""+pars[i]->par().name()+
"\" has "+
163 "array with the correct number of spectra.";
170 int nspec = spectra.
shape()[pars.
size()];
171 if (nebins != nspec) {
172 std::string msg =
"Spectra array has "+
gammalib::str(nspec)+
" energy "
174 " energy boundaries. Please specify a spectra "
175 "array with the correct number of energy bins.";
265 if (
this != &model) {
375 const GTime& srcTime,
376 const bool& gradients)
const
404 for (
int i = 0; i < dim; ++i) {
428 #if defined(G_NAN_CHECK)
430 std::cout <<
"*** ERROR: GModelSpectralTable::eval";
431 std::cout <<
"(srcEng=" << srcEng;
432 std::cout <<
", srcTime=" << srcTime <<
"):";
433 std::cout <<
" NaN/Inf encountered";
434 std::cout <<
" (value=" << value;
435 std::cout <<
", norm=" <<
norm();
436 std::cout <<
", func=" << func;
437 std::cout <<
")" << std::endl;
477 double e_min = emin.
MeV();
478 double e_max = emax.
MeV();
491 if (inx_emin == inx_emax) {
506 int i_start = (e_min <
m_lin_nodes[0]) ? inx_emin : inx_emin+1;
516 for (
int i = i_start; i < inx_emax; ++i) {
570 double e_min = emin.
MeV();
571 double e_max = emax.
MeV();
584 if (inx_emin == inx_emax) {
600 int i_start = (e_min <
m_lin_nodes[0]) ? inx_emin : inx_emin+1;
611 for (
int i = i_start; i < inx_emax; ++i) {
667 for (inx =
m_mc_cum.size()-1; inx > 0; --inx) {
679 double eng = (u > 0.0)
680 ? std::exp(std::log(u * (e_max - e_min) + e_min) /
m_mc_exp[inx])
688 double eng = std::exp(u * (e_max - e_min) + e_min);
841 const bool& clobber)
const
858 primary->
card(
"CONTENT",
"MODEL",
"Spectrum file");
860 primary->
card(
"ORIGIN", PACKAGE_NAME,
"Origin of FITS file");
861 primary->
card(
"MODLNAME",
"model",
"Model name");
862 primary->
card(
"MODLUNIT",
"photons/cm^2/s/MeV",
"Model units");
863 primary->
card(
"REDSHIFT",
false,
"If true then redshift will be included as a par");
864 primary->
card(
"ADDMODEL",
true,
"If true then this is an additive table model");
865 primary->
card(
"HDUCLASS",
"OGIP",
"Format conforms to OGIP standard");
866 primary->
card(
"HDUCLAS1",
"XSPEC TABLE MODEL",
"Model spectra for XSPEC");
867 primary->
card(
"HDUVERS",
"1.0.0",
"Version of format");
889 if (index < 0 || index >=
size()) {
908 if (index < 0 || index >=
size()) {
965 result.append(
"=== GModelSpectralTable ===");
972 for (
int i = 0; i <
size(); ++i) {
985 for (
int k = 0; k <
size; ++k) {
997 result.append(
" (linear)");
1000 result.append(
" (logarithmic)");
1007 result.append(
" [");
1009 result.append(
", ");
1194 for (
int i = 0; i < nebins; ++i) {
1196 double log10_energy_MeV = std::log10(energy_MeV);
1223 for (
int i = 0; i < nrows; ++i) {
1243 for (
int i = 0; i < nrows; ++i) {
1249 col_name(i) = par.
name();
1251 col_initial(i) = (float)par.
value();
1255 col_delta(i) = -1.0;
1265 double min_value = 0.0;
1266 double max_value = 0.0;
1272 if (value < min_value) {
1275 if (value > max_value) {
1278 col_value(i,k) = value;
1284 if (par.
min() < min_value) {
1285 col_minimum(i) = (float)par.
min();
1286 col_bottom(i) = min_value;
1289 col_minimum(i) = (float)par.
min();
1290 col_bottom(i) = (float)par.
min();
1294 col_minimum(i) = min_value;
1295 col_bottom(i) = min_value;
1300 if (par.
max() > max_value) {
1301 col_top(i) = max_value;
1302 col_maximum(i) = (float)par.
max();
1305 col_top(i) = (float)par.
max();
1306 col_maximum(i) = (float)par.
max();
1310 col_top(i) = max_value;
1311 col_maximum(i) = max_value;
1321 table.
append(col_method);
1322 table.
append(col_initial);
1324 table.
append(col_minimum);
1325 table.
append(col_bottom);
1327 table.
append(col_maximum);
1328 table.
append(col_numbvals);
1335 table.
card(
"HDUCLASS",
"OGIP",
"Format conforms to OGIP standard");
1336 table.
card(
"HDUCLAS1",
"XSPEC TABLE MODEL",
"Model spectra for XSPEC");
1337 table.
card(
"HDUCLAS2",
"PARAMETERS",
"Extension containing parameter info");
1338 table.
card(
"HDUVERS",
"1.0.0",
"Version of format");
1339 table.
card(
"NINTPARM", nrows,
"Number of interpolation parameters");
1340 table.
card(
"NADDPARM", 0,
"Number of additional parameters");
1382 table.
card(
"HDUCLASS",
"OGIP",
"Format conforms to OGIP standard");
1383 table.
card(
"HDUCLAS1",
"XSPEC TABLE MODEL",
"Model spectra for XSPEC");
1384 table.
card(
"HDUCLAS2",
"ENERGIES",
"Extension containing energy bin info");
1385 table.
card(
"HDUVERS",
"1.0.0",
"Version of format");
1419 col_spec.
unit(
"ph cm-2 s-1 MeV-1");
1422 std::vector<int> inx(npars+1,0);
1423 for (
int i = 0; i < nrows; ++i) {
1426 for (
int k = 0; k < npars; ++k) {
1431 std::vector<int> index = inx;
1434 for (
int k = 0; k < nebins; ++k, ++index[npars]) {
1449 }
while (ipar >= 0);
1464 table.
card(
"HDUCLASS",
"OGIP",
"Format conforms to OGIP standard");
1465 table.
card(
"HDUCLAS1",
"XSPEC TABLE MODEL",
"Model spectra for XSPEC");
1466 table.
card(
"HDUCLAS2",
"MODEL SPECTRA",
"Extension containing model spectra");
1467 table.
card(
"HDUVERS",
"1.0.0",
"Version of format");
1493 int npars = table.
integer(
"NAXIS2");
1496 for (
int i = 0; i < npars; ++i) {
1499 GModelPar par(table[
"NAME"]->
string(i), table[
"INITIAL"]->real(i));
1516 if (table[
"DELTA"]->real(i) < 0.0) {
1526 std::vector<double> values;
1527 for (
int k = 0; k < table[
"NUMBVALS"]->
integer(i); ++k) {
1530 double value = table[
"VALUE"]->
real(i,k);
1551 values.push_back(value);
1559 table_model_par.
method(table[
"METHOD"]->integer(i));
1589 int nebins = table.
integer(
"NAXIS2");
1595 std::string emin_unit = table[
"ENERG_LO"]->unit();
1596 std::string emax_unit = table[
"ENERG_HI"]->unit();
1597 if (emin_unit.empty()) {
1600 if (emax_unit.empty()) {
1605 for (
int i = 0; i < nebins; ++i) {
1606 GEnergy emin(table[
"ENERG_LO"]->real(i), emin_unit);
1607 GEnergy emax(table[
"ENERG_HI"]->real(i), emax_unit);
1630 int npars = fits.
table(
"PARAMETERS")->
integer(
"NAXIS2");
1631 int nebins = fits.
table(
"ENERGIES")->
integer(
"NAXIS2");
1634 std::vector<int> naxis(npars+1, 0);
1635 for (
int i = 0; i < npars; ++i) {
1636 naxis[i] = (*fits.
table(
"PARAMETERS"))[
"NUMBVALS"]->integer(i);
1638 naxis[npars] = nebins;
1647 int nrows = table.
integer(
"NAXIS2");
1650 for (
int i = 0; i < nrows; ++i) {
1653 std::vector<int> index(npars+1, 0);
1656 std::vector<double> parval(npars, 0.0);
1657 for (
int k = 0; k < npars; ++k) {
1663 nodes.
set_value(table[
"PARAMVAL"]->real(i,k));
1675 #if defined(G_DEBUG_LOAD_SPEC)
1676 std::cout << i <<
": (";
1677 for (
int k = 0; k < index.size(); ++k) {
1681 std::cout << index[k];
1683 std::cout <<
")" << std::endl;
1687 for (
int k = 0; k < nebins; ++k, ++index[npars]) {
1711 for (; index <
size(); ++index) {
1718 if (index >=
size()) {
1719 std::string msg =
"Parameter name \""+name+
"\" not found in spectral "
1720 "table. Please specify one of the following parameter "
1722 for (
int i = 0; i <
size(); ++i) {
1801 bool need_update =
false;
1812 for (
int i = 0; i <
m_npars; ++i) {
1824 #if defined(G_DEBUG_UPDATE)
1825 std::cout <<
"GModelSpectralTable::update() required" << std::endl;
1829 std::vector<double> weights(2*
m_npars, 0.0);
1830 std::vector<double> weight_gradients(2*
m_npars, 0.0);
1831 std::vector<int> indices(2*
m_npars, 0);
1835 for (
int i = 0; i <
m_npars; ++i) {
1842 double value = par->
value();
1878 #if defined(G_DEBUG_UPDATE)
1879 std::cout <<
" wgt_l=" << weights[il];
1880 std::cout <<
" wgt_r=" << weights[ir];
1881 std::cout <<
" wgt_grad_l=" << weight_gradients[il];
1882 std::cout <<
" wgt_grad_r=" << weight_gradients[ir];
1883 std::cout <<
" inx_l=" << indices[il];
1884 std::cout <<
" inx_r=" << indices[ir] << std::endl;
1890 int combinations = 1 <<
m_npars;
1897 #if defined(G_DEBUG_UPDATE)
1898 double weight_sum = 0.0;
1902 for (
int i = 0; i < combinations; ++i) {
1905 #if defined(G_DEBUG_UPDATE)
1906 std::cout <<
" " << i <<
": ";
1910 double weight = 1.0;
1911 std::vector<double> grad_weight(
m_npars, 1.0);
1914 std::vector<int> index_shape(
m_npars+1,0);
1917 for (
int k = 0, div = 1; k <
m_npars; ++k, div *= 2) {
1920 int index = i/div % 2 + k * 2;
1923 weight *= weights[index];
1926 index_shape[k] = indices[index];
1929 for (
int j = 0; j <
m_npars; ++j) {
1931 grad_weight[j] *= weight_gradients[index];
1934 grad_weight[j] *= weights[index];
1939 #if defined(G_DEBUG_UPDATE)
1941 std::cout <<
" (" << weights[index];
1942 std::cout <<
" @ " << indices[index] <<
") ";
1948 #if defined(G_DEBUG_UPDATE)
1949 std::cout <<
": wgt=" << weight;
1951 for (
int k = 0; k <
m_npars; ++k) {
1955 std::cout << grad_weight[k];
1957 std::cout <<
"]" << std::endl;
1958 weight_sum += weight;
1962 for (
int iebin = 0; iebin <
m_nebins; ++iebin) {
1974 for (
int j = 0; j <
m_npars; ++j) {
1995 #if defined(G_DEBUG_UPDATE)
1996 std::cout <<
" sum(wgt)=" << weight_sum << std::endl;
2029 double epivot = std::sqrt(emin*emax);
2032 double gamma = std::log(fmin/fmax) / std::log(emin/emax);
2035 double prefactor = fmin / std::pow(emin/epivot, gamma);
2038 double flux = prefactor *
2042 double eflux = prefactor *
2088 double e_min = emin.
MeV();
2089 double e_max = emax.
MeV();
2092 if (e_max > e_min) {
2107 if (inx_emin == inx_emax) {
2126 int i_start = (e_min <
m_lin_nodes[0]) ? inx_emin : inx_emin+1;
2140 for (
int i = i_start; i < inx_emax; ++i) {
2162 for (
int i = 1; i <
m_mc_cum.size(); ++i) {
2167 for (
int i = 0; i <
m_mc_cum.size(); ++i) {
2173 for (
int i = 0; i <
m_mc_cum.size(); ++i) {
2176 double exponent =
m_mc_exp[i] + 1.0;
2179 if (std::abs(exponent) > 1.0e-11) {
2183 if (exponent < -50.0) {
2191 else if (exponent > +50.0) {
Energy value class definition.
Exception handler interface definition.
FITS binary table class definition.
FITS table float column class interface definition.
FITS table long integer column class interface definition.
FITS table string column class interface definition.
FITS file class interface definition.
Mathematical function definitions.
Spectral model registry class definition.
Spectral table model parameter class definition.
Spectral table model parameter container class definition.
const GModelSpectralTable g_spectral_table_seed
Spectral table model class definition.
Random number generator class definition.
double norm(const GVector &vector)
Computes vector norm.
double min(const GVector &vector)
Computes minimum vector element.
double max(const GVector &vector)
Computes maximum vector element.
Energy boundaries container class.
const GEnergy & emax(void) const
Return maximum energy of all intervals.
void append(const GEnergy &emin, const GEnergy &emax)
Append energy interval.
int size(void) const
Return number of energy boundaries.
void clear(void)
Clear energy boundaries.
GEnergy elogmean(const int &index) const
Returns logarithmic mean energy for a given energy interval.
const GEnergy & emin(void) const
Return minimum energy of all intervals.
Class that handles energies in a unit independent way.
double keV(void) const
Return energy in keV.
double MeV(void) const
Return energy in MeV.
double log10MeV(void) const
Return log10 of energy in MeV.
std::string print(const GChatter &chatter=NORMAL) const
Print energy.
void clear(void)
Clear instance.
std::string url(void) const
Return Uniform Resource Locator (URL)
bool is_empty(void) const
Signal if filename is empty.
void clear(void)
Clear file name.
Abstract FITS extension base class.
const std::string & extname(void) const
Return extension name.
double real(const std::string &keyname) const
Return card value as double precision.
GFitsHeaderCard & card(const int &cardno)
Return header card.
int integer(const std::string &keyname) const
Return card value as integer.
void unit(const std::string &unit)
Set column unit.
FITS table long integer column.
FITS table string column.
Abstract interface for FITS table.
GFitsTableCol * append(const GFitsTableCol &column)
Append column to the table.
void saveto(const GFilename &filename, const bool &clobber=false)
Saves to specified FITS file.
GFitsHDU * append(const GFitsHDU &hdu)
Append HDU to FITS file.
void close(void)
Close FITS file.
GFitsTable * table(const int &extno)
Get pointer to table HDU.
void write(GXmlElement &xml) const
Set or update parameter attributes in XML element.
void read(const GXmlElement &xml)
Extract parameter attributes from XML element.
Interface definition for the spectral model registry class.
Spectral table model parameter class.
const int & method(void) const
Return reference to table model parameter interpolation method.
Spectral table model parameter container class.
GModelSpectralTablePar * append(const GModelSpectralTablePar &par)
Append table model parameter to container.
void clear(void)
Clear table model parameters.
int size(void) const
Return number of table model parameters in container.
Spectral table model class.
GModelSpectralTable(void)
Void constructor.
virtual void write(GXmlElement &xml) const
Write model into XML element.
std::vector< double > m_mc_cum
Cumulative distribution.
virtual double eval(const GEnergy &srcEng, const GTime &srcTime=GTime(), const bool &gradients=false) const
Evaluate spectral table model.
virtual void read(const GXmlElement &xml)
Read model from XML element.
std::vector< double > m_mc_min
Lower boundary for MC.
void save(const GFilename &filename, const bool &clobber=false) const
Save table into file.
GModelPar m_norm
Normalization factor.
void load_par(const GFits &fits)
Load data from PARAMETERS extension.
void load_spec(const GFits &fits)
Load data from SPECTRA extension.
GModelSpectralTablePar & table_par(const int &index)
Return reference to table parameter.
virtual std::string print(const GChatter &chatter=NORMAL) const
Print table model information.
GEnergy m_mc_emax
Maximum energy.
GEbounds m_ebounds
Energy boundaries.
GNodeArray m_log_nodes
log10(Energy) nodes of function
virtual GEnergy mc(const GEnergy &emin, const GEnergy &emax, const GTime &time, GRan &ran) const
Returns MC energy between [emin, emax].
std::vector< double > m_epivot
Power-law pivot energies.
virtual std::string type(void) const
Return model type.
void init_members(void)
Initialise class members.
GFilename m_filename
Filename of table.
void set_par_pointers(void)
Set parameter pointers.
int par_index(const std::string &name) const
Return index for parameter name.
void set_energy_nodes(void)
Set energy nodes from energy boundaries.
void free_members(void)
Delete class members.
std::vector< double > m_mc_exp
Exponent for MC.
virtual double eflux(const GEnergy &emin, const GEnergy &emax) const
Returns model energy flux between [emin, emax] (units: erg/cm2/s)
int m_nebins
Number of energy bins.
GNdarray m_log_values
log10(Function) values and grad's
GEnergy m_mc_emin
Minimum energy.
virtual ~GModelSpectralTable(void)
Destructor.
double norm(void) const
Return normalization factor.
GNodeArray m_lin_nodes
Energy nodes of function.
std::vector< double > m_flux
Photon fluxes.
void update_mc(const GEnergy &emin, const GEnergy &emax) const
Update MC cache.
void update(void) const
Update cache for spectral table model computation.
GNdarray m_lin_values
Function values and grad's.
GFitsBinTable create_eng_table(void) const
Create ENERGIES FITS table.
const GEbounds & ebounds(void) const
Return reference to energy boundaries.
virtual void clear(void)
Clear table model.
virtual GModelSpectralTable * clone(void) const
Clone table model.
void load(const GFilename &filename)
Load table from file.
std::vector< double > m_mc_max
Upper boundary for MC.
std::vector< double > m_last_values
Last parameter values.
void copy_members(const GModelSpectralTable &model)
Copy class members.
GFitsBinTable create_par_table(void) const
Create PARAMETERS FITS table.
void update_flux(void) const
Update flux cache.
std::vector< double > m_prefactor
Power-law normalisations.
int m_npars
Number of parameters.
void load_eng(const GFits &fits)
Load data from ENERGIES extension.
virtual GModelSpectralTable & operator=(const GModelSpectralTable &model)
Assignment operator.
std::vector< double > m_gamma
Power-law indices.
GModelSpectralTablePars m_table_pars
Table model parameters.
GNdarray m_spectra
Spectra.
virtual double flux(const GEnergy &emin, const GEnergy &emax) const
Returns model photon flux between [emin, emax] (units: ph/cm2/s)
const GFilename & filename(void) const
Return file name.
GFitsBinTable create_spec_table(void) const
Create SPECTRA FITS table.
std::vector< double > m_eflux
Energy fluxes.
Abstract spectral model base class.
void free_members(void)
Delete class members.
virtual GModelSpectral & operator=(const GModelSpectral &model)
Assignment operator.
std::vector< GModelPar * > m_pars
Parameter pointers.
int size(void) const
Return number of parameters.
void init_members(void)
Initialise class members.
N-dimensional array class.
void clear(void)
Clear array.
int dim(void) const
Return dimension of array.
const std::vector< int > & shape(void) const
Return shape of array.
void set_value(const double &value) const
Set indices and weighting factors for interpolation.
const double & wgt_grad_left(void) const
Returns left node weight gradient.
const int & inx_right(void) const
Returns right node index.
const int & inx_left(void) const
Returns left node index.
const double & wgt_grad_right(void) const
Returns right node weight gradient.
const double & wgt_right(void) const
Returns right node weight.
const double & wgt_left(void) const
Returns left node weight.
void clear(void)
Clear node array.
int size(void) const
Return number of nodes in node array.
void append(const double &node)
Append one node to array.
bool is_free(void) const
Signal if parameter is free.
void free(void)
Free a parameter.
const double & scale(void) const
Return parameter scale.
bool has_grad(void) const
Signal if parameter gradient is computed analytically.
void range(const double &min, const double &max)
Set minimum and maximum parameter boundaries.
bool has_min(void) const
Signal if parameter has minimum boundary.
double max(void) const
Return parameter maximum boundary.
bool is_fixed(void) const
Signal if parameter is fixed.
double min(void) const
Return parameter minimum boundary.
void remove_range(void)
Removes minimum and maximum boundary.
const double & factor_gradient(void) const
Return parameter factor gradient.
void fix(void)
Fix a parameter.
bool has_max(void) const
Signal if parameter has maximum boundary.
void factor_range(const double &min, const double &max)
Set minimum and maximum parameter boundary factors.
double gradient(void) const
Return parameter gradient.
void clear(void)
Clear parameter.
double value(void) const
Return parameter value.
const std::string & name(void) const
Return parameter name.
Random number generator class.
double uniform(void)
Returns random double precision floating value in range 0 to 1.
const GXmlAttribute * attribute(const int &index) const
Return attribute.
void check_energy_interval(const std::string &origin, const GEnergy &emin, const GEnergy &emax)
Checks energy interval.
std::string parformat(const std::string &s, const int &indent=0)
Convert string in parameter format.
bool is_infinite(const double &x)
Signal if argument is infinite.
double plaw_photon_flux(const double &emin, const double &emax, const double &epivot, const double &gamma)
Compute photon flux between two energies for a power law.
bool is_notanumber(const double &x)
Signal if argument is not a number.
std::string str(const unsigned short int &value)
Convert unsigned short integer value into string.
const GXmlElement * xml_get_par(const std::string &origin, const GXmlElement &xml, const std::string &name)
Return pointer to parameter with given name in XML element.
GXmlElement * xml_need_par(const std::string &origin, GXmlElement &xml, const std::string &name)
Return pointer to parameter with given name in XML element.
GFilename xml_file_reduce(const GXmlElement &xml, const std::string &filename)
Reduce file name provided for writing as XML attribute.
double plaw_energy_flux(const double &emin, const double &emax, const double &epivot, const double &gamma)
Compute energy flux between two energies for a power law.
bool xml_has_par(const GXmlElement &xml, const std::string &name)
Checks if parameter with given name in XML element exists.
GFilename xml_file_expand(const GXmlElement &xml, const std::string &filename)
Expand file name provided as XML attribute for loading.
void xml_check_type(const std::string &origin, GXmlElement &xml, const std::string &type)
Checks the model type.