GammaLib
2.0.0
|
Point source spatial model. More...
#include <GModelSpatialPointSource.hpp>
Public Member Functions | |
GModelSpatialPointSource (void) | |
Void constructor. More... | |
GModelSpatialPointSource (const bool &dummy, const std::string &type) | |
Model type constructor. More... | |
GModelSpatialPointSource (const GSkyDir &dir, const std::string &coordsys="CEL") | |
Sky direction constructor. More... | |
GModelSpatialPointSource (const double &lon, const double &lat, const std::string &coordsys="CEL") | |
Value constructor. More... | |
GModelSpatialPointSource (const GXmlElement &xml) | |
XML constructor. More... | |
GModelSpatialPointSource (const GModelSpatialPointSource &model) | |
Copy constructor. More... | |
virtual | ~GModelSpatialPointSource (void) |
Destructor. More... | |
virtual GModelSpatialPointSource & | operator= (const GModelSpatialPointSource &model) |
Assignment operator. More... | |
virtual void | clear (void) |
Clear point source model. More... | |
virtual GModelSpatialPointSource * | clone (void) const |
Clone point source model. More... | |
virtual std::string | classname (void) const |
Return class name. More... | |
virtual GClassCode | code (void) const |
Return class code. More... | |
virtual double | eval (const GPhoton &photon, const bool &gradients=false) const |
Evaluate function. More... | |
virtual GSkyDir | mc (const GEnergy &energy, const GTime &time, GRan &ran) const |
Returns MC sky direction. More... | |
virtual double | mc_norm (const GSkyDir &dir, const double &radius) const |
Return normalization of point source for Monte Carlo simulations. More... | |
virtual bool | contains (const GSkyDir &dir, const double &margin=0.0) const |
Checks where model contains specified sky direction. 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... | |
virtual double | flux (const GSkyRegion ®ion, const GEnergy &srcEng=GEnergy(), const GTime &srcTime=GTime()) const |
Returns model flux integrated in sky region. More... | |
std::string | coordsys (void) const |
Return coordinate system. More... | |
const GSkyDir & | dir (void) const |
Return position of point source. More... | |
void | dir (const GSkyDir &dir) |
Set position of point source. 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 GModelSpatial & | operator= (const GModelSpatial &model) |
Assignment operator. More... | |
virtual GModelPar & | operator[] (const int &index) |
Returns model parameter. More... | |
virtual const GModelPar & | operator[] (const int &index) const |
Returns model parameter (const version) More... | |
virtual GModelPar & | operator[] (const std::string &name) |
Returns model parameter. More... | |
virtual const GModelPar & | operator[] (const std::string &name) const |
Returns model parameter (const version) More... | |
std::string | type (void) const |
Return model type. More... | |
void | type (const std::string &type) |
Set model type. More... | |
GModelPar & | at (const int &index) |
Returns model parameter. More... | |
const GModelPar & | at (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 GSkyRegion * | region (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 GModelSpatialPointSource &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 |
Set boundary sky region. More... | |
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... | |
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... | |
Point source spatial model.
This class implements a point source as the spatial component of the factorised source model. The point source has two parameters: the longitude and latitude of the point source location. The parameters may either be specified in celestial or Galactic coordinates. For celestial coordinates the parameter names are "RA" and "DEC", for Galactic coordinates the parameter names are "GLON" and "GLAT".
Definition at line 57 of file GModelSpatialPointSource.hpp.
GModelSpatialPointSource::GModelSpatialPointSource | ( | void | ) |
Void constructor.
Constructs empty point source model.
Definition at line 75 of file GModelSpatialPointSource.cpp.
References init_members().
Referenced by clone().
GModelSpatialPointSource::GModelSpatialPointSource | ( | const bool & | dummy, |
const std::string & | type | ||
) |
Model type constructor.
[in] | dummy | Dummy flag. |
[in] | type | Model type. |
Constructs empty point source model by specifying a model type
.
Definition at line 93 of file GModelSpatialPointSource.cpp.
References init_members(), GModelSpatial::m_type, and GModelSpatial::type().
GModelSpatialPointSource::GModelSpatialPointSource | ( | const GSkyDir & | dir, |
const std::string & | coordsys = "CEL" |
||
) |
Sky direction constructor.
[in] | dir | Sky direction. |
[in] | coordsys | Coordinate system (either "CEL" or "GAL") |
GException::invalid_argument | Invalid coordsys argument specified. |
Construct a point source spatial model from a sky direction. The coordsys
parameter specifies whether the sky direction should be interpreted in the celestial or Galactic coordinate system.
Definition at line 121 of file GModelSpatialPointSource.cpp.
References dir(), G_CONSTRUCTOR1, init_members(), m_lat, m_lon, and GOptimizerPar::name().
GModelSpatialPointSource::GModelSpatialPointSource | ( | const double & | lon, |
const double & | lat, | ||
const std::string & | coordsys = "CEL" |
||
) |
Value constructor.
[in] | lon | Longitude of point source (deg). |
[in] | lat | Latitude of point source (deg). |
[in] | coordsys | Coordinate system (either "CEL" or "GAL") |
GException::invalid_argument | Invalid coordsys argument specified. |
Construct a point source spatial model from the longitude and latitude of a point source. The coordsys
parameter specifies the coordinate system in which lon
and lat
are provided. By default a celestial coordinate system is assumed, which means that lon
and lat
are taken as Right Ascension and Declination. If "GAL" is specified, lon
and lat
are taken as Galactic longitude and latitude.
Definition at line 171 of file GModelSpatialPointSource.cpp.
References G_CONSTRUCTOR2, init_members(), m_lat, m_lon, GOptimizerPar::name(), and GOptimizerPar::value().
|
explicit |
XML constructor.
[in] | xml | XML element. |
Construct a point source spatial model 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 215 of file GModelSpatialPointSource.cpp.
References init_members(), and read().
GModelSpatialPointSource::GModelSpatialPointSource | ( | const GModelSpatialPointSource & | model | ) |
Copy constructor.
[in] | model | Point source spatial model. |
Definition at line 234 of file GModelSpatialPointSource.cpp.
References copy_members(), and init_members().
|
virtual |
|
inlinevirtual |
Return class name.
Implements GModelSpatial.
Definition at line 128 of file GModelSpatialPointSource.hpp.
|
virtual |
Clear point source model.
Implements GModelSpatial.
Definition at line 306 of file GModelSpatialPointSource.cpp.
References free_members(), GModelSpatial::free_members(), init_members(), and GModelSpatial::init_members().
|
virtual |
Clone point source model.
Implements GModelSpatial.
Definition at line 326 of file GModelSpatialPointSource.cpp.
References GModelSpatialPointSource().
|
inlinevirtual |
Return class code.
Returns the code GModelSpatialPointSource of the class.
Implements GModelSpatial.
Definition at line 142 of file GModelSpatialPointSource.hpp.
References GMODEL_SPATIAL_POINT_SOURCE.
|
virtual |
Checks where model contains specified sky direction.
[in] | dir | Sky direction. |
[in] | margin | Margin to be added to sky direction (degrees) |
Signals whether a sky direction is contained in the point source model.
Implements GModelSpatial.
Definition at line 391 of file GModelSpatialPointSource.cpp.
References gammalib::deg2rad, dir(), and GSkyDir::dist().
|
inline |
Return coordinate system.
Returns "CEL" for a celestial coordinate system and "GAL" for a Galactic coordinate system.
Definition at line 177 of file GModelSpatialPointSource.hpp.
References is_celestial().
|
protected |
Copy class members.
[in] | model | Point source spatial model. |
Definition at line 670 of file GModelSpatialPointSource.cpp.
References m_dir, m_last_lat, m_last_lon, m_lat, m_lon, GModelSpatial::m_pars, and GModelSpatial::m_type.
Referenced by GModelSpatialPointSource(), and operator=().
const GSkyDir & GModelSpatialPointSource::dir | ( | void | ) | const |
Return position of point source.
Returns the sky direction of the point source.
Definition at line 565 of file GModelSpatialPointSource.cpp.
References is_celestial(), GSkyDir::lb_deg(), m_dir, m_last_lat, m_last_lon, m_lat, m_lon, GSkyDir::radec_deg(), and GOptimizerPar::value().
Referenced by contains(), eval(), GModelSpatialPointSource(), GResponse::irf_ptsrc(), GCTAResponseCube::irf_ptsrc(), GCTAResponseIrf::irf_ptsrc(), GLATResponse::irf_spatial_bin(), mc(), mc_norm(), GCTAResponseIrf::nroi_ptsrc(), and set_region().
void GModelSpatialPointSource::dir | ( | const GSkyDir & | dir | ) |
Set position of point source.
[in] | dir | Sky direction of point source. |
Sets the sky direction of the point source.
Definition at line 601 of file GModelSpatialPointSource.cpp.
References GSkyDir::b_deg(), GSkyDir::dec_deg(), is_celestial(), GSkyDir::l_deg(), m_lat, m_lon, GSkyDir::ra_deg(), and GOptimizerPar::value().
|
virtual |
Evaluate function.
[in] | photon | Incident photon. |
[in] | gradients | Compute gradients? |
Evaluates the spatial part for a point source model. It implements a delta function with respect to the coordinates of the source. For numerical reasons, a certain tolerance is accepted (typically 0.1 arcsec, i.e. well below the angular resolution of gamma-ray telescopes).
The method will not compute analytical parameter gradients, even if the gradients
argument is set to true. Point source parameter gradients need to be computed numerically.
Implements GModelSpatial.
Definition at line 349 of file GModelSpatialPointSource.cpp.
References GPhoton::dir(), dir(), GSkyDir::dist_deg(), and tolerance.
|
virtual |
Returns model flux integrated in sky region.
[in] | region | Sky region. |
[in] | srcEng | Energy. |
[in] | srcTime | Time. |
Returns point source flux within a sky region. If the point source is contained within the sky region the flux will be 1, otherwise the flux will be 0.
Reimplemented from GModelSpatial.
Definition at line 510 of file GModelSpatialPointSource.cpp.
References GModelSpatial::region().
|
protected |
Delete class members.
Definition at line 695 of file GModelSpatialPointSource.cpp.
Referenced by clear(), operator=(), and ~GModelSpatialPointSource().
|
protected |
Initialise class members.
Definition at line 627 of file GModelSpatialPointSource.cpp.
References GSkyDir::clear(), GOptimizerPar::clear(), GOptimizerPar::fix(), GOptimizerPar::gradient(), GOptimizerPar::has_grad(), m_dir, m_last_lat, m_last_lon, m_lat, m_lon, GModelSpatial::m_pars, GModelSpatial::m_type, GOptimizerPar::name(), GOptimizerPar::scale(), and GOptimizerPar::unit().
Referenced by clear(), GModelSpatialPointSource(), and operator=().
|
inlineprotected |
Check if model holds celestial coordinates.
Definition at line 189 of file GModelSpatialPointSource.hpp.
References m_lon, and GOptimizerPar::name().
Referenced by coordsys(), and dir().
|
virtual |
Returns MC sky direction.
[in] | energy | Photon energy. |
[in] | time | Photon arrival time. |
[in,out] | ran | Random number generator. |
Draws an arbitrary sky direction for the point source model. As the point source is a point in the sky, the methods always returns the directon of the point source.
Implements GModelSpatial.
Definition at line 372 of file GModelSpatialPointSource.cpp.
References dir().
|
inlinevirtual |
Return normalization of point source for Monte Carlo simulations.
[in] | dir | Centre of simulation cone. |
[in] | radius | Radius of simulation cone (degrees). |
Returns the normalization for a point source within a circular region. The normalization is 1 if the point source falls within the circle defined by dir
and radius
, 0 otherwise.
Implements GModelSpatial.
Definition at line 160 of file GModelSpatialPointSource.hpp.
References dir(), GSkyDir::dist_deg(), and norm().
|
virtual |
Assignment operator.
[in] | model | Point source spatial model. |
Definition at line 273 of file GModelSpatialPointSource.cpp.
References copy_members(), free_members(), init_members(), and GModelSpatial::operator=().
Print point source information.
[in] | chatter | Chattiness. |
Implements GModelSpatial.
Definition at line 533 of file GModelSpatialPointSource.cpp.
References GModelSpatial::m_pars, gammalib::parformat(), SILENT, GModelSpatial::size(), and gammalib::str().
|
virtual |
Read model from XML element.
[in] | xml | XML element containing point source model information. |
Read the point source information from an XML element with the following format
<spatialModel type="PointSource"> <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" /> </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" /> </spatialModel>
Implements GModelSpatial.
Definition at line 423 of file GModelSpatialPointSource.cpp.
References G_READ, m_lat, m_lon, GModelPar::read(), gammalib::xml_check_parnum(), gammalib::xml_get_par(), and gammalib::xml_has_par().
Referenced by GModelSpatialPointSource().
|
protectedvirtual |
Set boundary sky region.
Implements GModelSpatial.
Definition at line 705 of file GModelSpatialPointSource.cpp.
References dir(), GModelSpatial::m_region, and GModelSpatial::region().
|
virtual |
Write model into XML element.
[in] | xml | XML element into which model information is written. |
Depending on the coordinate system, write the point source information into an XML element with the following format
<spatialModel type="PointSource"> <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" /> </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" /> </spatialModel>
Implements GModelSpatial.
Definition at line 480 of file GModelSpatialPointSource.cpp.
References G_WRITE, m_lat, m_lon, GOptimizerPar::name(), GModelSpatial::type(), GModelPar::write(), gammalib::xml_check_type(), and gammalib::xml_need_par().
|
mutableprotected |
Sky direction representing parameters.
Definition at line 116 of file GModelSpatialPointSource.hpp.
Referenced by copy_members(), dir(), and init_members().
|
mutableprotected |
Last latitude.
Definition at line 118 of file GModelSpatialPointSource.hpp.
Referenced by copy_members(), dir(), and init_members().
|
mutableprotected |
Last longitude.
Definition at line 117 of file GModelSpatialPointSource.hpp.
Referenced by copy_members(), dir(), and init_members().
|
protected |
Declination or Galactic latitude (deg)
Definition at line 113 of file GModelSpatialPointSource.hpp.
Referenced by copy_members(), dir(), GModelSpatialPointSource(), init_members(), read(), and write().
|
protected |
Right Ascension or Galactic longitude (deg)
Definition at line 112 of file GModelSpatialPointSource.hpp.
Referenced by copy_members(), dir(), GModelSpatialPointSource(), init_members(), is_celestial(), read(), and write().