SatAndLight  2.2.2-hubble
Simulation toolkit for space telescopes
Telescope.h
Go to the documentation of this file.
1 // Author : florent robinet (IJCLab - Orsay): florent.robinet@ijclab.in2p3.fr
4 #ifndef __Telescope__
5 #define __Telescope__
6 
13 #include <TGraph.h>
14 #include <TH3.h>
15 #include "Camera.h"
16 #include "Source.h"
17 
18 using namespace std;
19 
123 class Telescope: public Camera {
124 
125  public:
126 
132 
146  Telescope(const int aNbits=14, const int aNy=256, const int aNz=256);
147 
149 
152  virtual ~Telescope(void);
158 
183  void TakePicture(const ULong64_t aTimeStart);
184 
186 
192  void StartObservation(TFile *aRootFile);
193 
195 
199  void StopObservation(void);
200 
202 
207  void SetName(const string aName);
208 
210 
213  inline string GetName(void){ return tel_name; };
214 
216 
227  void AddSource(Source *aSource, TGraph *aEffectiveArea, TH2F *aAcceptance);
228 
230 
238  TGraph* GetEffectiveArea(const UInt_t aSourceIndex){
239  return tel_effarea[aSourceIndex];
240  };
241 
243 
250  double GetEffectiveArea(const UInt_t aSourceIndex, const double aEnergy);
251 
253 
262  double GetEffectiveArea(const UInt_t aSourceIndex, const double aEnergyMin, const double aEnergyMax);
263 
265 
269  inline void SetFocalLength(const int aFocalLength=1000){ tel_foclen=aFocalLength; };
270 
272 
275  inline int GetFocalLength(void){ return tel_foclen; };
276 
278 
285  bool SetPsfFile(const string aFileName);
286 
288 
299  bool SetPsf(TH3F* aPsf);
300 
302 
315  bool SetPsf(TH2F* aPsf, const double aEnergyMin=0.0, const double aEnergyMax=1000000);
316 
318 
337  double SetPsfAiry(const double aOpeningAngle,
338  const int aYN, const int aZN,
339  const int aEnergyN, const double aEnergyMin, const double aEnergyMax);
340 
342 
346  inline TH3F* GetPsf(void){ return tel_psf; };
347 
349 
354  TH2F* GetPsf(const double aEnergy);
355 
357 
363  TH2F* GetPsf(const UInt_t aEnergyBinIndex);
364 
365  private:
366 
367  string tel_name;
370  TArrayD *tel_psf_ebins;
371  TH3F *tel_psf;
372  TH2F **tel_psfe;
373  vector <TGraph*> tel_effarea;
374  vector <TH2F*> tel_acc;
375  vector <Source*> tel_src;
376  TFile *tel_file;
377 
379 
392  void ApplyPsf(double& aYf, double& aZf, const double aY0, const double aZ0, const double aEnergy);
393 
394  ClassDef(Telescope,0)
395 };
396 
397 #endif
398 
399 
This module is used to describe a CCD camera.
This module is used to describe astrophysical sources.
CCD camera.
Definition: Camera.h:97
Describe astrophysical sources.
Definition: Source.h:62
Telescope.
Definition: Telescope.h:123
TFile * tel_file
current TFile.
Definition: Telescope.h:376
TArrayD * tel_psf_ebins
energy bins for the PSF.
Definition: Telescope.h:370
Telescope(const int aNbits=14, const int aNy=256, const int aNz=256)
Telescope class constuctor.
TH2F ** tel_psfe
telescope optical point spread function at one energy.
Definition: Telescope.h:372
string GetName(void)
Returns the telescope name.
Definition: Telescope.h:213
int GetFocalLength(void)
Returns the telescope focal length .
Definition: Telescope.h:275
vector< TGraph * > tel_effarea
telescope on-axis effective area to sources.
Definition: Telescope.h:373
vector< Source * > tel_src
vector of sources.
Definition: Telescope.h:375
int tel_psf_nebins
number of energy bins for the PSF.
Definition: Telescope.h:369
TGraph * GetEffectiveArea(const UInt_t aSourceIndex)
Returns the telescope on-axis effective area as a function of energy.
Definition: Telescope.h:238
vector< TH2F * > tel_acc
telescope acceptance to sources.
Definition: Telescope.h:374
int tel_foclen
telescope focal length .
Definition: Telescope.h:368
TH3F * tel_psf
telescope optical point spread function (PSF).
Definition: Telescope.h:371
TH3F * GetPsf(void)
Returns the point spread function.
Definition: Telescope.h:346
string tel_name
telescope name.
Definition: Telescope.h:367
void SetFocalLength(const int aFocalLength=1000)
Sets the telescope focal length.
Definition: Telescope.h:269