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

Abstract elliptical spatial model base class. More...

#include <GModelSpatialElliptical.hpp>

Inheritance diagram for GModelSpatialElliptical:
GModelSpatial GBase GModelSpatialEllipticalDisk GModelSpatialEllipticalGauss GModelSpatialEllipticalGeneralGauss

Public Member Functions

 GModelSpatialElliptical (void)
 Void constructor. More...
 
 GModelSpatialElliptical (const GModelSpatialElliptical &model)
 Copy constructor. More...
 
 GModelSpatialElliptical (const GXmlElement &xml)
 XML constructor. More...
 
virtual ~GModelSpatialElliptical (void)
 Destructor. More...
 
virtual GModelSpatialEllipticaloperator= (const GModelSpatialElliptical &model)
 Assignment operator. More...
 
virtual void clear (void)=0
 Clear object. More...
 
virtual GModelSpatialEllipticalclone (void) const =0
 Clones object. More...
 
virtual std::string classname (void) const =0
 Return class name. More...
 
virtual double eval (const double &theta, const double &posangle, const GEnergy &energy, const GTime &time, const bool &gradients=false) const =0
 
virtual GSkyDir mc (const GEnergy &energy, const GTime &time, GRan &ran) const =0
 
virtual bool contains (const GSkyDir &dir, const double &margin=0.0) const =0
 
virtual double theta_max (void) const =0
 
virtual std::string print (const GChatter &chatter=NORMAL) const =0
 Print content of object. More...
 
virtual GClassCode code (void) const
 Return class code. More...
 
virtual double eval (const GPhoton &photon, const bool &gradients=false) const
 Return model value. More...
 
virtual double mc_norm (const GSkyDir &dir, const double &radius) const
 Return normalization of elliptical source 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...
 
std::string coordsys (void) const
 Return coordinate system. More...
 
const GSkyDirdir (void) const
 Return position of elliptical spatial model. More...
 
void dir (const GSkyDir &dir)
 Set position of radial spatial model. More...
 
double posangle (void) const
 Return Position Angle of model. More...
 
void posangle (const double &posangle)
 Set Position Angle of model. More...
 
double semiminor (void) const
 Return semi-minor axis of ellipse. More...
 
double semimajor (void) const
 Return semi-major axis of ellipse. More...
 
void semiminor (const double &semiminor)
 Set semi-minor axis of ellipse. More...
 
void semimajor (const double &semimajor)
 Set semi-major axis of ellipse. More...
 
- Public Member Functions inherited from GModelSpatial
 GModelSpatial (void)
 Void constructor. More...
 
 GModelSpatial (const GModelSpatial &model)
 Copy constructor. More...
 
virtual ~GModelSpatial (void)
 Destructor. More...
 
virtual GModelSpatialoperator= (const GModelSpatial &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 model parameter. More...
 
virtual const GModelParoperator[] (const std::string &name) const
 Returns model parameter (const version) More...
 
virtual double flux (const GSkyRegion &region, const GEnergy &srcEng=GEnergy(), const GTime &srcTime=GTime()) const
 Returns model flux within sky region. More...
 
std::string type (void) const
 Return model type. More...
 
void type (const std::string &type)
 Set model type. More...
 
GModelParat (const int &index)
 Returns model parameter. More...
 
const GModelParat (const int &index) const
 Returns model parameter (const version) More...
 
bool has_par (const std::string &name) const
 Checks if parameter name exists. More...
 
bool has_free_pars (void) const
 Checks if the spatial model has free parameters. More...
 
int size (void) const
 Return number of parameters. More...
 
void autoscale (void)
 Autoscale parameters. More...
 
const GSkyRegionregion (void) const
 Return boundary sky region. 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 GModelSpatialElliptical &model)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
bool is_celestial (void) const
 Check if model holds celestial coordinates. More...
 
virtual void set_region (void) const =0
 
- Protected Member Functions inherited from GModelSpatial
void init_members (void)
 Initialise class members. More...
 
void copy_members (const GModelSpatial &model)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 

Protected Attributes

GModelPar m_lon
 Right Ascension or Galactic longitude (deg) More...
 
GModelPar m_lat
 Declination or Galactic latitude (deg) More...
 
GModelPar m_posangle
 Position angle from North, counterclockwise (deg) More...
 
GModelPar m_semiminor
 Semi-minor axis of ellipse (deg) More...
 
GModelPar m_semimajor
 Semi-major axis of ellipse (deg) More...
 
GSkyDir m_dir
 Sky direction representing parameters. More...
 
double m_last_lon
 Last longitude. More...
 
double m_last_lat
 Last latitude. More...
 
- Protected Attributes inherited from GModelSpatial
std::string m_type
 Spatial model type. More...
 
GSkyRegionCircle m_region
 Bounding circle. More...
 
std::vector< GModelPar * > m_pars
 Parameter pointers. More...
 

Detailed Description

Abstract elliptical spatial model base class.

This class defines the interface for an elliptical model as spatial component of the factorized source model. Typical examples of elliptical components are elliptical Disk, Gaussian or Shell shaped sources.

Definition at line 55 of file GModelSpatialElliptical.hpp.

Constructor & Destructor Documentation

GModelSpatialElliptical::GModelSpatialElliptical ( void  )

Void constructor.

Definition at line 54 of file GModelSpatialElliptical.cpp.

References init_members().

GModelSpatialElliptical::GModelSpatialElliptical ( const GModelSpatialElliptical model)

Copy constructor.

Parameters
[in]modelElliptical spatial model.

Definition at line 91 of file GModelSpatialElliptical.cpp.

References copy_members(), and init_members().

GModelSpatialElliptical::GModelSpatialElliptical ( const GXmlElement xml)
explicit

XML constructor.

Parameters
[in]xmlXML element.

Constructs an elliptical spatial model component by extracting information from an XML element. See the read() method for more information about the expected structure of the XML element.

Definition at line 73 of file GModelSpatialElliptical.cpp.

References init_members(), and read().

GModelSpatialElliptical::~GModelSpatialElliptical ( void  )
virtual

Destructor.

Definition at line 108 of file GModelSpatialElliptical.cpp.

References free_members().

Member Function Documentation

virtual std::string GModelSpatialElliptical::classname ( void  ) const
pure virtual

Return class name.

Returns
String containing the class name.

Returns the class name for non-abstract classes in a human readable way.

Implements GModelSpatial.

Implemented in GModelSpatialEllipticalGeneralGauss, GModelSpatialEllipticalGauss, and GModelSpatialEllipticalDisk.

virtual void GModelSpatialElliptical::clear ( void  )
pure virtual

Clear object.

Sets the object to a clean initial state. After calling the method the object will be in the same state as it were if an empty instance of the object would have been created.

Implements GModelSpatial.

Implemented in GModelSpatialEllipticalGeneralGauss, GModelSpatialEllipticalGauss, and GModelSpatialEllipticalDisk.

virtual GModelSpatialElliptical* GModelSpatialElliptical::clone ( void  ) const
pure virtual

Clones object.

Returns
Pointer to deep copy of object.

Creates a deep copy of the object and returns a pointer to the object.

Implements GModelSpatial.

Implemented in GModelSpatialEllipticalGeneralGauss, GModelSpatialEllipticalGauss, and GModelSpatialEllipticalDisk.

GClassCode GModelSpatialElliptical::code ( void  ) const
inlinevirtual

Return class code.

Returns
GModelSpatialElliptical.

Returns the code GModelSpatialElliptical of the class.

Implements GModelSpatial.

Definition at line 133 of file GModelSpatialElliptical.hpp.

References GMODEL_SPATIAL_ELLIPTICAL.

virtual bool GModelSpatialElliptical::contains ( const GSkyDir dir,
const double &  margin = 0.0 
) const
pure virtual
std::string GModelSpatialElliptical::coordsys ( void  ) const
inline

Return coordinate system.

Returns
Coordinate system of point source model.

Returns "CEL" for a celestial coordinate system and "GAL" for a Galactic coordinate system.

Definition at line 258 of file GModelSpatialElliptical.hpp.

References is_celestial().

void GModelSpatialElliptical::copy_members ( const GModelSpatialElliptical model)
protected

Copy class members.

Parameters
[in]modelRadial spatial model.

Definition at line 441 of file GModelSpatialElliptical.cpp.

References m_dir, m_last_lat, m_last_lon, m_lat, m_lon, GModelSpatial::m_pars, m_posangle, m_semimajor, and m_semiminor.

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

void GModelSpatialElliptical::dir ( const GSkyDir dir)

Set position of radial spatial model.

Parameters
[in]dirSky direction of radial spatial model.

Sets the sky direction of the radial spatial model.

Definition at line 341 of file GModelSpatialElliptical.cpp.

References GSkyDir::b_deg(), GSkyDir::dec_deg(), is_celestial(), GSkyDir::l_deg(), m_lat, m_lon, GSkyDir::ra_deg(), and GOptimizerPar::value().

virtual double GModelSpatialElliptical::eval ( const double &  theta,
const double &  posangle,
const GEnergy energy,
const GTime time,
const bool &  gradients = false 
) const
pure virtual
double GModelSpatialElliptical::eval ( const GPhoton photon,
const bool &  gradients = false 
) const
virtual

Return model value.

Parameters
[in]photonIncident Photon.
[in]gradientsCompute gradients?
Returns
Value of spatial elliptical model.

Evaluates the elliptical spatial model value for a specific incident photon.

If the gradients flag is true the method will also compute the parameter gradients for all model parameters.

Implements GModelSpatial.

Definition at line 173 of file GModelSpatialElliptical.cpp.

References GPhoton::dir(), dir(), GSkyDir::dist(), GPhoton::energy(), eval(), GSkyDir::posang(), and GPhoton::time().

void GModelSpatialElliptical::free_members ( void  )
protected
bool GModelSpatialElliptical::is_celestial ( void  ) const
inlineprotected

Check if model holds celestial coordinates.

Returns
True if model holds celestial coordinates, false otherwise.

Definition at line 270 of file GModelSpatialElliptical.hpp.

References m_lon, and GOptimizerPar::name().

Referenced by coordsys(), and dir().

virtual GSkyDir GModelSpatialElliptical::mc ( const GEnergy energy,
const GTime time,
GRan ran 
) const
pure virtual
double GModelSpatialElliptical::mc_norm ( const GSkyDir dir,
const double &  radius 
) const
inlinevirtual

Return normalization of elliptical source for Monte Carlo simulations.

Parameters
[in]dirCentre of simulation cone.
[in]radiusRadius of simulation cone (degrees).
Returns
Normalization.

Returns the normalization for an elliptical source within a circular region. The normalization is 1 if the elliptical source falls within the circle define by dir and radius, 0 otherwise.

Implements GModelSpatial.

Definition at line 241 of file GModelSpatialElliptical.hpp.

References dir(), GSkyDir::dist_deg(), norm(), and theta_max().

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

Assignment operator.

Parameters
[in]modelElliptical spatial model.
Returns
Elliptical spatial model.

Definition at line 130 of file GModelSpatialElliptical.cpp.

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

Referenced by GModelSpatialEllipticalDisk::operator=(), GModelSpatialEllipticalGauss::operator=(), and GModelSpatialEllipticalGeneralGauss::operator=().

void GModelSpatialElliptical::posangle ( const double &  posangle)
inline

Set Position Angle of model.

Parameters
[in]posanglePosition Angle of model (degrees).

Sets the Position Angle of model in degrees, measured counterclockwise from celestial North.

Definition at line 163 of file GModelSpatialElliptical.hpp.

References m_posangle, and GOptimizerPar::value().

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

Print content of object.

Parameters
[in]chatterChattiness (defaults to NORMAL).
Returns
String containing the content of the object.

Formats the content in a standard way and puts this content in a C++ string that is returned.

Implements GModelSpatial.

Implemented in GModelSpatialEllipticalGeneralGauss, GModelSpatialEllipticalGauss, and GModelSpatialEllipticalDisk.

void GModelSpatialElliptical::read ( const GXmlElement xml)
virtual

Read model from XML element.

Parameters
[in]xmlXML element.

Reads the elliptical source location and position angle information from an XML element in the following format

<spatialModel type="...">
  <parameter name="RA"  scale="1" value="83.63" min="-360" max="360" free="1"/>
  <parameter name="DEC" scale="1" value="22.01" min="-90"  max="90"  free="1"/>
  <parameter name="PA"  scale="1" value="45.0"  min="-360" max="360" free="1"/>
  ...
</spatialModel>

or

<spatialModel type="...">
  <parameter name="GLON" scale="1" value="83.63" min="-360" max="360" free="1"/>
  <parameter name="GLAT" scale="1" value="22.01" min="-90"  max="90"  free="1"/>
  <parameter name="PA"   scale="1" value="45.0"  min="-360" max="360" free="1"/>
  ...
</spatialModel>

Implements GModelSpatial.

Reimplemented in GModelSpatialEllipticalGeneralGauss, GModelSpatialEllipticalGauss, and GModelSpatialEllipticalDisk.

Definition at line 214 of file GModelSpatialElliptical.cpp.

References G_READ, m_lat, m_lon, m_posangle, GOptimizerPar::name(), GModelPar::read(), gammalib::xml_get_par(), and gammalib::xml_has_par().

Referenced by GModelSpatialElliptical(), GModelSpatialEllipticalDisk::read(), GModelSpatialEllipticalGauss::read(), and GModelSpatialEllipticalGeneralGauss::read().

void GModelSpatialElliptical::semimajor ( const double &  semimajor)
inline

Set semi-major axis of ellipse.

Parameters
[in]semimajorSemi-major axis of ellipse (degrees)

Sets the semi-major axis of the ellipse in degrees.

Definition at line 221 of file GModelSpatialElliptical.hpp.

References m_semimajor, and GOptimizerPar::value().

void GModelSpatialElliptical::semiminor ( const double &  semiminor)
inline

Set semi-minor axis of ellipse.

Parameters
[in]semiminorSemi-minor axis of ellipse (degrees)

Sets the semi-minor axis of the ellipse in degrees.

Definition at line 192 of file GModelSpatialElliptical.hpp.

References m_semiminor, and GOptimizerPar::value().

virtual void GModelSpatialElliptical::set_region ( void  ) const
protectedpure virtual
void GModelSpatialElliptical::write ( GXmlElement xml) const
virtual

Write model into XML element.

Parameters
[in]xmlXML element into which model information is written.

Depending on the coordinate system, write the elliptical source information into an XML element with the following format

<spatialModel type="...">
  <parameter name="RA"  scale="1" value="83.6331" min="-360" max="360" free="0" />
  <parameter name="DEC" scale="1" value="22.0145" min="-90"  max="90"  free="0" />
  <parameter name="PA"  scale="1" value="45.0"    min="-360" max="360" free="1" />
  ...
</spatialModel>

or

<spatialModel type="PointSource">
  <parameter name="GLON" scale="1" value="184.5575" min="-360" max="360" free="0" />
  <parameter name="GLAT" scale="1" value="-5.7843"  min="-90"  max="90"  free="0" />
  <parameter name="PA"   scale="1" value="45.0"     min="-360" max="360" free="1" />
  ...
</spatialModel>

Implements GModelSpatial.

Reimplemented in GModelSpatialEllipticalGeneralGauss, GModelSpatialEllipticalGauss, and GModelSpatialEllipticalDisk.

Definition at line 278 of file GModelSpatialElliptical.cpp.

References G_WRITE, m_lat, m_lon, m_posangle, GOptimizerPar::name(), GModelSpatial::type(), GModelPar::write(), gammalib::xml_check_type(), and gammalib::xml_need_par().

Referenced by GModelSpatialEllipticalDisk::write(), GModelSpatialEllipticalGauss::write(), and GModelSpatialEllipticalGeneralGauss::write().

Member Data Documentation

GSkyDir GModelSpatialElliptical::m_dir
mutableprotected

Sky direction representing parameters.

Definition at line 119 of file GModelSpatialElliptical.hpp.

Referenced by copy_members(), dir(), and init_members().

double GModelSpatialElliptical::m_last_lat
mutableprotected

Last latitude.

Definition at line 121 of file GModelSpatialElliptical.hpp.

Referenced by copy_members(), dir(), and init_members().

double GModelSpatialElliptical::m_last_lon
mutableprotected

Last longitude.

Definition at line 120 of file GModelSpatialElliptical.hpp.

Referenced by copy_members(), dir(), and init_members().

GModelPar GModelSpatialElliptical::m_posangle
protected

Position angle from North, counterclockwise (deg)

Definition at line 114 of file GModelSpatialElliptical.hpp.

Referenced by copy_members(), GModelSpatialEllipticalDisk::eval(), init_members(), posangle(), read(), and write().


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