GammaLib  2.0.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GModelSpatialDiffuseCube Class Reference

Spatial map cube model. More...

#include <GModelSpatialDiffuseCube.hpp>

Inheritance diagram for GModelSpatialDiffuseCube:
GModelSpatialDiffuse GModelSpatial GBase

Public Member Functions

 GModelSpatialDiffuseCube (void)
 Void constructor. More...
 
 GModelSpatialDiffuseCube (const bool &dummy, const std::string &type)
 Model type constructor. More...
 
 GModelSpatialDiffuseCube (const GXmlElement &xml)
 XML constructor. More...
 
 GModelSpatialDiffuseCube (const GFilename &filename, const double &value=1.0)
 Filename constructor. More...
 
 GModelSpatialDiffuseCube (const GSkyMap &cube, const GEnergies &energies, const double &value=1.0)
 Sky map constructor. More...
 
 GModelSpatialDiffuseCube (const GModelSpatialDiffuseCube &model)
 Copy constructor. More...
 
virtual ~GModelSpatialDiffuseCube (void)
 Destructor. More...
 
virtual GModelSpatialDiffuseCubeoperator= (const GModelSpatialDiffuseCube &model)
 Assignment operator. More...
 
virtual void clear (void)
 Clear map cube model. More...
 
virtual GModelSpatialDiffuseCubeclone (void) const
 Clone map cube model. More...
 
virtual std::string classname (void) const
 Return class name. 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 diffuse cube for Monte Carlo simulations. More...
 
virtual bool contains (const GSkyDir &dir, const double &margin=0.0) const
 Signals whether model contains 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 map cube information. More...
 
virtual double flux (const GSkyRegion &region, const GEnergy &srcEng=GEnergy(), const GTime &srcTime=GTime()) const
 Returns diffuse cube flux integrated in sky region. More...
 
double value (void) const
 Get model value. More...
 
void value (const double &value)
 Set model value. More...
 
const GFilenamefilename (void) const
 Get file name. More...
 
void filename (const GFilename &filename)
 Set file name. More...
 
const GSkyMapcube (void) const
 Get map cube. More...
 
void cube (const GSkyMap &cube)
 Set map cube. More...
 
GEnergies energies (void)
 Returns map cube energies. More...
 
void energies (const GEnergies &energies)
 Set energies for map cube. More...
 
const GModelSpectralNodesspectrum (void) const
 Get map cube spectrum. More...
 
void mc_cone (const GSkyRegionCircle &cone) const
 Set Monte Carlo simulation cone. More...
 
const GSkyRegionCirclemc_cone (void) const
 Get Monte Carlo simulation cone. More...
 
void load (const GFilename &filename)
 Load cube into the model class. More...
 
void save (const GFilename &filename, const bool &clobber=false) const
 Save cube into FITS file. More...
 
void read (const GFits &fits)
 Read diffuse cube from FITS file. More...
 
void write (GFits &fits) const
 Write diffuse cube into FITS file. More...
 
- Public Member Functions inherited from GModelSpatialDiffuse
 GModelSpatialDiffuse (void)
 Void constructor. More...
 
 GModelSpatialDiffuse (const GModelSpatialDiffuse &model)
 Copy constructor. More...
 
virtual ~GModelSpatialDiffuse (void)
 Destructor. More...
 
virtual GModelSpatialDiffuseoperator= (const GModelSpatialDiffuse &model)
 Assignment operator. More...
 
virtual GClassCode code (void) const
 Return class code. 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...
 
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 GModelSpatialDiffuseCube &model)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void fetch_cube (void) const
 Fetch cube. More...
 
void load_cube (const GFilename &filename)
 Load map cube. More...
 
void set_energy_boundaries (void)
 Set energy boundaries. More...
 
void update_mc_cache (void)
 Update Monte Carlo cache. More...
 
double cube_intensity (const GPhoton &photon) const
 Compute cube intensity by log-log interpolation. More...
 
virtual void set_region (void) const
 Set boundary sky region. More...
 
- Protected Member Functions inherited from GModelSpatialDiffuse
void init_members (void)
 Initialise class members. More...
 
void copy_members (const GModelSpatialDiffuse &model)
 Copy class members. More...
 
void free_members (void)
 Delete class members. 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_value
 Value. More...
 
GFilename m_filename
 Name of map cube. More...
 
bool m_loaded
 Signals if map cube has been loaded. More...
 
GSkyMap m_cube
 Map cube. More...
 
GNodeArray m_logE
 Log10(energy) values of the maps. More...
 
GEbounds m_ebounds
 Energy bounds of the maps. More...
 
GSkyRegionCircle m_mc_cone
 MC cone. More...
 
double m_mc_one_minus_cosrad
 1-cosine of radius More...
 
std::vector< double > m_mc_max
 Maximum values for MC. More...
 
GModelSpectralNodes m_mc_spectrum
 Map cube spectrum. 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

Spatial map cube model.

This class implements the spatial component of the factorised source model for a map cube. A map cube is a set of sky maps for different energies. It is assumed that the pixels of the sky map are given in the units ph/cm2/s/sr/MeV.

Definition at line 59 of file GModelSpatialDiffuseCube.hpp.

Constructor & Destructor Documentation

GModelSpatialDiffuseCube::GModelSpatialDiffuseCube ( void  )

Void constructor.

Constructs empty map cube model.

Definition at line 81 of file GModelSpatialDiffuseCube.cpp.

References init_members().

Referenced by clone().

GModelSpatialDiffuseCube::GModelSpatialDiffuseCube ( const bool &  dummy,
const std::string &  type 
)

Model type constructor.

Parameters
[in]dummyDummy flag.
[in]typeModel type.

Constructs empty map cube model by specifying a model type.

Definition at line 100 of file GModelSpatialDiffuseCube.cpp.

References init_members(), GModelSpatial::m_type, and GModelSpatial::type().

GModelSpatialDiffuseCube::GModelSpatialDiffuseCube ( const GXmlElement xml)
explicit

XML constructor.

Parameters
[in]xmlXML element.

Constructs map cube 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 124 of file GModelSpatialDiffuseCube.cpp.

References init_members(), and read().

GModelSpatialDiffuseCube::GModelSpatialDiffuseCube ( const GFilename filename,
const double &  value = 1.0 
)

Filename constructor.

Parameters
[in]filenameFile name.
[in]valueNormalization factor (defaults to 1).

Constructs map cube model by loading a map cube from filename and by assigning the normalization value.

Definition at line 147 of file GModelSpatialDiffuseCube.cpp.

References GModelSpatial::autoscale(), filename(), init_members(), m_filename, m_value, and GOptimizerPar::value().

GModelSpatialDiffuseCube::GModelSpatialDiffuseCube ( const GSkyMap cube,
const GEnergies energies,
const double &  value = 1.0 
)

Sky map constructor.

Parameters
[in]cubeSky map cube.
[in]energiesSky map energies.
[in]valueNormalization factor (defaults to 1).

Constructs map cube model by extracting a cube from a sky map. The constructor also assigns the energy values for all maps and sets the scaling value. The filename will remain blank.

Definition at line 179 of file GModelSpatialDiffuseCube.cpp.

References GModelSpatial::autoscale(), cube(), energies(), init_members(), m_value, and GOptimizerPar::value().

GModelSpatialDiffuseCube::GModelSpatialDiffuseCube ( const GModelSpatialDiffuseCube model)

Copy constructor.

Parameters
[in]modelMap cube model.

Definition at line 209 of file GModelSpatialDiffuseCube.cpp.

References copy_members(), and init_members().

GModelSpatialDiffuseCube::~GModelSpatialDiffuseCube ( void  )
virtual

Destructor.

Definition at line 226 of file GModelSpatialDiffuseCube.cpp.

References free_members().

Member Function Documentation

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

Return class name.

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

Implements GModelSpatialDiffuse.

Definition at line 151 of file GModelSpatialDiffuseCube.hpp.

void GModelSpatialDiffuseCube::clear ( void  )
virtual
GModelSpatialDiffuseCube * GModelSpatialDiffuseCube::clone ( void  ) const
virtual

Clone map cube model.

Returns
Pointer to deep copy of map cube model.

Implements GModelSpatialDiffuse.

Definition at line 303 of file GModelSpatialDiffuseCube.cpp.

References GModelSpatialDiffuseCube().

bool GModelSpatialDiffuseCube::contains ( const GSkyDir dir,
const double &  margin = 0.0 
) const
virtual

Signals whether model contains sky direction.

Parameters
[in]dirSky direction.
[in]marginMargin to be added to sky direction (degrees)
Returns
True if the model contains the sky direction.

Signals whether a sky direction falls within the bounding circle of the diffuse cube.

Implements GModelSpatialDiffuse.

Definition at line 517 of file GModelSpatialDiffuseCube.cpp.

References GSkyRegionCircle::centre(), GSkyDir::dist_deg(), fetch_cube(), GModelSpatial::m_region, and GSkyRegionCircle::radius().

void GModelSpatialDiffuseCube::copy_members ( const GModelSpatialDiffuseCube model)
protected

Copy class members.

Parameters
[in]modelSpatial map cube model.

Definition at line 1220 of file GModelSpatialDiffuseCube.cpp.

References m_cube, m_ebounds, m_filename, m_loaded, m_logE, m_mc_cone, m_mc_max, m_mc_one_minus_cosrad, m_mc_spectrum, GModelSpatial::m_pars, GModelSpatial::m_type, and m_value.

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

const GSkyMap & GModelSpatialDiffuseCube::cube ( void  ) const
inline

Get map cube.

Returns
Map cube.

Returns the map cube.

Definition at line 223 of file GModelSpatialDiffuseCube.hpp.

References fetch_cube(), and m_cube.

Referenced by cube(), GModelSpatialDiffuseCube(), GResponse::irf_diffuse(), mc(), mc_cone(), and gammalib::resolution().

void GModelSpatialDiffuseCube::cube ( const GSkyMap cube)

Set map cube.

Parameters
[in]cubeSky map.

Set the map cube of the spatial map cube model.

Definition at line 646 of file GModelSpatialDiffuseCube.cpp.

References GFilename::clear(), cube(), m_cube, m_filename, m_loaded, and update_mc_cache().

double GModelSpatialDiffuseCube::cube_intensity ( const GPhoton photon) const
protected

Compute cube intensity by log-log interpolation.

Parameters
[in]photonIncident photon.
Returns
Cube intensity value.

Definition at line 1415 of file GModelSpatialDiffuseCube.cpp.

References GPhoton::dir(), GPhoton::energy(), exp(), fetch_cube(), GNodeArray::inx_left(), GNodeArray::inx_right(), log(), GEnergy::log10MeV(), m_cube, m_logE, GNodeArray::set_value(), GNodeArray::size(), GNodeArray::wgt_left(), and GNodeArray::wgt_right().

Referenced by eval().

GEnergies GModelSpatialDiffuseCube::energies ( void  )

Returns map cube energies.

Returns
Map cube energies.

Returns the energies for the map cube in a vector.

Definition at line 717 of file GModelSpatialDiffuseCube.cpp.

References GEnergies::append(), fetch_cube(), GEnergy::log10MeV(), m_logE, GEnergies::reserve(), and GNodeArray::size().

Referenced by GModelSpatialDiffuseCube(), load_cube(), read(), and write().

void GModelSpatialDiffuseCube::energies ( const GEnergies energies)

Set energies for map cube.

Parameters
[in]energiesSky map energies.
Exceptions
GException::invalid_argumentSpecified sky map energies incompatible with map cube.

Sets the energies for the map cube.

Definition at line 673 of file GModelSpatialDiffuseCube.cpp.

References GNodeArray::append(), GNodeArray::clear(), fetch_cube(), G_ENERGIES, m_cube, m_logE, GSkyMap::nmaps(), set_energy_boundaries(), GEnergies::size(), gammalib::str(), and update_mc_cache().

double GModelSpatialDiffuseCube::eval ( const GPhoton photon,
const bool &  gradients = false 
) const
virtual

Evaluate function.

Parameters
[in]photonIncident photon.
[in]gradientsCompute gradients?
Returns
Sky map intensity ( \(\mbox{ph cm}^{-2}\mbox{sr}^{-1}\mbox{s}^{-1}\))

Computes the spatial diffuse model as function of photon parameters.

If the gradients flag is true the method will also evaluate the partial derivatives of the model.

Implements GModelSpatialDiffuse.

Definition at line 322 of file GModelSpatialDiffuseCube.cpp.

References cube_intensity(), GOptimizerPar::factor_gradient(), GOptimizerPar::is_free(), m_value, GOptimizerPar::scale(), value(), and GOptimizerPar::value().

Referenced by GResponse::irf_diffuse().

void GModelSpatialDiffuseCube::fetch_cube ( void  ) const
protected

Fetch cube.

Load diffuse cube if it is not yet loaded. The loading is thread save.

Definition at line 1261 of file GModelSpatialDiffuseCube.cpp.

References GFilename::is_empty(), load_cube(), m_filename, and m_loaded.

Referenced by contains(), cube(), cube_intensity(), energies(), flux(), mc(), mc_cone(), set_region(), spectrum(), and write().

const GFilename & GModelSpatialDiffuseCube::filename ( void  ) const
inline

Get file name.

Returns
File name.

Returns the file name of the spatial map cube model.

Definition at line 194 of file GModelSpatialDiffuseCube.hpp.

References m_filename.

Referenced by filename(), GModelSpatialDiffuseCube(), and load_cube().

void GModelSpatialDiffuseCube::filename ( const GFilename filename)
inline

Set file name.

Parameters
[in]filenameFile name.

Set the file name of the spatial map cube model.

Definition at line 208 of file GModelSpatialDiffuseCube.hpp.

References filename(), and m_filename.

double GModelSpatialDiffuseCube::flux ( const GSkyRegion region,
const GEnergy srcEng = GEnergy(),
const GTime srcTime = GTime() 
) const
virtual

Returns diffuse cube flux integrated in sky region.

Parameters
[in]regionSky region.
[in]srcEngEnergy.
[in]srcTimeTime.
Returns
Flux (adimensional or ph/cm2/s).

Returns diffuse cube flux within a sky region. The flux is computed by log-log interpolation of the flux of the sky maps that encompass the specified srcEng energy. Specifically

\[ F = {\tt value} \times \left \{ \begin{array}{l l} 0 & \mbox{if $F_l=0$ and $F_r=0$} \\ F_l & \mbox{if $F_l\ne0$ and $F_r=0$} \\ F_r & \mbox{if $F_l=0$ and $F_r\ne0$} \\ \exp \left( w_l \log F_l + w_r \log F_r \right) & \mbox{else} \end{array} \right . \]

where

  • \(F\) is the diffuse cube flux
  • \({\tt value}\) is the cube normalisation that is returned by the value() method,
  • \(w_l\) and \(w_r\) are the weighting factors from the logarithmic interpolation in srcEng, with \(w_l+w_r=1\), and
  • \(F_l\) and \(F_r\) are the fluxes in the sky maps that encompass the srcEng energy.

This method fetches the diffuse cube if it was not yet fetched before.

Reimplemented from GModelSpatial.

Definition at line 1042 of file GModelSpatialDiffuseCube.cpp.

References exp(), fetch_cube(), GSkyMap::flux(), GNodeArray::inx_left(), GNodeArray::inx_right(), log(), GEnergy::log10MeV(), m_cube, m_logE, m_value, GNodeArray::set_value(), GNodeArray::size(), GOptimizerPar::value(), GNodeArray::wgt_left(), and GNodeArray::wgt_right().

Referenced by mc_cone().

void GModelSpatialDiffuseCube::free_members ( void  )
protected

Delete class members.

Definition at line 1249 of file GModelSpatialDiffuseCube.cpp.

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

void GModelSpatialDiffuseCube::load ( const GFilename filename)

Load cube into the model class.

Parameters
[in]filenamecube file.

Loads cube into the model class.

Definition at line 879 of file GModelSpatialDiffuseCube.cpp.

References load_cube(), and update_mc_cache().

void GModelSpatialDiffuseCube::load_cube ( const GFilename filename)
protected

Load map cube.

Parameters
[in]filenameMap cube file.
Exceptions
GException::invalid_valueNumber of maps in cube mismatches number of energy bins.

Load diffuse map cube.

Definition at line 1289 of file GModelSpatialDiffuseCube.cpp.

References GNodeArray::append(), GSkyRegionCircle::clear(), GNodeArray::clear(), GSkyMap::clear(), energies(), filename(), G_LOAD_CUBE, GSkyMap::load(), m_cube, m_filename, m_loaded, m_logE, GModelSpatial::m_region, mc_cone(), GSkyMap::nmaps(), GSkyMap::region_circle(), set_energy_boundaries(), GEnergies::size(), and gammalib::str().

Referenced by fetch_cube(), and load().

GSkyDir GModelSpatialDiffuseCube::mc ( const GEnergy energy,
const GTime time,
GRan ran 
) const
virtual

Returns MC sky direction.

Parameters
[in]energyPhoton energy.
[in]timePhoton arrival time.
[in,out]ranRandom number generator.
Returns
Sky direction.
Exceptions
GException::invalid_valueNo map cube defined. No energy boundaries defined.
GException::invalid_return_valueSimulation cone not defined, does not overlap with map cube or map cube is empty for the specified energy.

Returns a random sky direction according to the intensity distribution of the model sky map and the specified energy. The method uses a rejection method to determine the sky direction. If no sky direction could be determined, the method throws an GException::invalid_return_value exception.

Todo:
Make sure that spatial model value is taken into account.

Implements GModelSpatialDiffuse.

Definition at line 379 of file GModelSpatialDiffuseCube.cpp.

References gammalib::acos(), GSkyRegionCircle::centre(), cube(), exp(), fetch_cube(), G_MC, GNodeArray::inx_left(), GNodeArray::inx_right(), log(), GEnergy::log10MeV(), m_cube, m_ebounds, m_logE, m_mc_cone, m_mc_max, m_mc_one_minus_cosrad, max(), GSkyMap::npix(), GEnergy::print(), gammalib::rad2deg, GSkyDir::rotate_deg(), GNodeArray::set_value(), GEbounds::size(), GRan::uniform(), value(), GNodeArray::wgt_left(), and GNodeArray::wgt_right().

void GModelSpatialDiffuseCube::mc_cone ( const GSkyRegionCircle cone) const

Set Monte Carlo simulation cone.

Parameters
[in]coneMonte Carlo simulation cone.

Sets the simulation cone that defines the directions that will be simulated using the mc() method and pre-computes the maximum intensity and the spatially integrated flux of each map within the simulation cone region.

Definition at line 758 of file GModelSpatialDiffuseCube.cpp.

References GModelSpectralNodes::append(), GSkyRegionCircle::centre(), GModelSpectralNodes::clear(), cos(), cube(), gammalib::deg2rad, GSkyDir::dist_deg(), GModelSpectralNodes::energy(), fetch_cube(), flux(), GSkyMap::flux(), GModelSpectralNodes::intensity(), GEnergy::log10MeV(), m_cube, m_logE, m_mc_cone, m_mc_max, m_mc_one_minus_cosrad, m_mc_spectrum, GSkyMap::nmaps(), GModelSpectralNodes::nodes(), GSkyMap::npix(), GSkyMap::pix2dir(), GSkyRegionCircle::radius(), GModelSpectralNodes::reserve(), GNodeArray::size(), and value().

Referenced by GModelSky::mc().

const GSkyRegionCircle & GModelSpatialDiffuseCube::mc_cone ( void  ) const
inline

Get Monte Carlo simulation cone.

Returns
Monte Carlo simulation sky region circle.

Returns the sky region circle used for Monte Carlo simulations.

Definition at line 271 of file GModelSpatialDiffuseCube.hpp.

References m_mc_cone.

Referenced by load_cube(), and update_mc_cache().

double GModelSpatialDiffuseCube::mc_norm ( const GSkyDir dir,
const double &  radius 
) const
inlinevirtual

Return normalization of diffuse cube for Monte Carlo simulations.

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

Returns the normalization of a diffuse cube. The normalization is given by the model value. The dir and radius parameters are not used.

Implements GModelSpatialDiffuse.

Definition at line 256 of file GModelSpatialDiffuseCube.hpp.

References value().

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

Assignment operator.

Parameters
[in]modelMap cube model.
Returns
Map cube model.

Definition at line 248 of file GModelSpatialDiffuseCube.cpp.

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

std::string GModelSpatialDiffuseCube::print ( const GChatter chatter = NORMAL) const
virtual
void GModelSpatialDiffuseCube::read ( const GXmlElement xml)
virtual

Read model from XML element.

Parameters
[in]xmlXML element.
Exceptions
GException::invalid_valueModel parameters not found in XML element.

Read the map cube information from an XML element. The XML element should have the format

<spatialModel type="DiffuseMapCube" file="test_file.fits">
  <parameter name="Normalization" ../>
</spatialModel>

Implements GModelSpatialDiffuse.

Definition at line 560 of file GModelSpatialDiffuseCube.cpp.

References GXmlElement::attribute(), G_READ, m_filename, m_value, GModelPar::read(), gammalib::xml_check_parnum(), gammalib::xml_file_expand(), gammalib::xml_get_par(), and gammalib::xml_has_par().

Referenced by GModelSpatialDiffuseCube().

void GModelSpatialDiffuseCube::read ( const GFits fits)

Read diffuse cube from FITS file.

Parameters
[in]fitsFITS file.
Exceptions
GException::invalid_valueENERGIES extension incompatible with map cube.

Reads the diffuse cube sky map from the first extension in the fits file and the energy extension from the "ENERGIES" extension.

Definition at line 928 of file GModelSpatialDiffuseCube.cpp.

References GNodeArray::append(), GFits::at(), GNodeArray::clear(), GSkyMap::clear(), energies(), gammalib::extname_energies, GFitsHDU::exttype(), G_READ_FITS, GFitsHDU::has_card(), GFitsHDU::HT_IMAGE, GFitsHDU::integer(), m_cube, m_loaded, m_logE, GSkyMap::nmaps(), GEnergies::read(), GSkyMap::read(), set_energy_boundaries(), GEnergies::size(), GFits::size(), gammalib::str(), and GFits::table().

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

Save cube into FITS file.

Parameters
[in]filenameCube FITS file name.
[in]clobberOverwrite existing FITS file (default: false).

Saves spatial cube model into a FITS file.

Definition at line 900 of file GModelSpatialDiffuseCube.cpp.

References GFits::saveto(), and write().

void GModelSpatialDiffuseCube::set_energy_boundaries ( void  )
protected

Set energy boundaries.

Computes energy boundaries from the energy nodes. The boundaries will be computed as the centre in log10(energy) between the nodes. If only a single map exists, no boundaries will be computed.

Definition at line 1350 of file GModelSpatialDiffuseCube.cpp.

References GEbounds::append(), GEbounds::clear(), GEnergy::log10MeV(), m_ebounds, m_logE, and GNodeArray::size().

Referenced by energies(), load_cube(), and read().

void GModelSpatialDiffuseCube::set_region ( void  ) const
protectedvirtual

Set boundary sky region.

Implements GModelSpatialDiffuse.

Definition at line 1456 of file GModelSpatialDiffuseCube.cpp.

References fetch_cube().

const GModelSpectralNodes & GModelSpatialDiffuseCube::spectrum ( void  ) const
inline

Get map cube spectrum.

Returns
Map cube spectrum.

Returns the map cube spectrum.

Definition at line 238 of file GModelSpatialDiffuseCube.hpp.

References fetch_cube(), and m_mc_spectrum.

Referenced by GModelSky::eflux(), GModelSky::flux(), and GModelSky::mc().

void GModelSpatialDiffuseCube::update_mc_cache ( void  )
protected

Update Monte Carlo cache.

Initialise the cache for Monte Carlo sampling of the map cube. See the mc_cone() for more information.

Definition at line 1396 of file GModelSpatialDiffuseCube.cpp.

References mc_cone().

Referenced by cube(), energies(), and load().

double GModelSpatialDiffuseCube::value ( void  ) const
inline

Get model value.

Returns
Model value.

Returns the value of the spatial map cube model.

Definition at line 165 of file GModelSpatialDiffuseCube.hpp.

References m_value, and GOptimizerPar::value().

Referenced by eval(), mc(), mc_cone(), mc_norm(), and write().

void GModelSpatialDiffuseCube::value ( const double &  value)
inline

Set model value.

Parameters
[in]valueModel value.

Set the value of the spatial map cube model.

Definition at line 179 of file GModelSpatialDiffuseCube.hpp.

References m_value, and GOptimizerPar::value().

void GModelSpatialDiffuseCube::write ( GXmlElement xml) const
virtual

Write model into XML element.

Parameters
[in]xmlXML element.

Write the map cube information into an XML element. The XML element will have either the format

<spatialModel type="DiffuseMapCube" file="test_file.fits">
  <parameter name="Normalization" scale="1" value="1" min="0.1" max="10" free="0"/>
</spatialModel>

or alternatively

<spatialModel type="DiffuseMapCube" file="test_file.fits">
  <parameter name="Value" scale="1" value="1" min="0.1" max="10" free="0"/>
</spatialModel>

The latter format is the default for newly written XML elements.

Implements GModelSpatialDiffuse.

Definition at line 620 of file GModelSpatialDiffuseCube.cpp.

References GXmlElement::attribute(), G_WRITE, m_filename, m_value, GOptimizerPar::name(), GModelSpatial::type(), value(), GModelPar::write(), gammalib::xml_check_type(), gammalib::xml_file_reduce(), and gammalib::xml_need_par().

Referenced by save().

void GModelSpatialDiffuseCube::write ( GFits fits) const

Write diffuse cube into FITS file.

Parameters
[in]fitsFITS file.

Writes the diffuse cube sky map and the energy intervals into fits file.

Definition at line 985 of file GModelSpatialDiffuseCube.cpp.

References GEnergies::append(), energies(), fetch_cube(), m_cube, m_logE, pow(), GNodeArray::size(), GEnergies::write(), and GSkyMap::write().

Member Data Documentation

GSkyMap GModelSpatialDiffuseCube::m_cube
protected
GEbounds GModelSpatialDiffuseCube::m_ebounds
protected

Energy bounds of the maps.

Definition at line 135 of file GModelSpatialDiffuseCube.hpp.

Referenced by copy_members(), init_members(), mc(), print(), and set_energy_boundaries().

GFilename GModelSpatialDiffuseCube::m_filename
protected
bool GModelSpatialDiffuseCube::m_loaded
protected

Signals if map cube has been loaded.

Definition at line 132 of file GModelSpatialDiffuseCube.hpp.

Referenced by copy_members(), cube(), fetch_cube(), init_members(), load_cube(), print(), and read().

GNodeArray GModelSpatialDiffuseCube::m_logE
protected

Log10(energy) values of the maps.

Definition at line 134 of file GModelSpatialDiffuseCube.hpp.

Referenced by copy_members(), cube_intensity(), energies(), flux(), init_members(), load_cube(), mc(), mc_cone(), print(), read(), set_energy_boundaries(), and write().

GSkyRegionCircle GModelSpatialDiffuseCube::m_mc_cone
mutableprotected

MC cone.

Definition at line 138 of file GModelSpatialDiffuseCube.hpp.

Referenced by copy_members(), init_members(), mc(), and mc_cone().

std::vector<double> GModelSpatialDiffuseCube::m_mc_max
mutableprotected

Maximum values for MC.

Definition at line 140 of file GModelSpatialDiffuseCube.hpp.

Referenced by copy_members(), init_members(), mc(), and mc_cone().

double GModelSpatialDiffuseCube::m_mc_one_minus_cosrad
mutableprotected

1-cosine of radius

Definition at line 139 of file GModelSpatialDiffuseCube.hpp.

Referenced by copy_members(), init_members(), mc(), and mc_cone().

GModelSpectralNodes GModelSpatialDiffuseCube::m_mc_spectrum
mutableprotected

Map cube spectrum.

Definition at line 141 of file GModelSpatialDiffuseCube.hpp.

Referenced by copy_members(), init_members(), mc_cone(), print(), and spectrum().

GModelPar GModelSpatialDiffuseCube::m_value
protected

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