source: trunk/MagicSoft/Mars/mhist/MHTimeDiffTheta.cc@ 1302

Last change on this file since 1302 was 1295, checked in by wittek, 23 years ago
*** empty log message ***
File size: 5.3 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): Thomas Bretz 1/2002 <mailto:tbretz@uni-sw.gwdg.de>
19! Author(s): Wolfgang Wittek 1/2002 <mailto:wittek@mppmu.mpg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2002
22!
23!
24\* ======================================================================== */
25
26//////////////////////////////////////////////////////////////////////////////
27// //
28// MHTimeDiffTheta //
29// //
30// calculates the 2D-histogram time-difference vs. Theta //
31// //
32//////////////////////////////////////////////////////////////////////////////
33
34#include "MHTimeDiffTheta.h"
35
36#include <TCanvas.h>
37
38#include "MTime.h"
39#include "MMcEvt.hxx"
40
41#include "MBinning.h"
42#include "MParList.h"
43
44#include "MLog.h"
45#include "MLogManip.h"
46
47ClassImp(MHTimeDiffTheta);
48
49
50// --------------------------------------------------------------------------
51//
52// Default Constructor. It sets name and title of the histogram.
53//
54MHTimeDiffTheta::MHTimeDiffTheta(const char *name, const char *title)
55 : fLastTime(0), fHist()
56{
57 //
58 // set the name and title of this object
59 //
60 fName = name ? name : "MHTimeDiffTheta";
61 fTitle = title ? title : "2-D histogram in Theta and time difference";
62
63 fHist.SetDirectory(NULL);
64
65 fHist.SetXTitle("\\Delta t [s]");
66 fHist.SetYTitle("\\Theta [\\circ]");
67}
68
69// --------------------------------------------------------------------------
70//
71// Set the binnings and prepare the filling of the histogram
72//
73Bool_t MHTimeDiffTheta::SetupFill(const MParList *plist)
74{
75 fTime = (MTime*)plist->FindObject("MTime");
76 if (!fTime)
77 {
78 *fLog << err << dbginf << "MTime not found... aborting." << endl;
79 return kFALSE;
80 }
81
82 fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
83 if (!fMcEvt)
84 {
85 *fLog << err << dbginf << "MMcEvt not found... aborting." << endl;
86 return kFALSE;
87 }
88
89 const MBinning* binsdtime = (MBinning*)plist->FindObject("BinningTimeDiff");
90 const MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
91 if (!binstheta || !binsdtime)
92 {
93 *fLog << err << dbginf << "At least one MBinning not found... aborting." << endl;
94 return kFALSE;
95 }
96
97 SetBinning(&fHist, binsdtime, binstheta);
98
99 return kTRUE;
100}
101
102// --------------------------------------------------------------------------
103//
104// Draw a copy of the histogram
105//
106TObject *MHTimeDiffTheta::DrawClone(Option_t *opt) const
107{
108 TCanvas &c = *MakeDefCanvas("DiffTimeTheta", "Distrib of \\Delta t, Theta");
109
110 c.Divide(2, 2);
111
112 gROOT->SetSelectedPad(NULL);
113
114 //
115 // FIXME: ProjectionX,Y is not const within root
116 //
117
118 TH1D *h;
119
120 c.cd(1);
121 h = ((TH2*)&fHist)->ProjectionX("ProjX-Theta", -1, 9999, "E");
122
123 h->SetTitle("Distribution of \\Delta t [s]");
124 h->SetXTitle("\\Delta t [s]");
125 h->SetYTitle("Counts");
126
127 h->Draw(opt);
128 h->SetBit(kCanDelete);;
129 gPad->SetLogy();
130
131 c.cd(2);
132 h = ((TH2*)&fHist)->ProjectionY("ProjY-timediff", -1, 9999, "E");
133
134 h->SetTitle("Distribution of \\Theta [\\circ]");
135 h->SetXTitle("\\Theta [\\circ]");
136 h->SetYTitle("Counts");
137
138 h->Draw(opt);
139 h->SetBit(kCanDelete);;
140
141 c.cd(3);
142 ((TH2*)&fHist)->DrawCopy(opt);
143
144 c.Modified();
145 c.Update();
146
147 return &c;
148}
149
150// --------------------------------------------------------------------------
151//
152// Draw the histogram
153//
154void MHTimeDiffTheta::Draw(Option_t *opt)
155{
156 if (!gPad)
157 MakeDefCanvas("DiffTimeTheta", "Distrib of Delta t, Theta");
158
159 TH1D *h;
160
161 gPad->Divide(2,2);
162
163 gPad->cd(1);
164 h = fHist.ProjectionX("ProjX_Theta", -1, 9999, "E");
165
166 h->SetTitle("Distribution of \\Delta t [s]");
167 h->SetXTitle("\\Delta t [s]");
168 h->SetYTitle("Counts");
169
170 h->Draw(opt);
171 h->SetBit(kCanDelete);;
172 gPad->SetLogy();
173
174 gPad->cd(2);
175 h = fHist.ProjectionY("ProjY_timediff", -1, 9999, "E");
176
177 h->SetTitle("Distribution of \\Theta [\\circ]");
178 h->SetXTitle("\\Theta [\\circ]");
179 h->SetYTitle("Counts");
180
181 h->Draw(opt);
182 h->SetBit(kCanDelete);;
183
184 gPad->cd(3);
185 fHist.DrawCopy(opt);
186
187 gPad->Modified();
188 gPad->Update();
189
190}
191
192// --------------------------------------------------------------------------
193//
194// Fill the histogram
195//
196Bool_t MHTimeDiffTheta::Fill(const MParContainer *par)
197{
198 const Int_t time = fTime->GetTimeLo();
199
200 fHist.Fill(0.0001*(time-fLastTime), fMcEvt->GetTheta()*kRad2Deg);
201 fLastTime = time;
202
203 return kTRUE;
204}
205
Note: See TracBrowser for help on using the repository browser.