GammaLib  2.1.0.dev
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GCTAModelSpatialLookup Class Reference

Spatial lookup table model class. More...

#include <GCTAModelSpatialLookup.hpp>

Inheritance diagram for GCTAModelSpatialLookup:
GCTAModelSpatial GBase

Public Member Functions

 GCTAModelSpatialLookup (void)
 Void constructor. More...
 
 GCTAModelSpatialLookup (const GFilename &filename)
 Lookup table file constructor. More...
 
 GCTAModelSpatialLookup (const GCTAResponseTable &table)
 Response table constructor. More...
 
 GCTAModelSpatialLookup (const GXmlElement &xml)
 XML constructor. More...
 
 GCTAModelSpatialLookup (const double &maxrad, const double &radbin, const GEbounds &ebds)
 Lookup table constructor. More...
 
 GCTAModelSpatialLookup (const GCTAModelSpatialLookup &model)
 Copy constructor. More...
 
virtual ~GCTAModelSpatialLookup (void)
 Destructor. More...
 
virtual GCTAModelSpatialLookupoperator= (const GCTAModelSpatialLookup &model)
 Assignment operator. More...
 
virtual void clear (void)
 Clear instance. More...
 
virtual GCTAModelSpatialLookupclone (void) const
 Clone instance. More...
 
virtual std::string classname (void) const
 Return class name. More...
 
virtual std::string type (void) const
 Return model type. More...
 
virtual double eval (const GCTAInstDir &dir, const GEnergy &energy, const GTime &time, const bool &gradients=false) const
 Evaluate function. More...
 
virtual double mc_max_value (const GCTAObservation &obs) const
 Return maximum function value for Monte Carlo simulations. More...
 
virtual void read (const GXmlElement &xml)
 Read model from XML element. More...
 
virtual void write (GXmlElement &xml) const
 Write model into XML element. More...
 
virtual std::string print (const GChatter &chatter=NORMAL) const
 Print lookup table information. More...
 
void fill (const GCTAObservation &obs)
 Fill lookup table with events from one CTA observation. More...
 
void fill (const GObservations &obs)
 Fill lookup table with CTA events in observation container. More...
 
const GCTAResponseTabletable (void) const
 Return lookup table. More...
 
void table (const GCTAResponseTable &table)
 Assign lookup table. More...
 
void read (const GFitsTable &table)
 Read lookup table from FITS table. More...
 
void write (GFitsBinTable &table) const
 Write lookup table into FITS binary table. More...
 
void load (const GFilename &filename)
 Load lookup table. More...
 
void save (const GFilename &filename, const bool &clobber=false) const
 Save lookup table into FITS file. More...
 
double norm (void) const
 Get lookup table model normalisation. More...
 
void norm (const double &norm)
 Set lookup table model normalisation. More...
 
- Public Member Functions inherited from GCTAModelSpatial
 GCTAModelSpatial (void)
 Void constructor. More...
 
 GCTAModelSpatial (const GCTAModelSpatial &model)
 Copy constructor. More...
 
virtual ~GCTAModelSpatial (void)
 Destructor. More...
 
virtual GCTAModelSpatialoperator= (const GCTAModelSpatial &model)
 Assignment operator. More...
 
virtual GModelParoperator[] (const int &index)
 Returns model parameter. More...
 
virtual const GModelParoperator[] (const int &index) const
 Returns model parameter (const version) More...
 
virtual GModelParoperator[] (const std::string &name)
 Returns reference to model parameter. More...
 
virtual const GModelParoperator[] (const std::string &name) const
 Returns reference to model parameter (const version) More...
 
virtual GCTAInstDir mc (const GEnergy &energy, const GTime &time, const GCTAObservation &obs, GRan &ran) const
 Returns MC instrument direction. More...
 
int size (void) const
 Return number of model parameters. More...
 
virtual double npred (const GEnergy &energy, const GTime &time, const GObservation &obs) const
 Return integral of spatial model component. More...
 
- Public Member Functions inherited from GBase
virtual ~GBase (void)
 Destructor. More...
 

Protected Member Functions

void init_members (void)
 Initialise class members. More...
 
void copy_members (const GCTAModelSpatialLookup &model)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void prepare_table (void)
 Prepare lookup table indices. More...
 
void normalise_table (void)
 Normalise lookup table. More...
 
int table_index (const int &ienergy, const int &itheta) const
 Return index of lookup table element. More...
 
void fill_buffer (const GCTAObservation &obs, std::vector< double > &buffer)
 Fill buffer from events in CTA observation. More...
 
void set_from_buffer (const std::vector< double > &buffer)
 Set lookup table from buffer. More...
 
- Protected Member Functions inherited from GCTAModelSpatial
void init_members (void)
 Initialise class members. More...
 
void copy_members (const GCTAModelSpatial &model)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 

Protected Attributes

GFilename m_filename
 Name of lookup table. More...
 
GCTAResponseTable m_lookup
 Lookup table. More...
 
GModelPar m_norm
 Normalization factor. More...
 
int m_inx_energy
 Energy index. More...
 
int m_inx_theta
 Theta index. More...
 
- Protected Attributes inherited from GCTAModelSpatial
std::vector< GModelPar * > m_pars
 Parameter pointers. More...
 

Detailed Description

Spatial lookup table model class.

Definition at line 58 of file GCTAModelSpatialLookup.hpp.

Constructor & Destructor Documentation

GCTAModelSpatialLookup::GCTAModelSpatialLookup ( void  )

Void constructor.

Definition at line 79 of file GCTAModelSpatialLookup.cpp.

References init_members().

Referenced by clone().

GCTAModelSpatialLookup::GCTAModelSpatialLookup ( const GFilename filename)
explicit

Lookup table file constructor.

Parameters
[in]filenameLookup table file name.

Creates instance of a spatial lookup table model from a lookup table FITS file.

Definition at line 97 of file GCTAModelSpatialLookup.cpp.

References init_members(), and load().

GCTAModelSpatialLookup::GCTAModelSpatialLookup ( const GCTAResponseTable table)
explicit

Response table constructor.

Parameters
[in]tableResponse table.

Creates instance of a spatial lookup table model from a response table.

Definition at line 118 of file GCTAModelSpatialLookup.cpp.

References init_members(), and table().

GCTAModelSpatialLookup::GCTAModelSpatialLookup ( const GXmlElement xml)
explicit

XML constructor.

Parameters
[in]xmlXML element.

Creates instance of a spatial lookup table model from an XML element. See the read() method for information about the expected structure of the XML element.

Definition at line 141 of file GCTAModelSpatialLookup.cpp.

References init_members(), and read().

GCTAModelSpatialLookup::GCTAModelSpatialLookup ( const double &  maxrad,
const double &  radbin,
const GEbounds ebds 
)

Lookup table constructor.

Parameters
[in]maxradMaximum radius (deg).
[in]radbinRadial binsize (deg).
[in]ebdsEnergy boundaries.
Exceptions
GException::invalid_argumentNon positive radial binsize specified.

Creates empty instance of a spatial lookup table model using the maximum lookup radius maxrad, the radial lookup binsize radbin, and some energy boundaries ebds. Following construction, the instance may be filled using the fill() methods.

The method throws an exception if a non-positive radial bin size is specified.

Definition at line 173 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::append_axis(), GCTAResponseTable::append_table(), GEbounds::emax(), GEbounds::emin(), G_CONSTRUCTOR, init_members(), m_lookup, prepare_table(), GEbounds::size(), gammalib::str(), and GEnergy::TeV().

GCTAModelSpatialLookup::GCTAModelSpatialLookup ( const GCTAModelSpatialLookup model)

Copy constructor.

Parameters
[in]modelLookup table spatial model.

Definition at line 223 of file GCTAModelSpatialLookup.cpp.

References copy_members(), and init_members().

GCTAModelSpatialLookup::~GCTAModelSpatialLookup ( void  )
virtual

Destructor.

Definition at line 240 of file GCTAModelSpatialLookup.cpp.

References free_members().

Member Function Documentation

std::string GCTAModelSpatialLookup::classname ( void  ) const
inlinevirtual

Return class name.

Returns
String containing the class name ("GCTAModelSpatialLookup").

Implements GCTAModelSpatial.

Definition at line 128 of file GCTAModelSpatialLookup.hpp.

void GCTAModelSpatialLookup::clear ( void  )
virtual

Clear instance.

Implements GCTAModelSpatial.

Definition at line 294 of file GCTAModelSpatialLookup.cpp.

References GCTAModelSpatial::free_members(), free_members(), GCTAModelSpatial::init_members(), and init_members().

Referenced by load(), and read().

GCTAModelSpatialLookup * GCTAModelSpatialLookup::clone ( void  ) const
virtual

Clone instance.

Implements GCTAModelSpatial.

Definition at line 312 of file GCTAModelSpatialLookup.cpp.

References GCTAModelSpatialLookup().

void GCTAModelSpatialLookup::copy_members ( const GCTAModelSpatialLookup model)
protected

Copy class members.

Parameters
[in]modelLookup table model.

Definition at line 836 of file GCTAModelSpatialLookup.cpp.

References m_filename, m_inx_energy, m_inx_theta, m_lookup, m_norm, and GCTAModelSpatial::m_pars.

Referenced by GCTAModelSpatialLookup(), and operator=().

double GCTAModelSpatialLookup::eval ( const GCTAInstDir dir,
const GEnergy energy,
const GTime time,
const bool &  gradients = false 
) const
virtual

Evaluate function.

Parameters
[in]dirEvent direction.
[in]energyEvent energy.
[in]timeEvent time.
[in]gradientsCompute gradients?
Returns
Function value

Evaluates the lookup table model for a given event direction and energy. The evaluation is done using bi-linear interpolation in the logarithm of the energy and the offset angle radius. The lookup table value is multiplied with a normalisation parameter. The method always returns a non-negative result.

If the gradients flag is true the method will also compute the partial derivatives of the normalisation parameter.

Implements GCTAModelSpatial.

Definition at line 336 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::elements(), GOptimizerPar::factor_gradient(), GOptimizerPar::is_free(), GEnergy::log10TeV(), m_lookup, m_norm, GOptimizerPar::scale(), GCTAResponseTable::tables(), GCTAInstDir::theta(), and GOptimizerPar::value().

void GCTAModelSpatialLookup::fill ( const GCTAObservation obs)

Fill lookup table with events from one CTA observation.

Parameters
[in]obsCTA observation.
Exceptions
GException::invalid_valueNo lookup table has been defined.

Fills lookup table with events from one CTA observation. The method needs that a lookup table is defined before. In case that no lookup table is defined the method throws an exception.

Definition at line 466 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::elements(), fill_buffer(), G_FILL_1, m_lookup, set_from_buffer(), and GCTAResponseTable::tables().

void GCTAModelSpatialLookup::fill ( const GObservations obs)

Fill lookup table with CTA events in observation container.

Parameters
[in]obsObservation container.
Exceptions
GException::invalid_valueNo lookup table has been defined.

Fills lookup table with CTA events in observation container. The method needs that a lookup table is defined before. In case that no lookup table is defined the method throws an exception.

Definition at line 509 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::elements(), GCTAObservation::eventtype(), fill_buffer(), G_FILL_2, m_lookup, set_from_buffer(), GObservations::size(), and GCTAResponseTable::tables().

void GCTAModelSpatialLookup::fill_buffer ( const GCTAObservation obs,
std::vector< double > &  buffer 
)
protected

Fill buffer from events in CTA observation.

Parameters
[in]obsCTA observation.
[in,out]bufferCounts buffer.

Definition at line 968 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::axis_bins(), GCTAResponseTable::axis_hi(), GCTAResponseTable::axis_lo(), GObservation::events(), G_FILL_BUFFER, m_inx_energy, m_inx_theta, m_lookup, gammalib::rad2deg, GCTAEventList::size(), table_index(), and GCTAInstDir::theta().

Referenced by fill().

void GCTAModelSpatialLookup::free_members ( void  )
protected

Delete class members.

Definition at line 857 of file GCTAModelSpatialLookup.cpp.

Referenced by clear(), operator=(), and ~GCTAModelSpatialLookup().

void GCTAModelSpatialLookup::load ( const GFilename filename)

Load lookup table.

Parameters
[in]filenameLoad lookup table file name.

Loads lookup table from FITS file.

Definition at line 646 of file GCTAModelSpatialLookup.cpp.

References clear(), GFits::close(), GFilename::extname(), GFilename::extno(), G_LOAD, GFilename::has_extname(), GFilename::has_extno(), GFilename::is_empty(), m_filename, read(), GFits::table(), table(), and GFilename::url().

Referenced by GCTAModelSpatialLookup(), and read().

double GCTAModelSpatialLookup::mc_max_value ( const GCTAObservation obs) const
inlinevirtual

Return maximum function value for Monte Carlo simulations.

Parameters
[in]obsCTA Observation.
Returns
Maximum function value for Monte Carlo simulations.

This method always returns 1.

Implements GCTAModelSpatial.

Definition at line 155 of file GCTAModelSpatialLookup.hpp.

double GCTAModelSpatialLookup::norm ( void  ) const
inline

Get lookup table model normalisation.

Returns
Lookup table model normalisation.

Returns the normalisation of the lookup table model.

Definition at line 183 of file GCTAModelSpatialLookup.hpp.

References m_norm, and GOptimizerPar::value().

Referenced by read(), and write().

void GCTAModelSpatialLookup::norm ( const double &  norm)
inline

Set lookup table model normalisation.

Parameters
[in]normLookup table model normalisation.

Set the normalisation of the lookup table model.

Definition at line 197 of file GCTAModelSpatialLookup.hpp.

References m_norm, and GOptimizerPar::value().

void GCTAModelSpatialLookup::normalise_table ( void  )
protected

Normalise lookup table.

Normalises lookup table so that each radial profile has a maximum value of one.

Definition at line 904 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::axis_bins(), m_inx_energy, m_inx_theta, m_lookup, and table_index().

Referenced by prepare_table(), and set_from_buffer().

GCTAModelSpatialLookup & GCTAModelSpatialLookup::operator= ( const GCTAModelSpatialLookup model)
virtual

Assignment operator.

Parameters
[in]modelLookup table spatial model.

Definition at line 261 of file GCTAModelSpatialLookup.cpp.

References copy_members(), free_members(), init_members(), and GCTAModelSpatial::operator=().

void GCTAModelSpatialLookup::prepare_table ( void  )
protected

Prepare lookup table indices.

Prepares a lookup table for usage.

The method sets the data members m_inx_energy and m_inx_theta which determine the axes order.

The method furthermore sets the energy axis to logarithmic scale and the offset angles to radians.

Finally, the method assures that the maximum values for each radial profile are one.

Definition at line 878 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::axis(), GCTAResponseTable::axis_log10(), GCTAResponseTable::axis_radians(), m_inx_energy, m_inx_theta, m_lookup, and normalise_table().

Referenced by GCTAModelSpatialLookup(), read(), and table().

std::string GCTAModelSpatialLookup::print ( const GChatter chatter = NORMAL) const
virtual

Print lookup table information.

Parameters
[in]chatterChattiness.
Returns
String containing lookup table information.

Implements GCTAModelSpatial.

Definition at line 748 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::axes(), GCTAResponseTable::axis_bins(), GCTAResponseTable::axis_hi(), GCTAResponseTable::axis_lo(), m_filename, m_inx_energy, m_inx_theta, m_lookup, gammalib::parformat(), SILENT, and gammalib::str().

void GCTAModelSpatialLookup::read ( const GXmlElement xml)
virtual

Read model from XML element.

Parameters
[in]xmlXML element.

Read the lookup table model information from an XML element. The XML element needs to be of the following format:

<spatialModel type="LookupTable" file="lookuptable.fits">
  <parameter name="Normalization" ../>
</spatialModel>

The file attribute provides the filename of the lookup table FITS file. The filename may be either an absolute filename (starting with '/') or a relative filename. If no access path is given, the file is expected to reside in the same location as the XML file.

Implements GCTAModelSpatial.

Definition at line 387 of file GCTAModelSpatialLookup.cpp.

References GXmlElement::attribute(), clear(), G_READ_XML, load(), m_norm, GOptimizerPar::name(), norm(), GModelPar::read(), gammalib::xml_file_expand(), and gammalib::xml_get_par().

Referenced by GCTAModelSpatialLookup(), and load().

void GCTAModelSpatialLookup::read ( const GFitsTable table)

Read lookup table from FITS table.

Parameters
[in]tableFITS table.
Exceptions
GException::invalid_valueResponse table is not two-dimensional.

Reads the lookup table form the FITS table. The following column names are mandatory:

ENERG_LO - Energy lower bin boundaries
ENERG_HI - Energy upper bin boundaries
THETA_LO - Offset angle lower bin boundaries
THETA_HI - Offset angle upper bin boundaries

After reading, the method assures that the lookup table is properly normalised, i.e. that for each energy vector the radial profile has a maximum value of one.

Definition at line 598 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::axes(), GCTAResponseTable::clear(), G_READ_TABLE, m_lookup, prepare_table(), GCTAResponseTable::read(), and gammalib::str().

void GCTAModelSpatialLookup::save ( const GFilename filename,
const bool &  clobber = false 
) const

Save lookup table into FITS file.

Parameters
[in]filenameFITS file name.
[in]clobberOverwrite existing file?

Saves lookup table into a FITS file. If a file with the given filename does not yet exist it will be created, otherwise the method opens the existing file. The method will create a (or replace an existing) effective area extension. The extension name can be specified as part of the filename, or if no extension name is given, is assumed to be RADIAL BACKGROUND LOOKUP.

An existing file will only be modified if the clobber flag is set to true.

Definition at line 708 of file GCTAModelSpatialLookup.cpp.

References GFitsHDU::extname(), GFilename::extname(), gammalib::extname_cta_spatial_lookup, GFits::remove(), table(), GFilename::url(), and write().

void GCTAModelSpatialLookup::set_from_buffer ( const std::vector< double > &  buffer)
protected

Set lookup table from buffer.

Parameters
[in]bufferCounts buffer.

Definition at line 1042 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::axis_bins(), GCTAResponseTable::axis_hi(), GCTAResponseTable::axis_lo(), gammalib::deg2rad, m_inx_energy, m_inx_theta, m_lookup, normalise_table(), gammalib::pi, and table_index().

Referenced by fill().

const GCTAResponseTable & GCTAModelSpatialLookup::table ( void  ) const
inline

Return lookup table.

Returns
Lookup table.

Returns the lookup table.

Definition at line 169 of file GCTAModelSpatialLookup.hpp.

References m_lookup.

Referenced by GCTAModelSpatialLookup(), load(), save(), and table().

void GCTAModelSpatialLookup::table ( const GCTAResponseTable table)

Assign lookup table.

Parameters
[in]tableLookup table.

Assigns lookup table.

Definition at line 565 of file GCTAModelSpatialLookup.cpp.

References m_lookup, prepare_table(), and table().

int GCTAModelSpatialLookup::table_index ( const int &  ienergy,
const int &  itheta 
) const
protected

Return index of lookup table element.

Parameters
[in]ienergyEnergy index.
[in]ithetaOffset index.
Returns
Index of lookup table element.

Definition at line 946 of file GCTAModelSpatialLookup.cpp.

References GCTAResponseTable::axis_bins(), m_inx_energy, m_inx_theta, and m_lookup.

Referenced by fill_buffer(), normalise_table(), and set_from_buffer().

std::string GCTAModelSpatialLookup::type ( void  ) const
inlinevirtual

Return model type.

Returns
Model type "LookupTable".

Implements GCTAModelSpatial.

Definition at line 140 of file GCTAModelSpatialLookup.hpp.

Referenced by write().

void GCTAModelSpatialLookup::write ( GXmlElement xml) const
virtual

Write model into XML element.

Parameters
[in]xmlXML element.
Exceptions
GException::invalid_valueSpatial model is not of valid type.

Write the lookup table model information into an XML element. The XML element will be of the following format:

<spatialModel type="LookupTable" file="lookuptable.fits">
  <parameter name="Normalization" ../>
</spatialModel>

The file attribute provides the filename of the lookup table FITS file. The filename may be either an absolute filename (starting with '/') or a relative filename. If no access path is given, the file is expected to reside in the same location as the XML file.

Implements GCTAModelSpatial.

Definition at line 426 of file GCTAModelSpatialLookup.cpp.

References GXmlElement::attribute(), G_WRITE_XML, m_filename, m_norm, GOptimizerPar::name(), norm(), type(), GModelPar::write(), gammalib::xml_file_reduce(), and gammalib::xml_need_par().

Referenced by save().

void GCTAModelSpatialLookup::write ( GFitsBinTable table) const

Write lookup table into FITS binary table.

Parameters
[in]tableFITS binary table.

Writes lookup table into a FITS binary table.

Definition at line 629 of file GCTAModelSpatialLookup.cpp.

References m_lookup, and GCTAResponseTable::write().

Member Data Documentation

GFilename GCTAModelSpatialLookup::m_filename
protected

Name of lookup table.

Definition at line 114 of file GCTAModelSpatialLookup.hpp.

Referenced by copy_members(), init_members(), load(), print(), and write().

int GCTAModelSpatialLookup::m_inx_energy
protected
int GCTAModelSpatialLookup::m_inx_theta
protected
GModelPar GCTAModelSpatialLookup::m_norm
protected

Normalization factor.

Definition at line 116 of file GCTAModelSpatialLookup.hpp.

Referenced by copy_members(), eval(), init_members(), norm(), read(), and write().


The documentation for this class was generated from the following files: