| 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 <TPad.h>
|
|---|
| 66 | #include <TCanvas.h>
|
|---|
| 67 | #include <TDataMember.h>
|
|---|
| 68 | #include <TMethodCall.h>
|
|---|
| 69 |
|
|---|
| 70 | #include "MLog.h"
|
|---|
| 71 | #include "MLogManip.h"
|
|---|
| 72 |
|
|---|
| 73 | #include "MParList.h"
|
|---|
| 74 |
|
|---|
| 75 | ClassImp(MH3);
|
|---|
| 76 |
|
|---|
| 77 | // --------------------------------------------------------------------------
|
|---|
| 78 | //
|
|---|
| 79 | // Creates an TH1F. memberx is filled into the X-bins. For a more detailed
|
|---|
| 80 | // description see the class description above.
|
|---|
| 81 | //
|
|---|
| 82 | MH3::MH3(const char *memberx) : fDimension(1)
|
|---|
| 83 | {
|
|---|
| 84 | fHist = new TH1F;
|
|---|
| 85 |
|
|---|
| 86 | fDataMember[0] = memberx;
|
|---|
| 87 |
|
|---|
| 88 | fName = "MH3";
|
|---|
| 89 | fTitle = "Container for a 1D Mars Histogram";
|
|---|
| 90 |
|
|---|
| 91 | fHist->SetDirectory(NULL);
|
|---|
| 92 | fHist->SetYTitle("Counts");
|
|---|
| 93 |
|
|---|
| 94 | fScale[0] = 1;
|
|---|
| 95 | fScale[1] = 1;
|
|---|
| 96 | fScale[2] = 1;
|
|---|
| 97 | }
|
|---|
| 98 |
|
|---|
| 99 | // --------------------------------------------------------------------------
|
|---|
| 100 | //
|
|---|
| 101 | // Creates an TH2F. memberx is filled into the X-bins. membery is filled
|
|---|
| 102 | // into the Y-bins. For a more detailed description see the class
|
|---|
| 103 | // description above.
|
|---|
| 104 | //
|
|---|
| 105 | MH3::MH3(const char *memberx, const char *membery) : fDimension(2)
|
|---|
| 106 | {
|
|---|
| 107 | fHist = new TH2F;
|
|---|
| 108 |
|
|---|
| 109 | fDataMember[0] = memberx;
|
|---|
| 110 | fDataMember[1] = membery;
|
|---|
| 111 |
|
|---|
| 112 | fName = "MH3";
|
|---|
| 113 | fTitle = "Container for a 2D Mars Histogram";
|
|---|
| 114 |
|
|---|
| 115 | fHist->SetDirectory(NULL);
|
|---|
| 116 | fHist->SetZTitle("Counts");
|
|---|
| 117 |
|
|---|
| 118 | fScale[0] = 1;
|
|---|
| 119 | fScale[1] = 1;
|
|---|
| 120 | fScale[2] = 1;
|
|---|
| 121 | }
|
|---|
| 122 |
|
|---|
| 123 | // --------------------------------------------------------------------------
|
|---|
| 124 | //
|
|---|
| 125 | // Creates an TH3F. memberx is filled into the X-bins. membery is filled
|
|---|
| 126 | // into the Y-bins. membery is filled into the Z-bins. For a more detailed
|
|---|
| 127 | // description see the class description above.
|
|---|
| 128 | //
|
|---|
| 129 | MH3::MH3(const char *memberx, const char *membery, const char *memberz) : fDimension(3)
|
|---|
| 130 | {
|
|---|
| 131 | fHist = new TH3F;
|
|---|
| 132 |
|
|---|
| 133 | fDataMember[0] = memberx;
|
|---|
| 134 | fDataMember[1] = membery;
|
|---|
| 135 | fDataMember[2] = memberz;
|
|---|
| 136 |
|
|---|
| 137 | fName = "MH3";
|
|---|
| 138 | fTitle = "Container for a 3D Mars Histogram";
|
|---|
| 139 |
|
|---|
| 140 | fScale[0] = 1;
|
|---|
| 141 | fScale[1] = 1;
|
|---|
| 142 | fScale[2] = 1;
|
|---|
| 143 | }
|
|---|
| 144 |
|
|---|
| 145 | // --------------------------------------------------------------------------
|
|---|
| 146 | //
|
|---|
| 147 | // Deletes the histogram
|
|---|
| 148 | //
|
|---|
| 149 | MH3::~MH3()
|
|---|
| 150 | {
|
|---|
| 151 | delete fHist;
|
|---|
| 152 | }
|
|---|
| 153 |
|
|---|
| 154 | // --------------------------------------------------------------------------
|
|---|
| 155 | //
|
|---|
| 156 | // Trys to determin the TMethodCall corresponding to the num-axis
|
|---|
| 157 | // corresponding data member.
|
|---|
| 158 | //
|
|---|
| 159 | Bool_t MH3::GetMethodCall(const MParList *plist, Int_t num)
|
|---|
| 160 | {
|
|---|
| 161 | TString cname(fDataMember[num]);
|
|---|
| 162 | TString mname(fDataMember[num]);
|
|---|
| 163 |
|
|---|
| 164 | const char *dot = strrchr(cname, '.');
|
|---|
| 165 |
|
|---|
| 166 | if (dot)
|
|---|
| 167 | {
|
|---|
| 168 | const int pos = dot-cname;
|
|---|
| 169 |
|
|---|
| 170 | cname.Remove(pos);
|
|---|
| 171 | mname.Remove(0, pos+1);
|
|---|
| 172 | }
|
|---|
| 173 |
|
|---|
| 174 | fObject[num] = (MParContainer*)plist->FindObject(cname);
|
|---|
| 175 | if (!fObject[num])
|
|---|
| 176 | {
|
|---|
| 177 | *fLog << err << "Object '" << cname << "' not in parameter list... aborting." << endl;
|
|---|
| 178 | return kFALSE;
|
|---|
| 179 | }
|
|---|
| 180 |
|
|---|
| 181 | fMethodCall[num] = fObject[num]->GetterMethod(mname);
|
|---|
| 182 |
|
|---|
| 183 | return fMethodCall[num] ? kTRUE : kFALSE;
|
|---|
| 184 | }
|
|---|
| 185 |
|
|---|
| 186 |
|
|---|
| 187 | // --------------------------------------------------------------------------
|
|---|
| 188 | //
|
|---|
| 189 | // Setup the Binning for the histograms automatically if the correct
|
|---|
| 190 | // instances of MBinning are found in the parameter list
|
|---|
| 191 | // For a more detailed description see class description above.
|
|---|
| 192 | //
|
|---|
| 193 | Bool_t MH3::SetupFill(const MParList *plist)
|
|---|
| 194 | {
|
|---|
| 195 | TString bname("Binning");
|
|---|
| 196 | bname += fName;
|
|---|
| 197 |
|
|---|
| 198 | MBinning *binsx = NULL;
|
|---|
| 199 | MBinning *binsy = NULL;
|
|---|
| 200 | MBinning *binsz = NULL;
|
|---|
| 201 | switch (fDimension)
|
|---|
| 202 | {
|
|---|
| 203 | case 3:
|
|---|
| 204 | binsz = (MBinning*)plist->FindObject(bname+"Z");
|
|---|
| 205 | if (!binsz)
|
|---|
| 206 | {
|
|---|
| 207 | *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl;
|
|---|
| 208 | return kFALSE;
|
|---|
| 209 | }
|
|---|
| 210 | fHist->SetZTitle(fName+"Z");
|
|---|
| 211 | if (!GetMethodCall(plist, 2))
|
|---|
| 212 | return kFALSE;
|
|---|
| 213 | case 2:
|
|---|
| 214 | binsy = (MBinning*)plist->FindObject(bname+"Y");
|
|---|
| 215 | if (!binsy)
|
|---|
| 216 | {
|
|---|
| 217 | *fLog << err << dbginf << "MBinning '" << bname << "Y' not found... aborting." << endl;
|
|---|
| 218 | return kFALSE;
|
|---|
| 219 | }
|
|---|
| 220 | fHist->SetYTitle(fName+"Y");
|
|---|
| 221 | if (!GetMethodCall(plist, 1))
|
|---|
| 222 | return kFALSE;
|
|---|
| 223 | case 1:
|
|---|
| 224 | binsx = (MBinning*)plist->FindObject(bname+"X");
|
|---|
| 225 | if (!binsx)
|
|---|
| 226 | {
|
|---|
| 227 | *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl;
|
|---|
| 228 | return kFALSE;
|
|---|
| 229 | }
|
|---|
| 230 | fHist->SetXTitle(fName+"X");
|
|---|
| 231 | if (!GetMethodCall(plist, 0))
|
|---|
| 232 | return kFALSE;
|
|---|
| 233 | }
|
|---|
| 234 |
|
|---|
| 235 | fHist->SetName(fName);
|
|---|
| 236 |
|
|---|
| 237 | TString title("Histogram for ");
|
|---|
| 238 | title += fName;
|
|---|
| 239 |
|
|---|
| 240 | switch (fDimension)
|
|---|
| 241 | {
|
|---|
| 242 | case 1:
|
|---|
| 243 | fHist->SetTitle(title+" (1D)");
|
|---|
| 244 | SetBinning(fHist, binsx);
|
|---|
| 245 | return kTRUE;
|
|---|
| 246 | case 2:
|
|---|
| 247 | fHist->SetTitle(title+" (2D)");
|
|---|
| 248 | SetBinning(fHist, binsx, binsy);
|
|---|
| 249 | return kTRUE;
|
|---|
| 250 | case 3:
|
|---|
| 251 | fHist->SetTitle(title+" (3D)");
|
|---|
| 252 | SetBinning(fHist, binsx, binsy, binsz);
|
|---|
| 253 | return kTRUE;
|
|---|
| 254 | }
|
|---|
| 255 |
|
|---|
| 256 | return kTRUE;
|
|---|
| 257 | }
|
|---|
| 258 |
|
|---|
| 259 | // --------------------------------------------------------------------------
|
|---|
| 260 | //
|
|---|
| 261 | // Set the name of the histogram ant the MH3 container
|
|---|
| 262 | //
|
|---|
| 263 | void MH3::SetName(const char *name)
|
|---|
| 264 | {
|
|---|
| 265 | fHist->SetName(name);
|
|---|
| 266 | MParContainer::SetName(name);
|
|---|
| 267 | }
|
|---|
| 268 |
|
|---|
| 269 | // --------------------------------------------------------------------------
|
|---|
| 270 | //
|
|---|
| 271 | // Set the title of the histogram ant the MH3 container
|
|---|
| 272 | //
|
|---|
| 273 | void MH3::SetTitle(const char *title)
|
|---|
| 274 | {
|
|---|
| 275 | fHist->SetTitle(title);
|
|---|
| 276 | MParContainer::SetTitle(title);
|
|---|
| 277 | }
|
|---|
| 278 |
|
|---|
| 279 | // --------------------------------------------------------------------------
|
|---|
| 280 | //
|
|---|
| 281 | // Returns the value corresponding to the data member of the given object
|
|---|
| 282 | //
|
|---|
| 283 | Bool_t MH3::GetValue(Int_t num, Double_t &v)
|
|---|
| 284 | {
|
|---|
| 285 | switch (fMethodCall[num]->ReturnType())
|
|---|
| 286 | {
|
|---|
| 287 | case TMethodCall::kLong:
|
|---|
| 288 | Long_t l;
|
|---|
| 289 | fMethodCall[num]->Execute(fObject[num], l); // FIXME: const, root
|
|---|
| 290 | v = fScale[num]*l;
|
|---|
| 291 | return kTRUE;
|
|---|
| 292 |
|
|---|
| 293 | case TMethodCall::kDouble:
|
|---|
| 294 | fMethodCall[num]->Execute(fObject[num], v); // FIXME: const, root
|
|---|
| 295 | v *= fScale[num];
|
|---|
| 296 | return kTRUE;
|
|---|
| 297 |
|
|---|
| 298 | default:
|
|---|
| 299 | *fLog << err << "DataMember " << fDataMember[num] << " of ";
|
|---|
| 300 | *fLog << fObject[num]->GetName() << " neither int nor float... abort." << endl;
|
|---|
| 301 | return kFALSE;
|
|---|
| 302 | }
|
|---|
| 303 | }
|
|---|
| 304 |
|
|---|
| 305 | // --------------------------------------------------------------------------
|
|---|
| 306 | //
|
|---|
| 307 | // Fills the one, two or three data members into our histogram
|
|---|
| 308 | //
|
|---|
| 309 | Bool_t MH3::Fill(const MParContainer *par)
|
|---|
| 310 | {
|
|---|
| 311 | Double_t x, y, z;
|
|---|
| 312 |
|
|---|
| 313 | switch (fDimension)
|
|---|
| 314 | {
|
|---|
| 315 | case 3:
|
|---|
| 316 | if (!GetValue(2, z))
|
|---|
| 317 | return kFALSE;
|
|---|
| 318 | case 2:
|
|---|
| 319 | if (!GetValue(1, y))
|
|---|
| 320 | return kFALSE;
|
|---|
| 321 | case 1:
|
|---|
| 322 | if (!GetValue(0, x))
|
|---|
| 323 | return kFALSE;
|
|---|
| 324 | }
|
|---|
| 325 |
|
|---|
| 326 | switch (fDimension)
|
|---|
| 327 | {
|
|---|
| 328 | case 3:
|
|---|
| 329 | ((TH3*)fHist)->Fill(x, y, z);
|
|---|
| 330 | return kTRUE;
|
|---|
| 331 | case 2:
|
|---|
| 332 | ((TH2*)fHist)->Fill(x, y);
|
|---|
| 333 | return kTRUE;
|
|---|
| 334 | case 1:
|
|---|
| 335 | fHist->Fill(x);
|
|---|
| 336 | return kTRUE;
|
|---|
| 337 | }
|
|---|
| 338 |
|
|---|
| 339 | return kFALSE;
|
|---|
| 340 | }
|
|---|
| 341 |
|
|---|
| 342 | // --------------------------------------------------------------------------
|
|---|
| 343 | //
|
|---|
| 344 | // Draw clone of histogram. So that the object can be deleted
|
|---|
| 345 | // and the histogram is still visible in the canvas.
|
|---|
| 346 | // The cloned object are deleted together with the canvas if the canvas is
|
|---|
| 347 | // destroyed. If you want to handle destroying the canvas you can get a
|
|---|
| 348 | // pointer to it from this function
|
|---|
| 349 | //
|
|---|
| 350 | TObject *MH3::DrawClone(Option_t *opt) const
|
|---|
| 351 | {
|
|---|
| 352 | TCanvas &c = *MH::MakeDefCanvas(fHist);
|
|---|
| 353 |
|
|---|
| 354 | //
|
|---|
| 355 | // This is necessary to get the expected bahviour of DrawClone
|
|---|
| 356 | //
|
|---|
| 357 | gROOT->SetSelectedPad(NULL);
|
|---|
| 358 |
|
|---|
| 359 | fHist->DrawCopy(opt);
|
|---|
| 360 |
|
|---|
| 361 | c.Modified();
|
|---|
| 362 | c.Update();
|
|---|
| 363 |
|
|---|
| 364 | return &c;
|
|---|
| 365 | }
|
|---|
| 366 |
|
|---|
| 367 | // --------------------------------------------------------------------------
|
|---|
| 368 | //
|
|---|
| 369 | // Creates a new canvas and draws the histogram into it.
|
|---|
| 370 | // Be careful: The histogram belongs to this object and won't get deleted
|
|---|
| 371 | // together with the canvas.
|
|---|
| 372 | //
|
|---|
| 373 | void MH3::Draw(Option_t *opt)
|
|---|
| 374 | {
|
|---|
| 375 | if (!gPad)
|
|---|
| 376 | MH::MakeDefCanvas(fHist);
|
|---|
| 377 |
|
|---|
| 378 | fHist->Draw(opt);
|
|---|
| 379 |
|
|---|
| 380 | gPad->Modified();
|
|---|
| 381 | gPad->Update();
|
|---|
| 382 | }
|
|---|