source: trunk/MagicSoft/Mars/mmontecarlo/MMcEnerHisto.cc@ 853

Last change on this file since 853 was 834, checked in by jlopez, 23 years ago
Tobject to store the energy threshold information
File size: 3.0 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Javier Lopez 05/2001 (jlopez@ifae.es)
19!
20! Copyright: MAGIC Software Development, 2000-2001
21!
22!
23\* ======================================================================== */
24
25#include "MMcEnerHisto.h"
26
27#include <MLog.h>
28#include <TH1.h>
29#include <TF1.h>
30#include <TPaveLabel.h>
31
32ClassImp(MMcEnerHisto)
33
34MMcEnerHisto::MMcEnerHisto(const int index)
35{
36 //
37 // default constructor
38 //
39
40 // - we initialize the histogram and the gaus function
41 // - we have to give diferent names for the diferent histograms because
42 // root don't allow us to have diferent histograms with the same name
43
44 char aux[15];
45 sprintf (aux,"hLogEner%i",index);
46 hLogEner = new TH1F(aux,"",100,0.5,4.5);
47 sprintf (aux,"fLogEner%i",index);
48 fLogEner = new TF1(aux,"gaus",1.,3.);
49
50}
51
52MMcEnerHisto::~MMcEnerHisto()
53{
54 // we can not delete the histogram because it dissappear
55 // when we display it using root
56
57 // delete hLogEner ;
58 delete fLogEner ;
59}
60
61void MMcEnerHisto::SetBins(Int_t nbins, Float_t xmin, Float_t xmax)
62{
63 hLogEner->SetBins(nbins,xmin,xmax);
64}
65
66void MMcEnerHisto::Fill(Float_t log10E, Float_t w)
67{
68 hLogEner->Fill(log10E,w) ;
69}
70
71void MMcEnerHisto::Fit(const char *fname, Option_t *option, Option_t *goption, Axis_t xxmin, Axis_t xxmax)
72{
73 hLogEner->Fit(fname,option,goption,xxmin,xxmax);
74}
75
76void MMcEnerHisto::Draw(Option_t* option)
77{
78 char text[50];
79 sprintf(text,"Energy Threshold = %4.1f +- %4.1f GeV",GetEnerThre(),GetEnerThreErr());
80
81 TPaveLabel* label = new TPaveLabel(2.2,0.75*(hLogEner->GetMinimum()+hLogEner->GetMaximum()),4.4,0.90*(hLogEner->GetMinimum()+hLogEner->GetMaximum()),text);
82
83 hLogEner->SetYTitle("dN/dE") ;
84 hLogEner->SetXTitle("Log(E) [GeV]") ;
85 hLogEner->Draw(option) ;
86 label->SetFillColor(10);
87 label->SetTextSize(0.3);
88 label->Draw();
89}
90
91void MMcEnerHisto::Print()
92{
93 cout << "Energy Threshold: " << GetEnerThre() << " +- " << GetEnerThreErr() << endl;
94}
95
96Float_t MMcEnerHisto::GetEnerThre()
97{
98 return pow(10,fLogEner->GetParameter(1));
99}
100
101Float_t MMcEnerHisto::GetEnerThreErr()
102{
103 return pow(10,fLogEner->GetParameter(1))*log(10)*fLogEner->GetParError(1);
104}
105
106Float_t MMcEnerHisto::GetPeakAtLogE()
107{
108 return fLogEner->GetParameter(1);
109}
110
111Float_t MMcEnerHisto::GetSigmaAtLogE()
112{
113 return fLogEner->GetParameter(2);
114}
115
116
117
Note: See TracBrowser for help on using the repository browser.