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

Spatial energy dependent Gaussian model class. More...

#include <GCTAModelSpatialGaussSpectrum.hpp>

Inheritance diagram for GCTAModelSpatialGaussSpectrum:
GCTAModelSpatial GBase

Public Member Functions

 GCTAModelSpatialGaussSpectrum (void)
 Void constructor. More...
 
 GCTAModelSpatialGaussSpectrum (const double &sigma)
 Energy-independent Gaussian constructor. More...
 
 GCTAModelSpatialGaussSpectrum (const GModelSpectral &sigma)
 Gaussian spectrum constructor. More...
 
 GCTAModelSpatialGaussSpectrum (const GXmlElement &xml)
 XML constructor. More...
 
 GCTAModelSpatialGaussSpectrum (const GCTAModelSpatialGaussSpectrum &model)
 Copy constructor. More...
 
virtual ~GCTAModelSpatialGaussSpectrum (void)
 Destructor. More...
 
virtual
GCTAModelSpatialGaussSpectrum
operator= (const GCTAModelSpatialGaussSpectrum &model)
 Assignment operator. More...
 
virtual void clear (void)
 Clear instance. More...
 
virtual
GCTAModelSpatialGaussSpectrum
clone (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 point source information. More...
 
const GModelSpectralsigma (void) const
 Return pointer to sigma spectrum. More...
 
void sigma (const double &sigma)
 Set sigma spectrum from value. More...
 
void sigma (const GModelSpectral &sigma)
 Set sigma spectrum. 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 GCTAModelSpatialGaussSpectrum &model)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void set_pointers (void)
 Set pointers. 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

GModelSpectralm_sigma
 
- Protected Attributes inherited from GCTAModelSpatial
std::vector< GModelPar * > m_pars
 Parameter pointers. More...
 

Detailed Description

Spatial energy dependent Gaussian model class.

Definition at line 48 of file GCTAModelSpatialGaussSpectrum.hpp.

Constructor & Destructor Documentation

GCTAModelSpatialGaussSpectrum::GCTAModelSpatialGaussSpectrum ( void  )

Void constructor.

Definition at line 68 of file GCTAModelSpatialGaussSpectrum.cpp.

References init_members().

Referenced by clone().

GCTAModelSpatialGaussSpectrum::GCTAModelSpatialGaussSpectrum ( const double &  sigma)
explicit

Energy-independent Gaussian constructor.

Parameters
[in]sigmaGaussian width (degrees \(^2\)).

Definition at line 84 of file GCTAModelSpatialGaussSpectrum.cpp.

References init_members(), and sigma().

GCTAModelSpatialGaussSpectrum::GCTAModelSpatialGaussSpectrum ( const GModelSpectral sigma)
explicit

Gaussian spectrum constructor.

Parameters
[in]sigmaSpectrum defining the energy-deependent Gaussian width (degrees \(^2\)).

Definition at line 104 of file GCTAModelSpatialGaussSpectrum.cpp.

References init_members(), and sigma().

GCTAModelSpatialGaussSpectrum::GCTAModelSpatialGaussSpectrum ( const GXmlElement xml)
explicit

XML constructor.

Parameters
[in]xmlXML element.

Creates instance of a spatial energy-dependent Gaussian model by extracting information from an XML element. See GCTAModelSpatialGaussSpectrum::read() for more information about the expected structure of the XML element.

Definition at line 128 of file GCTAModelSpatialGaussSpectrum.cpp.

References init_members(), and read().

GCTAModelSpatialGaussSpectrum::GCTAModelSpatialGaussSpectrum ( const GCTAModelSpatialGaussSpectrum model)

Copy constructor.

Parameters
[in]modelEnergy-dependent Gaussian model.

Definition at line 147 of file GCTAModelSpatialGaussSpectrum.cpp.

References copy_members(), and init_members().

GCTAModelSpatialGaussSpectrum::~GCTAModelSpatialGaussSpectrum ( void  )
virtual

Destructor.

Definition at line 164 of file GCTAModelSpatialGaussSpectrum.cpp.

References free_members().

Member Function Documentation

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

Return class name.

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

Implements GCTAModelSpatial.

Definition at line 100 of file GCTAModelSpatialGaussSpectrum.hpp.

void GCTAModelSpatialGaussSpectrum::clear ( void  )
virtual

Clear instance.

Implements GCTAModelSpatial.

Definition at line 218 of file GCTAModelSpatialGaussSpectrum.cpp.

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

Referenced by read(), and sigma().

GCTAModelSpatialGaussSpectrum * GCTAModelSpatialGaussSpectrum::clone ( void  ) const
virtual

Clone instance.

Implements GCTAModelSpatial.

Definition at line 236 of file GCTAModelSpatialGaussSpectrum.cpp.

References GCTAModelSpatialGaussSpectrum().

void GCTAModelSpatialGaussSpectrum::copy_members ( const GCTAModelSpatialGaussSpectrum model)
protected

Copy class members.

Parameters
[in]modelEnergy-dependent Gaussian model.

Definition at line 530 of file GCTAModelSpatialGaussSpectrum.cpp.

References GModelSpectral::clone(), m_sigma, and set_pointers().

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

double GCTAModelSpatialGaussSpectrum::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 energy-dependent Gaussian model for a given event direction. The Gaussian model is defined by

\[f(\theta,E) = \exp \left(-\frac{1}{2} \left( \frac{\theta^2}{\sigma(E)} \right)^2 \right)\]

where \(\theta\) is the offset angle (in degrees), and \(\sigma(E)\) is the energy-dependent Gaussian width (in degrees \(^2\)).

If the gradients flag is true the method will also compute the partial derivatives of the parameters. The partial derivative of the Gaussian width is given by

\[\frac{df}{d\sigma_v} = f(\theta) \frac{\theta^4}{\sigma^3} \sigma_s\]

where \(\sigma_v\) is the value part, \(\sigma_s\) is the scaling part, and \(\sigma = \sigma_v \sigma_s\).

Note that this method implements a function which is unity for \(\theta=0\).

Implements GCTAModelSpatial.

Definition at line 275 of file GCTAModelSpatialGaussSpectrum.cpp.

References GModelSpectral::eval(), exp(), GOptimizerPar::gradient(), gammalib::is_infinite(), gammalib::is_notanumber(), m_sigma, gammalib::rad2deg, sigma(), GModelSpectral::size(), and GCTAInstDir::theta().

void GCTAModelSpatialGaussSpectrum::free_members ( void  )
protected

Delete class members.

Definition at line 548 of file GCTAModelSpatialGaussSpectrum.cpp.

References m_sigma.

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

void GCTAModelSpatialGaussSpectrum::init_members ( void  )
protected

Initialise class members.

Definition at line 512 of file GCTAModelSpatialGaussSpectrum.cpp.

References GCTAModelSpatial::m_pars, and m_sigma.

Referenced by clear(), GCTAModelSpatialGaussSpectrum(), and operator=().

double GCTAModelSpatialGaussSpectrum::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 127 of file GCTAModelSpatialGaussSpectrum.hpp.

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

Assignment operator.

Parameters
[in]modelEnergy-dependent Gaussian model.

Definition at line 185 of file GCTAModelSpatialGaussSpectrum.cpp.

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

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

Print point source information.

Parameters
[in]chatterChattiness.
Returns
String containing point source information.

Implements GCTAModelSpatial.

Definition at line 471 of file GCTAModelSpatialGaussSpectrum.cpp.

References m_sigma, gammalib::parformat(), SILENT, GModelSpectral::size(), and GModelSpectral::type().

void GCTAModelSpatialGaussSpectrum::read ( const GXmlElement xml)
virtual

Read model from XML element.

Parameters
[in]xmlXML element.

Read the energy-dependent Gaussian spatial model information from an XML element. The XML element needs to be of the following format:

<spatialModel type="EnergyDependentGaussian">
  <sigma type="...">
    ...
  </sigma>
</spatialModel>

where any spectral model type can be provided for the sigma tag.

Implements GCTAModelSpatial.

Definition at line 351 of file GCTAModelSpatialGaussSpectrum.cpp.

References GModelSpectralRegistry::alloc(), clear(), GXmlNode::element(), m_sigma, set_pointers(), and sigma().

Referenced by GCTAModelSpatialGaussSpectrum().

void GCTAModelSpatialGaussSpectrum::set_pointers ( void  )
protected

Set pointers.

Set pointers to all model parameters. The pointers are stored in a vector that is member of the GModelData base class.

Definition at line 567 of file GCTAModelSpatialGaussSpectrum.cpp.

References GCTAModelSpatial::m_pars, m_sigma, and GModelSpectral::size().

Referenced by copy_members(), read(), and sigma().

const GModelSpectral * GCTAModelSpatialGaussSpectrum::sigma ( void  ) const
inline

Return pointer to sigma spectrum.

Returns
Pointer to sigma spectrum.

Definition at line 139 of file GCTAModelSpatialGaussSpectrum.hpp.

References m_sigma.

Referenced by eval(), GCTAModelSpatialGaussSpectrum(), read(), and write().

void GCTAModelSpatialGaussSpectrum::sigma ( const double &  sigma)

Set sigma spectrum from value.

Parameters
[in]sigmaConstant sigma value.

Definition at line 428 of file GCTAModelSpatialGaussSpectrum.cpp.

References clear(), m_sigma, and set_pointers().

void GCTAModelSpatialGaussSpectrum::sigma ( const GModelSpectral sigma)

Set sigma spectrum.

Parameters
[in]sigmaSigma spectrum.

Definition at line 449 of file GCTAModelSpatialGaussSpectrum.cpp.

References clear(), GModelSpectral::clone(), m_sigma, and set_pointers().

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

Return model type.

Returns
Model type "EnergyDependentGaussian".

Implements GCTAModelSpatial.

Definition at line 112 of file GCTAModelSpatialGaussSpectrum.hpp.

Referenced by write().

void GCTAModelSpatialGaussSpectrum::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 energy-dependent Gaussian spatial model information into an XML element. The XML element will be of the following format:

<spatialModel type="EnergyDependentGaussian">
  <sigma type="...">
    ...
  </sigma>
</spatialModel>

where the spectral sigma model will be written into the sigma tag.

Implements GCTAModelSpatial.

Definition at line 390 of file GCTAModelSpatialGaussSpectrum.cpp.

References GXmlNode::append(), GXmlElement::attribute(), GXmlNode::element(), GXmlNode::elements(), G_WRITE, m_sigma, sigma(), type(), and GModelSpectral::write().

Member Data Documentation

GModelSpectral* GCTAModelSpatialGaussSpectrum::m_sigma
protected

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