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 time saving 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): Wolfgang Wittek 6/2002 <mailto:wittek@mppmu.mpg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2002
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | //////////////////////////////////////////////////////////////////////////////
|
---|
26 | // //
|
---|
27 | // MHGamma //
|
---|
28 | // //
|
---|
29 | // manipulation of alpha distributions //
|
---|
30 | // //
|
---|
31 | //////////////////////////////////////////////////////////////////////////////
|
---|
32 |
|
---|
33 | #include "MHGamma.h"
|
---|
34 |
|
---|
35 | #include <TCanvas.h>
|
---|
36 | #include <TPad.h>
|
---|
37 |
|
---|
38 | #include <math.h>
|
---|
39 |
|
---|
40 | #include <TH1.h>
|
---|
41 | #include <TH2.h>
|
---|
42 | #include <TH3.h>
|
---|
43 |
|
---|
44 | #include "MLog.h"
|
---|
45 | #include "MLogManip.h"
|
---|
46 |
|
---|
47 | ClassImp(MHGamma);
|
---|
48 |
|
---|
49 |
|
---|
50 | // --------------------------------------------------------------------------
|
---|
51 | //
|
---|
52 | // Default Constructor.
|
---|
53 | //
|
---|
54 | MHGamma::MHGamma()
|
---|
55 | {
|
---|
56 | }
|
---|
57 |
|
---|
58 | // --------------------------------------------------------------------------
|
---|
59 | //
|
---|
60 | // Dummy Fill
|
---|
61 | //
|
---|
62 | Bool_t MHGamma::Fill(const MParContainer *par)
|
---|
63 | {
|
---|
64 | return kTRUE;
|
---|
65 | }
|
---|
66 |
|
---|
67 | // --------------------------------------------------------------------------
|
---|
68 | //
|
---|
69 | // Calculate the histogram as the difference of two histograms :
|
---|
70 | // fHist(gamma) = h1(source) - h2(antisource)
|
---|
71 | //
|
---|
72 | TH3D *MHGamma::Subtract(const TH3D *h1, const TH3D *h2,
|
---|
73 | const char *name, const char *title,
|
---|
74 | Bool_t Draw)
|
---|
75 | {
|
---|
76 | TH3D *fHist;
|
---|
77 | fHist = new TH3D();
|
---|
78 | fHist->SetName(name);
|
---|
79 | fHist->SetTitle(title);
|
---|
80 |
|
---|
81 | fHist->SetDirectory(NULL);
|
---|
82 |
|
---|
83 | // SetBinning((TH3D*)fHist, (TH3D*)h1);
|
---|
84 | SetBinning((TH1*)fHist, (TH1*)h1);
|
---|
85 |
|
---|
86 | TString strg1 = (((TH1*)h1)->GetXaxis())->GetTitle();
|
---|
87 | TString strg2 = (((TH1*)h1)->GetYaxis())->GetTitle();
|
---|
88 | TString strg3 = (((TH1*)h1)->GetZaxis())->GetTitle();
|
---|
89 | fHist->SetXTitle(strg1);
|
---|
90 | fHist->SetYTitle(strg2);
|
---|
91 | fHist->SetZTitle(strg3);
|
---|
92 |
|
---|
93 | fHist->Add((TH1*)h1, (TH1*)h2, 1, -1); // ROOT: FIXME!
|
---|
94 |
|
---|
95 | //...........................................................
|
---|
96 | // draw histogram
|
---|
97 | if (Draw == kTRUE)
|
---|
98 | {
|
---|
99 | TString strg7 = "3D-plot of ";
|
---|
100 | strg7 += strg1;
|
---|
101 | strg7 += ",";
|
---|
102 | strg7 += strg2;
|
---|
103 | strg7 += ",";
|
---|
104 | strg7 += strg3;
|
---|
105 | strg7 += " for \'gamma\' sample";
|
---|
106 |
|
---|
107 | TCanvas &c = *MakeDefCanvas("Alpha", strg7);
|
---|
108 |
|
---|
109 | c.Divide(2, 2);
|
---|
110 |
|
---|
111 | gROOT->SetSelectedPad(NULL);
|
---|
112 |
|
---|
113 | TH1 *h;
|
---|
114 |
|
---|
115 | c.cd(1);
|
---|
116 | h = ((TH3D*)(fHist))->Project3D("ex");
|
---|
117 |
|
---|
118 | TString strg0 = "SRC-ASRC : ";
|
---|
119 | TString strg4 = strg0 + strg1;
|
---|
120 | h->SetTitle(strg4);
|
---|
121 | h->SetXTitle(strg1);
|
---|
122 | h->SetYTitle("Counts");
|
---|
123 |
|
---|
124 | h->Draw();
|
---|
125 | h->SetBit(kCanDelete);
|
---|
126 |
|
---|
127 | c.cd(2);
|
---|
128 | h = ((TH3D*)(fHist))->Project3D("ey");
|
---|
129 |
|
---|
130 | TString strg5 = strg0 + strg2;
|
---|
131 | h->SetTitle(strg5);
|
---|
132 | h->SetXTitle(strg2);
|
---|
133 | h->SetYTitle("Counts");
|
---|
134 |
|
---|
135 | h->Draw();
|
---|
136 | h->SetBit(kCanDelete);
|
---|
137 | gPad->SetLogx();
|
---|
138 |
|
---|
139 | c.cd(3);
|
---|
140 | h = ((TH3D*)(fHist))->Project3D("ez");
|
---|
141 |
|
---|
142 | TString strg6 = strg0 + strg3;
|
---|
143 | h->SetTitle(strg6);
|
---|
144 | h->SetXTitle(strg3);
|
---|
145 | h->SetYTitle("Counts");
|
---|
146 |
|
---|
147 | h->Draw();
|
---|
148 | h->SetBit(kCanDelete);
|
---|
149 |
|
---|
150 | c.cd(4);
|
---|
151 | ((TH3D*)fHist)->DrawCopy();
|
---|
152 |
|
---|
153 | c.Modified();
|
---|
154 | c.Update();
|
---|
155 | }
|
---|
156 |
|
---|
157 | return fHist;
|
---|
158 | }
|
---|
159 |
|
---|
160 | // --------------------------------------------------------------------------
|
---|
161 | //
|
---|
162 | // Integrate fHist(gamma) in the alpha range (lo, up)
|
---|
163 | //
|
---|
164 | TH2D *MHGamma::GetAlphaProjection(TH3D *fHist, Axis_t lo, Axis_t up,
|
---|
165 | Bool_t Drawp)
|
---|
166 | {
|
---|
167 | if (up < lo)
|
---|
168 | {
|
---|
169 | *fLog << err << fName << "MHGamma : Alpha projection not possible: lo=" << lo << " up=" << up << endl;
|
---|
170 | return NULL;
|
---|
171 | }
|
---|
172 |
|
---|
173 | TAxis &axe = *fHist->GetXaxis();
|
---|
174 |
|
---|
175 | Int_t ilo = axe.FindFixBin(lo);
|
---|
176 | Int_t iup = axe.FindFixBin(up);
|
---|
177 |
|
---|
178 | const Double_t epslo1 = lo-axe.GetBinLowEdge(ilo);
|
---|
179 | const Double_t epslo2 = axe.GetBinUpEdge(ilo)-lo;
|
---|
180 |
|
---|
181 | const Double_t epsup1 = up-axe.GetBinLowEdge(iup);
|
---|
182 | const Double_t epsup2 = axe.GetBinUpEdge(iup)-up;
|
---|
183 |
|
---|
184 | const Double_t epslo = epslo1<epslo2 ? epslo1 : epslo2;
|
---|
185 | const Double_t epsup = epsup1<epsup2 ? epsup1 : epsup2;
|
---|
186 |
|
---|
187 | if (epslo1>epslo2)
|
---|
188 | ilo++;
|
---|
189 |
|
---|
190 | if (epsup1<epsup2)
|
---|
191 | iup--;
|
---|
192 |
|
---|
193 | if (epslo>0.01*axe.GetBinWidth(ilo) || epsup>0.01*axe.GetBinWidth(iup))
|
---|
194 | {
|
---|
195 | *fLog << err << fName << "MHGamma : binning is not adequate for the requested projection:" << endl;
|
---|
196 | *fLog << "Please specify a lower or upper limit which is not more than 1% away from a bin edge" << endl;
|
---|
197 | *fLog << " epslo = " << epslo << endl;
|
---|
198 | *fLog << " epsup = " << epsup << endl;
|
---|
199 | *fLog << " dwl = " << axe.GetBinWidth(ilo) << endl;
|
---|
200 | *fLog << " dwu = " << axe.GetBinWidth(iup) << endl;
|
---|
201 | return NULL;
|
---|
202 | }
|
---|
203 |
|
---|
204 | axe.SetRange(ilo, iup);
|
---|
205 |
|
---|
206 | TH2D &h2D = *(TH2D *)fHist->Project3D("ezy");
|
---|
207 |
|
---|
208 | TString strg0 = "2D-plot of ";
|
---|
209 | TString strg1 = (fHist->GetYaxis())->GetTitle();
|
---|
210 | TString strg2 = (fHist->GetZaxis())->GetTitle();
|
---|
211 | strg0 += strg2;
|
---|
212 | strg0 += " vs. ";
|
---|
213 | strg0 += strg1;
|
---|
214 | h2D.SetTitle(strg0);
|
---|
215 | h2D.SetXTitle(strg1);
|
---|
216 | h2D.SetYTitle(strg2);
|
---|
217 |
|
---|
218 |
|
---|
219 | //...........................................................
|
---|
220 | // draw histogram
|
---|
221 | if (Drawp == kTRUE)
|
---|
222 | {
|
---|
223 | char txt[100];
|
---|
224 | TString strg3 = "No.of Gammas vs. ";
|
---|
225 | strg3 += strg1;
|
---|
226 | strg3 += " and ";
|
---|
227 | strg3 += strg2;
|
---|
228 | sprintf(txt, " (%.1f < alpha < %.1f deg)", lo, up);
|
---|
229 | strg3 += txt;
|
---|
230 |
|
---|
231 | TCanvas &c = *MakeDefCanvas("Gamma", strg3);
|
---|
232 |
|
---|
233 | c.Divide(2, 2);
|
---|
234 |
|
---|
235 | gROOT->SetSelectedPad(NULL);
|
---|
236 |
|
---|
237 | TH1 *h;
|
---|
238 |
|
---|
239 | c.cd(1);
|
---|
240 | h = h2D.ProjectionX("xpro", -1, 9999, "E");
|
---|
241 | TString strg0 = "No.of gammas : ";
|
---|
242 | TString strg7 = strg0 + strg1;
|
---|
243 | h->SetTitle(strg7);
|
---|
244 | h->SetXTitle(strg1);
|
---|
245 | h->SetYTitle("Counts");
|
---|
246 |
|
---|
247 | h->Draw();
|
---|
248 | h->SetBit(kCanDelete);
|
---|
249 | gPad->SetLogx();
|
---|
250 |
|
---|
251 | c.cd(2);
|
---|
252 | h = h2D.ProjectionY("ypro", -1, 9999, "E");
|
---|
253 | TString strg8 = strg0 + strg2;
|
---|
254 | h->SetTitle(strg8);
|
---|
255 | h->SetXTitle(strg2);
|
---|
256 | h->SetYTitle("Counts");
|
---|
257 |
|
---|
258 | h->Draw();
|
---|
259 | h->SetBit(kCanDelete);
|
---|
260 |
|
---|
261 | c.cd(3);
|
---|
262 |
|
---|
263 | h2D.DrawCopy();
|
---|
264 | gPad->SetLogx();
|
---|
265 |
|
---|
266 | c.Modified();
|
---|
267 | c.Update();
|
---|
268 | }
|
---|
269 | //...........................................................
|
---|
270 |
|
---|
271 | return &h2D;
|
---|
272 | }
|
---|