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 | // MHAlphaEnergyTheta //
|
---|
29 | // //
|
---|
30 | // //
|
---|
31 | //////////////////////////////////////////////////////////////////////////////
|
---|
32 |
|
---|
33 | #include "MHAlphaEnergyTheta.h"
|
---|
34 |
|
---|
35 | #include <TCanvas.h>
|
---|
36 |
|
---|
37 | #include "MMcEvt.hxx"
|
---|
38 | #include "MHillasSrc.h"
|
---|
39 | #include "MEnergyEst.h"
|
---|
40 |
|
---|
41 | #include "MBinning.h"
|
---|
42 | #include "MParList.h"
|
---|
43 |
|
---|
44 | #include "MLog.h"
|
---|
45 | #include "MLogManip.h"
|
---|
46 |
|
---|
47 | ClassImp(MHAlphaEnergyTheta);
|
---|
48 |
|
---|
49 | // --------------------------------------------------------------------------
|
---|
50 | //
|
---|
51 | // Default Constructor. It sets name and title only. Typically you won't
|
---|
52 | // need to change this.
|
---|
53 | //
|
---|
54 | MHAlphaEnergyTheta::MHAlphaEnergyTheta(const char *name, const char *title)
|
---|
55 | {
|
---|
56 | //
|
---|
57 | // set the name and title of this object
|
---|
58 | //
|
---|
59 | fName = name ? name : "MHAlphaEnergyTheta";
|
---|
60 | fTitle = title ? title : "3-D histogram in alpha, energy and theta";
|
---|
61 |
|
---|
62 | fHist.SetDirectory(NULL);
|
---|
63 |
|
---|
64 | fHist.GetXaxis()->SetTitle("\\alpha [\\circ]");
|
---|
65 | fHist.GetYaxis()->SetTitle("E_{est} [GeV]");
|
---|
66 | fHist.GetZaxis()->SetTitle("\\Theta [\\circ]");
|
---|
67 | }
|
---|
68 |
|
---|
69 | Bool_t MHAlphaEnergyTheta::SetupFill(const MParList *plist)
|
---|
70 | {
|
---|
71 | fEnergy = (MEnergyEst*)plist->FindObject("MEnergyEst");
|
---|
72 | if (!fEnergy)
|
---|
73 | {
|
---|
74 | *fLog << err << dbginf << "MEnergyEst not found... aborting." << endl;
|
---|
75 | return kFALSE;
|
---|
76 | }
|
---|
77 |
|
---|
78 | fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
|
---|
79 | if (!fMcEvt)
|
---|
80 | {
|
---|
81 | *fLog << err << dbginf << "MMcEvt not found... aborting." << endl;
|
---|
82 | return kFALSE;
|
---|
83 | }
|
---|
84 |
|
---|
85 | MBinning* binsenergy = (MBinning*)plist->FindObject("BinningE");
|
---|
86 | MBinning* binsalpha = (MBinning*)plist->FindObject("BinningAlpha");
|
---|
87 | MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
|
---|
88 | if (!binsenergy || !binsalpha || !binstheta)
|
---|
89 | {
|
---|
90 | *fLog << err << dbginf << "At least one MBinning not found... aborting." << endl;
|
---|
91 | return kFALSE;
|
---|
92 | }
|
---|
93 |
|
---|
94 | SetBinning(&fHist, binsalpha, binsenergy, binstheta);
|
---|
95 |
|
---|
96 | return kTRUE;
|
---|
97 | }
|
---|
98 |
|
---|
99 | Bool_t MHAlphaEnergyTheta::Fill(const MParContainer *par)
|
---|
100 | {
|
---|
101 | MHillasSrc &hil = *(MHillasSrc*)par;
|
---|
102 |
|
---|
103 | fHist.Fill(hil.GetAlpha(), fEnergy->GetEnergy(), fMcEvt->GetTheta()*kRad2Deg);
|
---|
104 | return kTRUE;
|
---|
105 | }
|
---|
106 |
|
---|
107 | void MHAlphaEnergyTheta::Draw(Option_t *opt)
|
---|
108 | {
|
---|
109 | if (!gPad)
|
---|
110 | MakeDefCanvas("AlphaEnergyTheta", "Distrib of \\alpha, E, \\Theta");
|
---|
111 |
|
---|
112 | gPad->Divide(2,2);
|
---|
113 |
|
---|
114 | gPad->cd(1);
|
---|
115 | fHist.Project3D("x")->Draw(opt);
|
---|
116 |
|
---|
117 | gPad->cd(2);
|
---|
118 | fHist.Project3D("y")->Draw(opt);
|
---|
119 |
|
---|
120 | gPad->cd(3);
|
---|
121 | fHist.Project3D("z")->Draw(opt);
|
---|
122 |
|
---|
123 | gPad->cd(4);
|
---|
124 | fHist.Draw(opt);
|
---|
125 |
|
---|
126 | gPad->Modified();
|
---|
127 | gPad->Update();
|
---|
128 | }
|
---|
129 |
|
---|
130 | TObject *MHAlphaEnergyTheta::DrawClone(Option_t *opt) const
|
---|
131 | {
|
---|
132 | TCanvas *c = MakeDefCanvas("DiffTimeTheta", "Distrib of \\Delta t, \\Theta");
|
---|
133 | c->Divide(2, 2);
|
---|
134 |
|
---|
135 | gROOT->SetSelectedPad(NULL);
|
---|
136 |
|
---|
137 | //
|
---|
138 | // FIXME: ProjectionX,Y is not const within root
|
---|
139 | //
|
---|
140 | c->cd(1);
|
---|
141 | ((TH3*)(&fHist))->Project3D("x")->DrawCopy(opt);
|
---|
142 |
|
---|
143 | c->cd(2);
|
---|
144 | ((TH3*)(&fHist))->Project3D("y")->DrawCopy(opt);
|
---|
145 |
|
---|
146 | c->cd(3);
|
---|
147 | ((TH3*)(&fHist))->Project3D("z")->DrawCopy(opt);
|
---|
148 |
|
---|
149 | c->cd(4);
|
---|
150 | ((TH3*)(&fHist))->DrawCopy(opt);
|
---|
151 |
|
---|
152 | c->Modified();
|
---|
153 | c->Update();
|
---|
154 |
|
---|
155 | return c;
|
---|
156 | }
|
---|
157 |
|
---|
158 | void MHAlphaEnergyTheta::Substract(const TH3D *h1, const TH3D *h2)
|
---|
159 | {
|
---|
160 | MH::SetBinning(&fHist, h1);
|
---|
161 |
|
---|
162 | fHist.Sumw2();
|
---|
163 | fHist.Add((TH1*)h1, (TH1*)h2, 1, -1); // Root: FIXME
|
---|
164 | }
|
---|
165 |
|
---|
166 | void MHAlphaEnergyTheta::SetAlphaRange(Axis_t lo, Axis_t up)
|
---|
167 | {
|
---|
168 | TAxis &axe = *fHist.GetXaxis();
|
---|
169 |
|
---|
170 | //
|
---|
171 | // FIXME: ROOT Binning??? of projection
|
---|
172 | // root 3.02: SetRangeUser
|
---|
173 | axe.SetRange(axe.FindFixBin(lo), axe.FindFixBin(up));
|
---|
174 | }
|
---|
175 |
|
---|
176 | TH2D *MHAlphaEnergyTheta::GetAlphaProjection(Axis_t lo, Axis_t up)
|
---|
177 | {
|
---|
178 | SetAlphaRange(lo, up);
|
---|
179 | return (TH2D*)fHist.Project3D("yz");
|
---|
180 | }
|
---|