source: trunk/MagicSoft/Mars/mhist/MHMcEnergy.cc@ 861

Last change on this file since 861 was 861, checked in by jlopez, 23 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 5.2 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/////////////////////////////////////////////////////////////////////////////
26//
27// MHMcEnergy
28//
29// This class holds the information ( histogram and fit function )
30// about the energy threshold for a particular trigger condition.
31//
32////////////////////////////////////////////////////////////////////////////
33#include "MHMcEnergy.h"
34
35#include <iostream.h>
36
37#include <TH1.h>
38#include <TF1.h>
39#include <TCanvas.h>
40#include <TPaveLabel.h>
41
42ClassImp(MHMcEnergy);
43
44// -------------------------------------------------------------------------
45//
46// Default Constructor.
47//
48MHMcEnergy::MHMcEnergy(const UInt_t idx, const char *name, const char *title)
49{
50 char aux[15]="MHMcEnergy";
51
52 if (idx>0)
53 sprintf(aux+10, ";%i", idx);
54 *fName = name ? name : aux;
55 *fTitle = title ? title : "Container for a MC energy histogram" ;
56
57 // - we initialize the histogram and the gaus function
58 // - we have to give diferent names for the diferent histograms because
59 // root don't allow us to have diferent histograms with the same name
60
61 strcpy(aux, "fLogEnergy");
62 if (idx>0)
63 sprintf(aux+10, ";%i", idx);
64 fLogEner = new TF1(aux, "gaus", 1., 3.);
65
66 strcpy(aux, "hLogEnergy");
67 if (idx>0)
68 sprintf(aux+10, ";%i", idx);
69 hLogEner = new TH1F(aux, "", 40, 0.5, 4.5);
70 hLogEner->SetXTitle("log(E) [GeV]");
71 hLogEner->SetYTitle("dN/dE");
72}
73
74//-------------------------------------------------------------------------
75//
76// Defualt Destructor
77//
78MHMcEnergy::~MHMcEnergy()
79{
80 delete hLogEner;
81 delete fLogEner;
82}
83
84//--------------------------------------------------------------------------
85//
86// Fill the histogram with the log10 of the energy for triggered events.
87//
88void MHMcEnergy::Fill(Float_t log10E, Float_t w)
89{
90 hLogEner->Fill(log10E, w);
91}
92
93// -------------------------------------------------------------------------
94//
95// Fitting function
96//
97void MHMcEnergy::Fit(Axis_t xxmin, Axis_t xxmax)
98{
99 //
100 // 0: don't draw the function (it is drawn together with the histogram)
101 // Q: quiet mode
102 //
103 hLogEner->Fit(fLogEner->GetName(), "Q0", "", xxmin, xxmax);
104}
105
106// ------------------------------------------------------------------------
107//
108// Drawing function. It creates its own canvas.
109//
110void MHMcEnergy::Draw(Option_t *option)
111{
112 char text[50];
113 sprintf(text, "Energy Threshold = %4.1f +- %4.1f GeV",
114 GetThreshold(), GetThresholdErr());
115
116 const Float_t min = hLogEner->GetMinimum();
117 const Float_t max = hLogEner->GetMaximum();
118 const Float_t sum = min+max;
119
120 TCanvas *c=new TCanvas("Energy Distribution","Energy distribution for triggered events");
121
122 hLogEner->Draw(option);
123
124 TPaveLabel* label = new TPaveLabel(2.2, 0.75*sum, 4.4, 0.90*sum, text);
125 label->SetFillColor(10);
126 label->SetTextSize(0.3);
127 label->Draw();
128
129 c->Modified();
130 c->Update();
131}
132
133// --------------------------------------------------------------------------
134//
135// Set the number of bins in the histogran.
136//
137void MHMcEnergy::SetBins(Int_t nbins)
138{
139 hLogEner->SetBins(nbins,0.5,4.5);
140}
141// --------------------------------------------------------------------------
142//
143// Write the threshold and its error in the standard output
144//
145void MHMcEnergy::Print(Option_t*)
146{
147 cout << "Threshold: " << GetThreshold() << " +- " << GetThresholdErr() << endl;
148}
149
150// -------------------------------------------------------------------------
151//
152// Return the threshold
153//
154Float_t MHMcEnergy::GetThreshold() const
155{
156 const Float_t p1 = fLogEner->GetParameter(1);
157
158 return pow(10, p1);
159}
160
161// -------------------------------------------------------------------------
162//
163// Return the error of the threshold.
164//
165Float_t MHMcEnergy::GetThresholdErr() const
166{
167 const Float_t lg10 = log(10);
168 const Float_t p1 = fLogEner->GetParameter(1);
169 const Float_t p1err = fLogEner->GetParError(1);
170
171 // The error has into accuont the error in the fit
172 return pow(10, p1) * p1err * lg10;
173}
174
175// -------------------------------------------------------------------------
176//
177// Return the peak of the fitted gaussan function.
178//
179Float_t MHMcEnergy::GetGaussPeak() const
180{
181 return fLogEner->GetParameter(1);
182}
183
184// -------------------------------------------------------------------------
185//
186// Return the sigma of the fitted gaussan function.
187//
188Float_t MHMcEnergy::GetGaussSigma() const
189{
190 return fLogEner->GetParameter(2);
191}
192
193
194
195
Note: See TracBrowser for help on using the repository browser.