ctools  2.0.0
 All Classes Namespaces Files Functions Variables Macros Pages
ctskymap Class Reference

Sky mapping tool. More...

#include <ctskymap.hpp>

Inheritance diagram for ctskymap:
ctobservation ctool GApplication

Public Member Functions

 ctskymap (void)
 Void constructor. More...
 
 ctskymap (const GObservations &obs)
 Observations constructor. More...
 
 ctskymap (int argc, char *argv[])
 Command line constructor. More...
 
 ctskymap (const ctskymap &app)
 Copy constructor. More...
 
virtual ~ctskymap (void)
 Destructor. More...
 
ctskymapoperator= (const ctskymap &app)
 Assignment operator. More...
 
void clear (void)
 Clear sky mapping tool. More...
 
void process (void)
 Process the sky mapping tool. More...
 
void save (void)
 Save sky map. More...
 
void publish (const std::string &name="")
 Publish sky map. More...
 
const GSkyMap & skymap (void) const
 Return observation container. More...
 
const GFits & fits (void) const
 Return fits container. More...
 
void exclusion_map (const GSkyRegionMap &exclusion_map)
 Set exclusion region map. More...
 
GSkyRegionMap exclusion_map (void) const
 Return exclusion region map. More...
 
- Public Member Functions inherited from ctobservation
 ctobservation (const std::string &name, const std::string &version)
 Name constructor. More...
 
 ctobservation (const std::string &name, const std::string &version, const GApplicationPars &pars)
 Application parameters constructor. More...
 
 ctobservation (const std::string &name, const std::string &version, int argc, char *argv[])
 Command line constructor. More...
 
 ctobservation (const std::string &name, const std::string &version, const GObservations &obs)
 Observations constructor. More...
 
 ctobservation (const ctobservation &app)
 Copy constructor. More...
 
virtual ~ctobservation (void)
 Destructor. More...
 
ctobservationoperator= (const ctobservation &app)
 Assignment operator. More...
 
void obs (const GObservations &obs)
 Set observation container. More...
 
const GObservations & obs (void) const
 Return observation container. More...
 
- Public Member Functions inherited from ctool
 ctool (const std::string &name, const std::string &version)
 Name constructor. More...
 
 ctool (const std::string &name, const std::string &version, const GApplicationPars &pars)
 Application parameter constructor. More...
 
 ctool (const std::string &name, const std::string &version, int argc, char *argv[])
 Command line constructor. More...
 
 ctool (const ctool &app)
 Copy constructor. More...
 
virtual ~ctool (void)
 Destructor. More...
 
ctooloperator= (const ctool &app)
 Assignment operator. More...
 
virtual void run (void)
 Run ctool. More...
 
virtual void execute (void)
 Execute ctool. More...
 

Protected Member Functions

void init_members (void)
 Initialise class members. More...
 
void copy_members (const ctskymap &app)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void get_parameters (void)
 Get application parameters. More...
 
void setup_maps (void)
 Setup maps. More...
 
void setup_exclusion_map (void)
 Generates map of pixel exclusions. More...
 
void setup_exclusion_map_fits (const GFilename &filename)
 Fills exclusions map from FITS image. More...
 
void setup_exclusion_map_region (const GFilename &filename)
 Fills exclusions map from DS9 region file. More...
 
void adjust_exclusion_map (void)
 Verifys that exclusion map points in fov of counts. More...
 
void fill_maps (void)
 Fill maps from observation container. More...
 
void fill_maps_counts (GCTAObservation *obs)
 Fill events into counts map. More...
 
void fill_maps_acceptance (GCTAObservation *obs)
 Compute background acceptance sky map based on IRF template. More...
 
void compute_maps (void)
 Compute sky map, background map and significance map. More...
 
void compute_maps_ring_fft (void)
 Compute the maps for RING background using a FFT. More...
 
void compute_maps_ring_direct (void)
 Compute the pixel significance for RING background. More...
 
void compute_ring_values (const int &ipixel, const GSkyMap &counts, const GSkyMap &background, double &non, double &noff, double &alpha)
 Computes Non, Noff and alpha for a counts map and sensitivity map. More...
 
void construct_fits (void)
 Construct GFits object consisting of all maps. More...
 
void ring_bounding_box (const int &ipixel, int &ix1, int &ix2, int &iy1, int &iy2) const
 Computes bounding box for RING background computation. More...
 
GSkyMap ring_convolve (const GSkyMap &map, const double &rmin, const double &rmax) const
 Return FFT kernel for background ring. More...
 
GNdarray ring_kernel (const double &rmin, const double &rmax) const
 Return FFT kernel for background ring. More...
 
double sigma_li_ma (const double &n_on, const double &n_off, const double &alpha) const
 Compute significance following Li & Ma. More...
 
void write_map (GFits &fits, const GSkyMap &map, const std::string &extname) const
 Write sky map into FITS file. More...
 
void write_hdu_keywords (GFitsHDU *hdu) const
 Write keywords in FITS HDU. More...
 
- Protected Member Functions inherited from ctobservation
GCTAObservation * first_unbinned_observation (void)
 Return first unbinned CTA observation. More...
 
GCTAObservation * next_unbinned_observation (void)
 Return next unbinned CTA observation. More...
 
const GCTAObservation * first_unbinned_observation (void) const
 Return first unbinned CTA observation (const version) More...
 
const GCTAObservation * next_unbinned_observation (void) const
 Return next unbinned CTA observation (const version) More...
 
void read_ogip_keywords (GFitsHDU *hdu) const
 Read OGIP keywords from FITS HDU. More...
 
void write_ogip_keywords (GFitsHDU *hdu) const
 Write OGIP keywords in FITS HDU. More...
 
void set_obs_statistic (const std::string &statistic)
 Set fit statistic for CTA observations. More...
 
void set_obs_bounds ()
 Set observation boundaries for CTA observations. More...
 
void save_events_fits (void)
 Save event list in FITS format. More...
 
void save_events_xml (void)
 Save event list(s) in XML format. More...
 
void init_members (void)
 Initialise class members. More...
 
void copy_members (const ctobservation &app)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
- Protected Member Functions inherited from ctool
void init_members (void)
 Initialise class members. More...
 
void copy_members (const ctool &app)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void sync_pfiles (void)
 Synchronise parameter files. More...
 
const bool & read_ahead (void) const
 Signal whether parameters should be read ahead. More...
 
void setup_observations (GObservations &obs, const bool &response=true, const bool &list=true, const bool &cube=true)
 Setup observation container. More...
 
void setup_models (GObservations &obs, const std::string &name="")
 Setup model container. More...
 
GEbounds create_ebounds (void)
 Create energy boundaries from user parameters. More...
 
GEnergies create_energies (void)
 Create energies from user parameters. More...
 
GSkyMap create_map (const GObservations &obs)
 Create a skymap from user parameters. More...
 
GCTAEventCube create_cube (const GObservations &obs)
 Create a CTA event cube from user parameters. More...
 
GCTAObservation create_cta_obs (void)
 Create a CTA observation from User parameters. More...
 
void require_inobs (const std::string &method)
 Throws exception if inobs parameter is not valid. More...
 
void require_inobs_nolist (const std::string &method)
 
void require_inobs_nocube (const std::string &method)
 Throws exception if inobs parameter is a counts cube. More...
 
GCTARoi get_roi (const GCTAPointing &pnt=GCTAPointing())
 Return RoI from User parameters. More...
 
GEbounds get_ebounds (void)
 Return energy boundaries from User parameters. More...
 
GGti get_gti (const GTimeReference &ref)
 Return Good Time Intervals from User parameter. More...
 
GCTAPointing get_pointing (void)
 Return CTA pointing from User parameters. More...
 
GSkyDir get_skydir (void)
 Return sky direction from User parameters. More...
 
std::string set_outfile_name (const std::string &filename)
 Set output file name. More...
 
bool is_stacked (void)
 Query user parameters for stacked analysis. More...
 
bool is_onoff (void)
 Query user parameters for On/Off analysis. More...
 
void log_parameters (const GChatter &chatter)
 Log application parameters. More...
 
void log_observations (const GChatter &chatter, const GObservations &obs, const std::string &what="Observation")
 Log observation container. More...
 
void log_models (const GChatter &chatter, const GModels &models, const std::string &what="Model")
 Log model container. More...
 
void set_response (GObservations &obs)
 Set response for all CTA observations in container. More...
 
std::vector< bool > set_edisp (GObservations &obs, const bool &edisp) const
 Set energy dispersion to CTA observations. More...
 
void restore_edisp (GObservations &obs, const std::vector< bool > &edisp) const
 Restore energy dispersion flags of CTA observations. More...
 
void set_obs_response (GCTAObservation *obs)
 Set response for CTA observation. More...
 
GObservations get_observations (const bool &get_response=true)
 Get observation container. More...
 
GSkyDir get_mean_pointing (const GObservations &obs)
 Derives mean pointing from CTA observations. More...
 
size_t get_current_rss (void)
 Get current resident set size (physical memory use) in Bytes. More...
 
std::string get_obs_header (const GObservation *obs) const
 Return observation header string. More...
 
GEnergies insert_energy_boundaries (const GEnergies &energies, const GCTAObservation &obs)
 Insert observation energy boundaries into list of energies. More...
 
std::vector< bool > cube_layer_usage (const GEbounds &cube_ebounds, const GEbounds &list_ebounds) const
 Determine the counts cube layer usage. More...
 
std::string get_gtiname (const std::string &filename, const std::string &evtname) const
 Get Good Time Intervals extension name. More...
 
void save_event_list (const GCTAObservation *obs, const std::string &infile, const std::string &evtname, const std::string &gtiname, const std::string &outfile) const
 Save event list into FITS file. More...
 
std::string warn_too_few_energies (const GEnergies &energies) const
 Set warning string if there are too few energies. More...
 
std::string warn_xml_suffix (const GFilename &filename) const
 Set warning string if file has no .xml suffix. More...
 
void provide_help (void) const
 Dump help text in the console. More...
 

Protected Attributes

GFilename m_outmap
 Output file name. More...
 
GFilename m_inexclusion
 Exclusion map file name. More...
 
double m_emin
 Minimum energy (TeV) More...
 
double m_emax
 Maximum energy (TeV) More...
 
std::string m_bkgsubtract
 Background subtraction method. More...
 
double m_roiradius
 Region of interest radius for RING bkg. More...
 
double m_inradius
 Inner ring radius for RING background. More...
 
double m_outradius
 Outer ring radius for RING background. More...
 
int m_iterations
 Number of iterations for RING background. More...
 
double m_threshold
 Threshold for RING background. More...
 
bool m_usefft
 Use FFT for RING background. More...
 
bool m_publish
 Publish sky map? More...
 
GChatter m_chatter
 Chattiness. More...
 
bool m_has_inmap
 Has valid input map. More...
 
GSkyMap m_skymap
 Sky map. More...
 
GSkyMap m_bkgmap
 Background map. More...
 
GSkyMap m_sigmap
 Significance map. More...
 
GSkyMap m_counts
 Counts map. More...
 
GSkyMap m_acceptance
 Acceptance map. More...
 
GSkyMap m_exclmap
 Exclusion map for RING background. More...
 
GFits m_fits
 Output GFits object. More...
 
std::vector< double > m_solidangle
 Cached pixel solid angles. More...
 
std::vector< GSkyDir > m_dirs
 Cached pixel directions. More...
 
double m_cos_roiradius
 Cosine of RoI radius. More...
 
double m_cos_inradius
 Cosine of inner ring radius. More...
 
double m_cos_outradius
 Cosine of outer ring radius. More...
 
- Protected Attributes inherited from ctobservation
GObservations m_obs
 Observation container. More...
 
- Protected Attributes inherited from ctool
bool m_read_ahead
 Read ahead output parameters. More...
 
bool m_use_xml
 Use XML file instead of FITS file for observations. More...
 

Detailed Description

Sky mapping tool.

This class creates a sky map from a CTA event list.

Definition at line 45 of file ctskymap.hpp.

Constructor & Destructor Documentation

ctskymap::ctskymap ( void  )

Void constructor.

Constructs an empty sky mapping tool.

Definition at line 59 of file ctskymap.cpp.

References init_members().

ctskymap::ctskymap ( const GObservations &  obs)
explicit

Observations constructor.

param[in] obs Observation container.

Constructs sky mapping tool from an observation container.

Definition at line 76 of file ctskymap.cpp.

References init_members().

ctskymap::ctskymap ( int  argc,
char *  argv[] 
)

Command line constructor.

Parameters
[in]argcNumber of arguments in command line.
[in]argvArray of command line arguments.

Constructs sky mapping tool using command line arguments for user parameter setting.

Definition at line 97 of file ctskymap.cpp.

References init_members().

ctskymap::ctskymap ( const ctskymap app)

Copy constructor.

Parameters
[in]appSky mapping tool.

Definition at line 113 of file ctskymap.cpp.

References copy_members(), and init_members().

ctskymap::~ctskymap ( void  )
virtual

Destructor.

Definition at line 129 of file ctskymap.cpp.

References free_members().

Member Function Documentation

void ctskymap::adjust_exclusion_map ( void  )
protected

Verifys that exclusion map points in fov of counts.

Adjust exclusion map. Reinitialise exclusion map via the counts object and transfer the masked areas provided by the current exclusion map to the updated fov.

Definition at line 721 of file ctskymap.cpp.

References m_counts, and m_exclmap.

Referenced by setup_exclusion_map_fits(), and setup_maps().

void ctskymap::clear ( void  )
virtual

Clear sky mapping tool.

Clears sky mapping tool.

Implements ctobservation.

Definition at line 186 of file ctskymap.cpp.

References free_members(), ctool::free_members(), ctobservation::free_members(), init_members(), ctool::init_members(), and ctobservation::init_members().

void ctskymap::compute_maps ( void  )
protected

Compute sky map, background map and significance map.

Computes the sky map, background map and significance map. The following background subtraction methods are supported:

NONE

The sky map is simply the binned counts map, no background and significance maps are computed.

IRF

The background map is the acceptance map, and this background map is subtracted from the counts map and assigned to the sky map. The significance \(\sigma_i\) is computed assuming Poisson statistic in the Gaussian limit, using

\[\sigma_i = \frac{N_i - B_i}{\sqrt{N_i}}\]

where \(N_i\) is the number of observed counts and \(B_i\) is the estimated number of background counts.

RING

For the RING method the Li & Ma significance is computed for each pixel, using the roiradius, inradius and outradius parameters to specify the radius of the On region and the Off background ring, respectively. Depending on the usefft parameter, either an FFT is used to compute the number of events and the acceptance of the On and Off regions, or a direct computation is performed. While the former is faster it is less accurate since it assumes Euclidean distances in the sky map. The latter is slower but compute exact distance between pixels of the sky map. Use FFT if the sky map is close to a cartesian grid, and use the direct method if the sky map shows important distortions.

Definition at line 1070 of file ctskymap.cpp.

References compute_maps_ring_direct(), compute_maps_ring_fft(), m_acceptance, m_bkgmap, m_bkgsubtract, m_counts, m_exclmap, m_iterations, m_sigmap, m_skymap, m_threshold, and m_usefft.

Referenced by process().

void ctskymap::compute_maps_ring_direct ( void  )
protected

Compute the pixel significance for RING background.

Computes the Li & Ma significance for each sky map pixel and replaces the sky and background maps by the On- and Off-count maps. The computation is done by computing the exact distances between pixels in the sky map.

Definition at line 1223 of file ctskymap.cpp.

References compute_ring_values(), m_acceptance, m_bkgmap, m_counts, m_sigmap, m_skymap, and sigma_li_ma().

Referenced by compute_maps().

void ctskymap::compute_maps_ring_fft ( void  )
protected

Compute the maps for RING background using a FFT.

Computes the Li & Ma significance for each sky map pixel and replaces the sky and background maps by the On- and Off-count maps. The computation is done using a FFT.

Definition at line 1137 of file ctskymap.cpp.

References m_acceptance, m_bkgmap, m_counts, m_exclmap, m_inradius, m_outradius, m_roiradius, m_sigmap, m_skymap, ring_convolve(), and sigma_li_ma().

Referenced by compute_maps().

void ctskymap::compute_ring_values ( const int &  ipixel,
const GSkyMap &  counts,
const GSkyMap &  background,
double &  non,
double &  noff,
double &  alpha 
)
protected

Computes Non, Noff and alpha for a counts map and sensitivity map.

Parameters
[in]ipixelSky map pixel to consider.
[in]countsCounts map.
[in]backgroundBackground map.
[out]nonReturned estimate of On-counts.
[out]noffReturned estimate of Off-counts.
[out]alphaReturned estimate of alpha.

Computes the On- and Off-counts values at a given position from in counts map and the alpha parameter from the background map.

To speed-up the computations the method considers only the pixels within a bounding box that comprises the outer background ring radius. The method considers wrapping around of pixel indices in the Right Ascension / Galactic longitude direction.

If the alpha of the ring region is zero the values of non, noff, and alpha will be set to zero.

Definition at line 1323 of file ctskymap.cpp.

References m_cos_inradius, m_cos_outradius, m_cos_roiradius, m_dirs, m_exclmap, and ring_bounding_box().

Referenced by compute_maps_ring_direct().

void ctskymap::construct_fits ( void  )
protected

Construct GFits object consisting of all maps.

Assembles all of the skymaps into a GFits object. This object is saved when ctskymap::save() is called. the object can also be accessed by calling ctskymap::fits().

Definition at line 287 of file ctskymap.cpp.

References m_acceptance, m_bkgmap, m_bkgsubtract, m_counts, m_exclmap, m_fits, m_sigmap, m_skymap, and write_map().

Referenced by process().

void ctskymap::copy_members ( const ctskymap app)
protected
void ctskymap::exclusion_map ( const GSkyRegionMap &  exclusion_map)
inline

Set exclusion region map.

Parameters
[in]exclusion_mapExclusion region map.

Definition at line 169 of file ctskymap.hpp.

References m_exclmap.

GSkyRegionMap ctskymap::exclusion_map ( void  ) const
inline

Return exclusion region map.

Returns
Exclusion region map

Definition at line 185 of file ctskymap.hpp.

References m_exclmap.

void ctskymap::fill_maps ( void  )
protected

Fill maps from observation container.

Fill counts and optionally acceptance maps from data in the observation container.

Definition at line 759 of file ctskymap.cpp.

References fill_maps_acceptance(), fill_maps_counts(), ctobservation::first_unbinned_observation(), ctool::log_observations(), m_bkgsubtract, ctobservation::m_obs, ctobservation::next_unbinned_observation(), and ctobservation::obs().

Referenced by process().

void ctskymap::fill_maps_acceptance ( GCTAObservation *  obs)
protected

Compute background acceptance sky map based on IRF template.

Parameters
[in]obsCTA observation.
Exceptions
GException::invalid_valueNo response information available for observation. No background template available in instrument response function.

Computes a background acceptance sky map for a given observation and adds the estimate to the acceptance sky map.

If the response contains a background template, that template is used to compute the background acceptance map.

If no background template is available, the outcome depends on the background subtraction method. If IRF background subtraction is requested, an exception is thrown. For RING background subtraction, a constant background rate is assumed.

Definition at line 930 of file ctskymap.cpp.

References G_FILL_MAPS_ACCEPTANCE, ctool::get_obs_header(), m_acceptance, m_bkgsubtract, m_dirs, m_emax, m_emin, and m_solidangle.

Referenced by fill_maps().

void ctskymap::fill_maps_counts ( GCTAObservation *  obs)
protected

Fill events into counts map.

Parameters
[in]obsCTA observation.
Exceptions
GException::invalid_valueObservation does not contain an event list.

Fills the events found in a CTA events list into a sky map. The method adds the events to the m_counts member.

Definition at line 822 of file ctskymap.cpp.

References G_FILL_MAPS_COUNTS, ctool::get_obs_header(), m_chatter, m_counts, m_emax, and m_emin.

Referenced by fill_maps().

const GFits & ctskymap::fits ( void  ) const
inline

Return fits container.

Returns
Reference to fits container

Definition at line 157 of file ctskymap.hpp.

References m_fits.

Referenced by get_parameters().

void ctskymap::free_members ( void  )
protected

Delete class members.

Definition at line 450 of file ctskymap.cpp.

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

void ctskymap::get_parameters ( void  )
protected

Get application parameters.

Get all task parameters from parameter file or (if required) by querying the user. Most parameters are only required if no observation exists so far in the observation container. In this case, a single CTA observation will be added to the container, using the definition provided in the parameter file.

Definition at line 466 of file ctskymap.cpp.

References ctool::create_map(), fits(), G_GET_PARAMETERS, ctool::log_parameters(), m_acceptance, m_bkgsubtract, m_chatter, m_counts, m_emax, m_emin, m_has_inmap, m_inradius, m_iterations, ctobservation::m_obs, m_outradius, m_publish, m_roiradius, m_threshold, m_usefft, ctool::read_ahead(), ctobservation::read_ogip_keywords(), ctool::set_response(), and ctool::setup_observations().

Referenced by process().

ctskymap & ctskymap::operator= ( const ctskymap app)

Assignment operator.

Parameters
[in]appSky mapping tool.
Returns
Sky mapping tool.

Definition at line 151 of file ctskymap.cpp.

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

void ctskymap::process ( void  )
virtual

Process the sky mapping tool.

Generates a sky map from event list by looping over all unbinned CTA observation in the observation container and filling all events into a sky map.

Implements ctobservation.

Definition at line 216 of file ctskymap.cpp.

References compute_maps(), construct_fits(), fill_maps(), get_parameters(), m_has_inmap, m_publish, publish(), and setup_maps().

void ctskymap::publish ( const std::string &  name = "")

Publish sky map.

Parameters
[in]nameSky map name.

Definition at line 333 of file ctskymap.cpp.

References CTSKYMAP_NAME, and m_skymap.

Referenced by process().

void ctskymap::ring_bounding_box ( const int &  ipixel,
int &  ix1,
int &  ix2,
int &  iy1,
int &  iy2 
) const
protected

Computes bounding box for RING background computation.

Parameters
[in]ipixelSky map pixel to consider.
[out]ix1Index of first pixel in x.
[out]ix2Index after last pixel in x.
[out]iy1Index of first pixel in y.
[out]iy2Index after last pixel in y.

Computes the bounding box that contained the background ring for a specific pixel for the RING background method.

The method determines the local pixel scale for the requested pixel and draws a bounding box with 1.5 times the outer background ring radius around the pixel. In the x direction the pixel indices are unconstrained, and the client has to assure that the pixel value is comprised within the validity range (this is required to handle the longitude wrap around). In the y direction the pixel indices are constrained to [0,ny], where ny is the number of y pixels in the sky map.

Definition at line 1444 of file ctskymap.cpp.

References G_RING_BOUNDING_BOX, m_counts, m_dirs, and m_outradius.

Referenced by compute_ring_values().

GSkyMap ctskymap::ring_convolve ( const GSkyMap &  map,
const double &  rmin,
const double &  rmax 
) const
protected

Return FFT kernel for background ring.

Parameters
[in]mapSky map to be convolved with a background ring.
[in]rminMinimum ring radius (degrees).
[in]rmaxMaximum ring radius (degrees).
Returns
FFT kernel.

Computes the FFT kernel for a backgrund ring.

Definition at line 1535 of file ctskymap.cpp.

References ring_kernel().

Referenced by compute_maps_ring_fft().

GNdarray ctskymap::ring_kernel ( const double &  rmin,
const double &  rmax 
) const
protected

Return FFT kernel for background ring.

Parameters
[in]rminMinimum ring radius (degrees).
[in]rmaxMaximum ring radius (degrees).
Returns
FFT kernel.
Exceptions
GException::invalid_valueSky map is not a WCS projection.

Computes the FFT kernel for a background ring. The computation is done assuming that the sky map represents a cartesian grid, and distances are computed in that grid using Euclidean distances. The rmin and rmax arguments refer to a ring radius in these Euclidean distances.

Definition at line 1586 of file ctskymap.cpp.

References G_RING_KERNEL, and m_counts.

Referenced by ring_convolve().

void ctskymap::save ( void  )
virtual

Save sky map.

Saves the sky map into a FITS file. The FITS file name is specified by the outname parameter.

Implements ctobservation.

Definition at line 252 of file ctskymap.cpp.

References m_fits, m_outmap, and m_skymap.

void ctskymap::setup_exclusion_map ( void  )
protected

Generates map of pixel exclusions.

Generates a sky map of the pixels that are to be excluded from the background estimation. Pixels with values different from 0 will be excluded.

Definition at line 633 of file ctskymap.cpp.

References m_counts, m_exclmap, m_inexclusion, setup_exclusion_map_fits(), and setup_exclusion_map_region().

Referenced by setup_maps().

void ctskymap::setup_exclusion_map_fits ( const GFilename &  filename)
protected

Fills exclusions map from FITS image.

Parameters
[in]filenameFITS image file name.

Sets all exclusion map pixels to 1 that correspond to non-zero pixels in the exclusion sky map FITS file.

Definition at line 670 of file ctskymap.cpp.

References adjust_exclusion_map(), and m_exclmap.

Referenced by setup_exclusion_map().

void ctskymap::setup_exclusion_map_region ( const GFilename &  filename)
protected

Fills exclusions map from DS9 region file.

Parameters
[in]filenameDS9 region file name.

Sets all exclusion map pixels to 1 that are contained in any of the DS9 regions.

Definition at line 691 of file ctskymap.cpp.

References m_exclmap.

Referenced by setup_exclusion_map().

void ctskymap::setup_maps ( void  )
protected

Setup maps.

Setup maps for sky map generation.

Definition at line 586 of file ctskymap.cpp.

References adjust_exclusion_map(), m_bkgsubtract, m_cos_inradius, m_cos_outradius, m_cos_roiradius, m_counts, m_dirs, m_exclmap, m_inradius, m_outradius, m_roiradius, m_solidangle, and setup_exclusion_map().

Referenced by process().

double ctskymap::sigma_li_ma ( const double &  n_on,
const double &  n_off,
const double &  alpha 
) const
protected

Compute significance following Li & Ma.

Parameters
[in]n_onNumber of On-counts.
[in]n_offNumber of Off-counts.
[in]alphaAlpha parameters.
Returns
Significance.

Computes the significance following Li & Ma, Equation (17).

Definition at line 1640 of file ctskymap.cpp.

Referenced by compute_maps_ring_direct(), and compute_maps_ring_fft().

const GSkyMap & ctskymap::skymap ( void  ) const
inline

Return observation container.

Returns
Reference to observation container

Definition at line 145 of file ctskymap.hpp.

References m_skymap.

void ctskymap::write_hdu_keywords ( GFitsHDU *  hdu) const
protected

Write keywords in FITS HDU.

Parameters
[in,out]hduPointer to FITS HDU.

Writes keywords in FITS HDU.

Definition at line 1704 of file ctskymap.cpp.

References m_bkgsubtract, m_emax, m_emin, m_inexclusion, m_inradius, m_iterations, m_outradius, m_roiradius, m_threshold, and m_usefft.

Referenced by write_map().

void ctskymap::write_map ( GFits &  fits,
const GSkyMap &  map,
const std::string &  extname 
) const
protected

Write sky map into FITS file.

Parameters
[in,out]fitsFITS file.
[in]mapSky map.
[in]extnameExtension name.

Write one sky map with the extension name and keywords into the FITS file.

Definition at line 1678 of file ctskymap.cpp.

References write_hdu_keywords(), and ctobservation::write_ogip_keywords().

Referenced by construct_fits().

Member Data Documentation

GSkyMap ctskymap::m_acceptance
protected
GSkyMap ctskymap::m_bkgmap
protected
std::string ctskymap::m_bkgsubtract
protected
GChatter ctskymap::m_chatter
protected

Chattiness.

Definition at line 118 of file ctskymap.hpp.

Referenced by copy_members(), fill_maps_counts(), get_parameters(), and init_members().

double ctskymap::m_cos_inradius
protected

Cosine of inner ring radius.

Definition at line 134 of file ctskymap.hpp.

Referenced by compute_ring_values(), copy_members(), init_members(), and setup_maps().

double ctskymap::m_cos_outradius
protected

Cosine of outer ring radius.

Definition at line 135 of file ctskymap.hpp.

Referenced by compute_ring_values(), copy_members(), init_members(), and setup_maps().

double ctskymap::m_cos_roiradius
protected

Cosine of RoI radius.

Definition at line 133 of file ctskymap.hpp.

Referenced by compute_ring_values(), copy_members(), init_members(), and setup_maps().

std::vector<GSkyDir> ctskymap::m_dirs
protected

Cached pixel directions.

Definition at line 132 of file ctskymap.hpp.

Referenced by compute_ring_values(), copy_members(), fill_maps_acceptance(), init_members(), ring_bounding_box(), and setup_maps().

double ctskymap::m_emax
protected

Maximum energy (TeV)

Definition at line 109 of file ctskymap.hpp.

Referenced by copy_members(), fill_maps_acceptance(), fill_maps_counts(), get_parameters(), init_members(), and write_hdu_keywords().

double ctskymap::m_emin
protected

Minimum energy (TeV)

Definition at line 108 of file ctskymap.hpp.

Referenced by copy_members(), fill_maps_acceptance(), fill_maps_counts(), get_parameters(), init_members(), and write_hdu_keywords().

GFits ctskymap::m_fits
protected

Output GFits object.

Definition at line 128 of file ctskymap.hpp.

Referenced by construct_fits(), copy_members(), fits(), init_members(), and save().

bool ctskymap::m_has_inmap
protected

Has valid input map.

Definition at line 121 of file ctskymap.hpp.

Referenced by copy_members(), get_parameters(), init_members(), and process().

GFilename ctskymap::m_inexclusion
protected

Exclusion map file name.

Definition at line 107 of file ctskymap.hpp.

Referenced by copy_members(), init_members(), setup_exclusion_map(), and write_hdu_keywords().

double ctskymap::m_inradius
protected

Inner ring radius for RING background.

Definition at line 112 of file ctskymap.hpp.

Referenced by compute_maps_ring_fft(), copy_members(), get_parameters(), init_members(), setup_maps(), and write_hdu_keywords().

int ctskymap::m_iterations
protected

Number of iterations for RING background.

Definition at line 114 of file ctskymap.hpp.

Referenced by compute_maps(), copy_members(), get_parameters(), init_members(), and write_hdu_keywords().

GFilename ctskymap::m_outmap
protected

Output file name.

Definition at line 106 of file ctskymap.hpp.

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

double ctskymap::m_outradius
protected

Outer ring radius for RING background.

Definition at line 113 of file ctskymap.hpp.

Referenced by compute_maps_ring_fft(), copy_members(), get_parameters(), init_members(), ring_bounding_box(), setup_maps(), and write_hdu_keywords().

bool ctskymap::m_publish
protected

Publish sky map?

Definition at line 117 of file ctskymap.hpp.

Referenced by copy_members(), get_parameters(), init_members(), and process().

double ctskymap::m_roiradius
protected

Region of interest radius for RING bkg.

Definition at line 111 of file ctskymap.hpp.

Referenced by compute_maps_ring_fft(), copy_members(), get_parameters(), init_members(), setup_maps(), and write_hdu_keywords().

GSkyMap ctskymap::m_sigmap
protected
GSkyMap ctskymap::m_skymap
protected
std::vector<double> ctskymap::m_solidangle
protected

Cached pixel solid angles.

Definition at line 131 of file ctskymap.hpp.

Referenced by copy_members(), fill_maps_acceptance(), init_members(), and setup_maps().

double ctskymap::m_threshold
protected

Threshold for RING background.

Definition at line 115 of file ctskymap.hpp.

Referenced by compute_maps(), copy_members(), get_parameters(), init_members(), and write_hdu_keywords().

bool ctskymap::m_usefft
protected

Use FFT for RING background.

Definition at line 116 of file ctskymap.hpp.

Referenced by compute_maps(), copy_members(), get_parameters(), init_members(), and write_hdu_keywords().


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