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 2002 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2002
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | /////////////////////////////////////////////////////////////////////////////
|
---|
26 | //
|
---|
27 | // MH3
|
---|
28 | //
|
---|
29 | // With this histogram you can fill a histogram with up to three
|
---|
30 | // variables from Mars parameter containers in an eventloop.
|
---|
31 | //
|
---|
32 | // In the constructor you can give up to three variables which should be
|
---|
33 | // filled in the histogram. Dependend on the number of given variables
|
---|
34 | // (data members) a TH1F, TH2F or TH3F is created.
|
---|
35 | // Specify the data mamber like the following:
|
---|
36 | // "MHillas.fLength"
|
---|
37 | // Where MHillas is the name of the parameter container in the parameter
|
---|
38 | // list and fLength is the name of the data member which should be filled
|
---|
39 | // in the histogram. Assuming that your MHillas container has a different
|
---|
40 | // name (MyHillas) the name to give would be:
|
---|
41 | // "MyHillas.fLength"
|
---|
42 | //
|
---|
43 | // The axis binning is retrieved from the parameter list, too. Create a
|
---|
44 | // MBinning with the name "Binning" plus the name of your MH3 container
|
---|
45 | // plus the axis name ("X", "Y" or "Z") and add it to the parameter list.
|
---|
46 | //
|
---|
47 | // If you want to use a different unit for histogramming use SetScaleX,
|
---|
48 | // SetScaleY and SetScaleZ.
|
---|
49 | //
|
---|
50 | // For example:
|
---|
51 | // MH3 myhist("MHillas.fLength");
|
---|
52 | // myhist.SetName("MyHist");
|
---|
53 | // myhist.SetScaleX(geomcam.GetConvMm2Deg()); //convert length to degree
|
---|
54 | // MBinning bins("BinningMyHistX");
|
---|
55 | // bins.SetEdges(10, 0, 150);
|
---|
56 | // plist.AddToList(&myhist);
|
---|
57 | // plist.AddToList(&bins);
|
---|
58 | //
|
---|
59 | /////////////////////////////////////////////////////////////////////////////
|
---|
60 |
|
---|
61 | #include "MH3.h"
|
---|
62 |
|
---|
63 | #include <TH2.h>
|
---|
64 | #include <TH3.h>
|
---|
65 | #include <TProfile.h>
|
---|
66 | #include <TProfile2D.h>
|
---|
67 | #include <TPad.h>
|
---|
68 | #include <TCanvas.h>
|
---|
69 | #include <TDataMember.h>
|
---|
70 | #include <TMethodCall.h>
|
---|
71 |
|
---|
72 | #include "MLog.h"
|
---|
73 | #include "MLogManip.h"
|
---|
74 |
|
---|
75 | #include "MParList.h"
|
---|
76 |
|
---|
77 | #include "MDataChain.h"
|
---|
78 |
|
---|
79 | ClassImp(MH3);
|
---|
80 |
|
---|
81 | // --------------------------------------------------------------------------
|
---|
82 | //
|
---|
83 | // Creates an TH1F. memberx is filled into the X-bins. For a more detailed
|
---|
84 | // description see the class description above.
|
---|
85 | //
|
---|
86 | MH3::MH3(const char *memberx)
|
---|
87 | : fDimension(1)
|
---|
88 | {
|
---|
89 | fHist = new TH1F;
|
---|
90 |
|
---|
91 | fData[0] = new MDataChain(memberx);
|
---|
92 | fData[1] = NULL;
|
---|
93 | fData[2] = NULL;
|
---|
94 |
|
---|
95 | fName = "MH3";
|
---|
96 | fTitle = "Container for a 1D Mars Histogram";
|
---|
97 |
|
---|
98 | fHist->SetDirectory(NULL);
|
---|
99 | fHist->SetYTitle("Counts");
|
---|
100 |
|
---|
101 | fScale[0] = 1;
|
---|
102 | fScale[1] = 1;
|
---|
103 | fScale[2] = 1;
|
---|
104 | }
|
---|
105 |
|
---|
106 | // --------------------------------------------------------------------------
|
---|
107 | //
|
---|
108 | // Creates an TH2F. memberx is filled into the X-bins. membery is filled
|
---|
109 | // into the Y-bins. For a more detailed description see the class
|
---|
110 | // description above.
|
---|
111 | //
|
---|
112 | MH3::MH3(const char *memberx, const char *membery)
|
---|
113 | : fDimension(2)
|
---|
114 | {
|
---|
115 | fHist = new TH2F;
|
---|
116 |
|
---|
117 | fData[0] = new MDataChain(memberx);
|
---|
118 | fData[1] = new MDataChain(membery);
|
---|
119 | fData[2] = NULL;
|
---|
120 |
|
---|
121 | fName = "MH3";
|
---|
122 | fTitle = "Container for a 2D Mars Histogram";
|
---|
123 |
|
---|
124 | fHist->SetDirectory(NULL);
|
---|
125 | fHist->SetZTitle("Counts");
|
---|
126 |
|
---|
127 | fScale[0] = 1;
|
---|
128 | fScale[1] = 1;
|
---|
129 | fScale[2] = 1;
|
---|
130 | }
|
---|
131 |
|
---|
132 | // --------------------------------------------------------------------------
|
---|
133 | //
|
---|
134 | // Creates an TH3F. memberx is filled into the X-bins. membery is filled
|
---|
135 | // into the Y-bins. membery is filled into the Z-bins. For a more detailed
|
---|
136 | // description see the class description above.
|
---|
137 | //
|
---|
138 | MH3::MH3(const char *memberx, const char *membery, const char *memberz)
|
---|
139 | : fDimension(3)
|
---|
140 | {
|
---|
141 | fHist = new TH3F;
|
---|
142 |
|
---|
143 | fData[0] = new MDataChain(memberx);
|
---|
144 | fData[1] = new MDataChain(membery);
|
---|
145 | fData[2] = new MDataChain(memberz);
|
---|
146 |
|
---|
147 | fName = "MH3";
|
---|
148 | fTitle = "Container for a 3D Mars Histogram";
|
---|
149 |
|
---|
150 | fHist->SetDirectory(NULL);
|
---|
151 |
|
---|
152 | fScale[0] = 1;
|
---|
153 | fScale[1] = 1;
|
---|
154 | fScale[2] = 1;
|
---|
155 | }
|
---|
156 |
|
---|
157 | // --------------------------------------------------------------------------
|
---|
158 | //
|
---|
159 | // Deletes the histogram
|
---|
160 | //
|
---|
161 | MH3::~MH3()
|
---|
162 | {
|
---|
163 | delete fHist;
|
---|
164 |
|
---|
165 | for (int i=0; i<3; i++)
|
---|
166 | if (fData[i])
|
---|
167 | delete fData[i];
|
---|
168 | }
|
---|
169 |
|
---|
170 | // --------------------------------------------------------------------------
|
---|
171 | //
|
---|
172 | // Setup the Binning for the histograms automatically if the correct
|
---|
173 | // instances of MBinning are found in the parameter list
|
---|
174 | // For a more detailed description see class description above.
|
---|
175 | //
|
---|
176 | Bool_t MH3::SetupFill(const MParList *plist)
|
---|
177 | {
|
---|
178 | TString bname("Binning");
|
---|
179 | bname += fName;
|
---|
180 |
|
---|
181 | MBinning *binsx = NULL;
|
---|
182 | MBinning *binsy = NULL;
|
---|
183 | MBinning *binsz = NULL;
|
---|
184 | switch (fDimension)
|
---|
185 | {
|
---|
186 | case 3:
|
---|
187 | binsz = (MBinning*)plist->FindObject(bname+"Z");
|
---|
188 | if (!binsz)
|
---|
189 | {
|
---|
190 | *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl;
|
---|
191 | return kFALSE;
|
---|
192 | }
|
---|
193 | fHist->SetZTitle(fData[2]->GetTitle());
|
---|
194 | if (!fData[2]->PreProcess(plist))
|
---|
195 | return kFALSE;
|
---|
196 | case 2:
|
---|
197 | binsy = (MBinning*)plist->FindObject(bname+"Y");
|
---|
198 | if (!binsy)
|
---|
199 | {
|
---|
200 | *fLog << err << dbginf << "MBinning '" << bname << "Y' not found... aborting." << endl;
|
---|
201 | return kFALSE;
|
---|
202 | }
|
---|
203 | fHist->SetYTitle(fData[1]->GetTitle());
|
---|
204 | if (!fData[1]->PreProcess(plist))
|
---|
205 | return kFALSE;
|
---|
206 | case 1:
|
---|
207 | binsx = (MBinning*)plist->FindObject(bname+"X");
|
---|
208 | if (!binsx)
|
---|
209 | {
|
---|
210 | *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl;
|
---|
211 | return kFALSE;
|
---|
212 | }
|
---|
213 | fHist->SetXTitle(fData[0]->GetTitle());
|
---|
214 | if (!fData[0]->PreProcess(plist))
|
---|
215 | return kFALSE;
|
---|
216 | }
|
---|
217 |
|
---|
218 | fHist->SetName(fName);
|
---|
219 |
|
---|
220 | TString title("Histogram for ");
|
---|
221 | title += fName;
|
---|
222 |
|
---|
223 | switch (fDimension)
|
---|
224 | {
|
---|
225 | case 1:
|
---|
226 | fHist->SetTitle(title+" (1D)");
|
---|
227 | SetBinning(fHist, binsx);
|
---|
228 | return kTRUE;
|
---|
229 | case 2:
|
---|
230 | fHist->SetTitle(title+" (2D)");
|
---|
231 | SetBinning((TH2*)fHist, binsx, binsy);
|
---|
232 | return kTRUE;
|
---|
233 | case 3:
|
---|
234 | fHist->SetTitle(title+" (3D)");
|
---|
235 | SetBinning((TH3*)fHist, binsx, binsy, binsz);
|
---|
236 | return kTRUE;
|
---|
237 | }
|
---|
238 |
|
---|
239 | return kTRUE;
|
---|
240 | }
|
---|
241 |
|
---|
242 | // --------------------------------------------------------------------------
|
---|
243 | //
|
---|
244 | // Set the name of the histogram ant the MH3 container
|
---|
245 | //
|
---|
246 | void MH3::SetName(const char *name)
|
---|
247 | {
|
---|
248 | fHist->SetName(name);
|
---|
249 | MParContainer::SetName(name);
|
---|
250 | }
|
---|
251 |
|
---|
252 | // --------------------------------------------------------------------------
|
---|
253 | //
|
---|
254 | // Set the title of the histogram ant the MH3 container
|
---|
255 | //
|
---|
256 | void MH3::SetTitle(const char *title)
|
---|
257 | {
|
---|
258 | fHist->SetTitle(title);
|
---|
259 | MParContainer::SetTitle(title);
|
---|
260 | }
|
---|
261 |
|
---|
262 | // --------------------------------------------------------------------------
|
---|
263 | //
|
---|
264 | // Fills the one, two or three data members into our histogram
|
---|
265 | //
|
---|
266 | Bool_t MH3::Fill(const MParContainer *par)
|
---|
267 | {
|
---|
268 | Double_t x=0;
|
---|
269 | Double_t y=0;
|
---|
270 | Double_t z=0;
|
---|
271 |
|
---|
272 | switch (fDimension)
|
---|
273 | {
|
---|
274 | case 3:
|
---|
275 | z = fData[2]->GetValue()*fScale[2];
|
---|
276 | case 2:
|
---|
277 | y = fData[1]->GetValue()*fScale[1];
|
---|
278 | case 1:
|
---|
279 | x = fData[0]->GetValue()*fScale[0];
|
---|
280 | }
|
---|
281 |
|
---|
282 | switch (fDimension)
|
---|
283 | {
|
---|
284 | case 3:
|
---|
285 | ((TH3*)fHist)->Fill(x, y, z);
|
---|
286 | return kTRUE;
|
---|
287 | case 2:
|
---|
288 | ((TH2*)fHist)->Fill(x, y);
|
---|
289 | return kTRUE;
|
---|
290 | case 1:
|
---|
291 | fHist->Fill(x);
|
---|
292 | return kTRUE;
|
---|
293 | }
|
---|
294 |
|
---|
295 | return kFALSE;
|
---|
296 | }
|
---|
297 |
|
---|
298 | // --------------------------------------------------------------------------
|
---|
299 | //
|
---|
300 | // Draw clone of histogram. So that the object can be deleted
|
---|
301 | // and the histogram is still visible in the canvas.
|
---|
302 | // The cloned object are deleted together with the canvas if the canvas is
|
---|
303 | // destroyed. If you want to handle destroying the canvas you can get a
|
---|
304 | // pointer to it from this function
|
---|
305 | //
|
---|
306 | TObject *MH3::DrawClone(Option_t *opt) const
|
---|
307 | {
|
---|
308 | TCanvas &c = *MH::MakeDefCanvas(fHist);
|
---|
309 |
|
---|
310 | //
|
---|
311 | // This is necessary to get the expected bahviour of DrawClone
|
---|
312 | //
|
---|
313 | gROOT->SetSelectedPad(NULL);
|
---|
314 |
|
---|
315 | fHist->DrawCopy(opt);
|
---|
316 |
|
---|
317 | TString str(opt);
|
---|
318 | if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
|
---|
319 | {
|
---|
320 | TProfile *p = ((TH2*)fHist)->ProfileX();
|
---|
321 | p->Draw("same");
|
---|
322 | p->SetBit(kCanDelete);
|
---|
323 | }
|
---|
324 | if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
|
---|
325 | {
|
---|
326 | TProfile *p = ((TH2*)fHist)->ProfileY();
|
---|
327 | p->Draw("same");
|
---|
328 | p->SetBit(kCanDelete);
|
---|
329 | }
|
---|
330 |
|
---|
331 | c.Modified();
|
---|
332 | c.Update();
|
---|
333 |
|
---|
334 | return &c;
|
---|
335 | }
|
---|
336 |
|
---|
337 | // --------------------------------------------------------------------------
|
---|
338 | //
|
---|
339 | // Creates a new canvas and draws the histogram into it.
|
---|
340 | // Be careful: The histogram belongs to this object and won't get deleted
|
---|
341 | // together with the canvas.
|
---|
342 | //
|
---|
343 | void MH3::Draw(Option_t *opt)
|
---|
344 | {
|
---|
345 | if (!gPad)
|
---|
346 | MH::MakeDefCanvas(fHist);
|
---|
347 |
|
---|
348 | fHist->Draw(opt);
|
---|
349 |
|
---|
350 | TString str(opt);
|
---|
351 | if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
|
---|
352 | {
|
---|
353 | TProfile *p = ((TH2*)fHist)->ProfileX();
|
---|
354 | p->Draw("same");
|
---|
355 | p->SetBit(kCanDelete);
|
---|
356 | }
|
---|
357 | if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
|
---|
358 | {
|
---|
359 | TProfile *p = ((TH2*)fHist)->ProfileY();
|
---|
360 | p->Draw("same");
|
---|
361 | p->SetBit(kCanDelete);
|
---|
362 | }
|
---|
363 |
|
---|
364 | gPad->Modified();
|
---|
365 | gPad->Update();
|
---|
366 | }
|
---|