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

Last change on this file since 2482 was 2461, checked in by tbretz, 21 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
49using namespace std;
50
51// --------------------------------------------------------------------------
52//
53// Default Constructor. It sets name and title of the histogram.
54//
55MHTimeDiffTheta::MHTimeDiffTheta(const char *name, const char *title)
56 : fLastTime(0), fHist()
57{
58 //
59 // set the name and title of this object
60 //
61 fName = name ? name : "MHTimeDiffTheta";
62 fTitle = title ? title : "2-D histogram in Theta and time difference";
63
64 fHist.SetDirectory(NULL);
65
66 fHist.SetXTitle("\\Delta t [s]");
67 fHist.SetYTitle("\\Theta [\\circ]");
68}
69
70// --------------------------------------------------------------------------
71//
72// Set the binnings and prepare the filling of the histogram
73//
74Bool_t MHTimeDiffTheta::SetupFill(const MParList *plist)
75{
76 fTime = (MTime*)plist->FindObject("MTime");
77 if (!fTime)
78 {
79 *fLog << err << dbginf << "MTime not found... aborting." << endl;
80 return kFALSE;
81 }
82
83 fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
84 if (!fMcEvt)
85 {
86 *fLog << err << dbginf << "MMcEvt not found... aborting." << endl;
87 return kFALSE;
88 }
89
90 const MBinning* binsdtime = (MBinning*)plist->FindObject("BinningTimeDiff");
91 const MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
92 if (!binstheta || !binsdtime)
93 {
94 *fLog << err << dbginf << "At least one MBinning not found... aborting." << endl;
95 return kFALSE;
96 }
97
98 SetBinning(&fHist, binsdtime, binstheta);
99
100 return kTRUE;
101}
102
103// --------------------------------------------------------------------------
104//
105// Draw the histogram
106//
107void MHTimeDiffTheta::Draw(Option_t *opt)
108{
109 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
110 pad->SetBorderMode(0);
111
112 AppendPad("");
113
114 TH1 *h;
115
116 pad->Divide(2,2);
117
118 pad->cd(1);
119 gPad->SetLogy();
120 h = fHist.ProjectionX("ProjX_Theta", -1, 9999, "E");
121 h->SetTitle("Distribution of \\Delta t [s]");
122 h->SetXTitle("\\Delta t [s]");
123 h->SetYTitle("Counts");
124 h->Draw(opt);
125 h->SetBit(kCanDelete);;
126
127 pad->cd(2);
128 h = fHist.ProjectionY("ProjY_timediff", -1, 9999, "E");
129 h->SetTitle("Distribution of \\Theta [\\circ]");
130 h->SetXTitle("\\Theta [\\circ]");
131 h->SetYTitle("Counts");
132 h->Draw(opt);
133 h->SetBit(kCanDelete);;
134
135 pad->cd(3);
136 fHist.Draw(opt);
137
138 pad->Modified();
139 pad->Update();
140}
141
142// --------------------------------------------------------------------------
143//
144// Fill the histogram
145//
146Bool_t MHTimeDiffTheta::Fill(const MParContainer *par, const Stat_t w)
147{
148 const Double_t time = *fTime;
149
150 fHist.Fill(time-fLastTime, fMcEvt->GetTelescopeTheta()*kRad2Deg, w);
151 fLastTime = time;
152
153 return kTRUE;
154}
155
Note: See TracBrowser for help on using the repository browser.