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

Last change on this file since 1975 was 1966, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 4.4 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 the histogram
105//
106void MHTimeDiffTheta::Draw(Option_t *opt)
107{
108 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
109 pad->SetBorderMode(0);
110
111 AppendPad("");
112
113 TH1 *h;
114
115 pad->Divide(2,2);
116
117 pad->cd(1);
118 gPad->SetLogy();
119 h = fHist.ProjectionX("ProjX_Theta", -1, 9999, "E");
120 h->SetTitle("Distribution of \\Delta t [s]");
121 h->SetXTitle("\\Delta t [s]");
122 h->SetYTitle("Counts");
123 h->Draw(opt);
124 h->SetBit(kCanDelete);;
125
126 pad->cd(2);
127 h = fHist.ProjectionY("ProjY_timediff", -1, 9999, "E");
128 h->SetTitle("Distribution of \\Theta [\\circ]");
129 h->SetXTitle("\\Theta [\\circ]");
130 h->SetYTitle("Counts");
131 h->Draw(opt);
132 h->SetBit(kCanDelete);;
133
134 pad->cd(3);
135 fHist.Draw(opt);
136
137 pad->Modified();
138 pad->Update();
139}
140
141// --------------------------------------------------------------------------
142//
143// Fill the histogram
144//
145Bool_t MHTimeDiffTheta::Fill(const MParContainer *par)
146{
147 const Double_t time = 200e-9*fTime->GetTimeLo() + fTime->GetTimeHi();
148
149 fHist.Fill(time-fLastTime, fMcEvt->GetTelescopeTheta()*kRad2Deg);
150 fLastTime = time;
151
152 return kTRUE;
153}
154
Note: See TracBrowser for help on using the repository browser.