/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz 2002 ! ! Copyright: MAGIC Software Development, 2000-2002 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MH3 // // With this histogram you can fill a histogram with up to three // variables from Mars parameter containers in an eventloop. // // In the constructor you can give up to three variables which should be // filled in the histogram. Dependend on the number of given variables // (data members) a TH1F, TH2F or TH3F is created. // Specify the data mamber like the following: // "MHillas.fLength" // Where MHillas is the name of the parameter container in the parameter // list and fLength is the name of the data member which should be filled // in the histogram. Assuming that your MHillas container has a different // name (MyHillas) the name to give would be: // "MyHillas.fLength" // // The axis binning is retrieved from the parameter list, too. Create a // MBinning with the name "Binning" plus the name of your MH3 container // plus the axis name ("X", "Y" or "Z") and add it to the parameter list. // // If you want to use a different unit for histogramming use SetScaleX, // SetScaleY and SetScaleZ. // // For example: // MH3 myhist("MHillas.fLength"); // myhist.SetName("MyHist"); // myhist.SetScaleX(geomcam.GetConvMm2Deg()); //convert length to degree // MBinning bins("BinningMyHistX"); // bins.SetEdges(10, 0, 150); // plist.AddToList(&myhist); // plist.AddToList(&bins); // ///////////////////////////////////////////////////////////////////////////// #include "MH3.h" #include #include #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MParList.h" #include "MDataChain.h" ClassImp(MH3); static const TString gsDefName = "MH3"; static const TString gsDefTitle = "Container for a %dD Mars Histogram"; MH3::MH3() : fDimension(0), fHist(NULL) { fName = gsDefName; fTitle = Form(gsDefTitle.Data(), 0); fData[0] = fData[1] = fData[2] = NULL; fScale[0] = fScale[1] = fScale[2] = 1; } // -------------------------------------------------------------------------- // // Creates an TH1F. memberx is filled into the X-bins. For a more detailed // description see the class description above. // MH3::MH3(const char *memberx) : fDimension(1) { fHist = new TH1F; fData[0] = new MDataChain(memberx); fData[1] = NULL; fData[2] = NULL; fName = gsDefName; fTitle = Form(gsDefTitle.Data(), 1); fHist->SetDirectory(NULL); fHist->SetYTitle("Counts"); fScale[0] = 1; fScale[1] = 1; fScale[2] = 1; } // -------------------------------------------------------------------------- // // Creates an TH2F. memberx is filled into the X-bins. membery is filled // into the Y-bins. For a more detailed description see the class // description above. // MH3::MH3(const char *memberx, const char *membery) : fDimension(2) { fHist = new TH2F; fData[0] = new MDataChain(memberx); fData[1] = new MDataChain(membery); fData[2] = NULL; fName = gsDefName; fTitle = Form(gsDefTitle.Data(), 2); fHist->SetDirectory(NULL); fHist->SetZTitle("Counts"); fScale[0] = 1; fScale[1] = 1; fScale[2] = 1; } // -------------------------------------------------------------------------- // // Creates an TH3F. memberx is filled into the X-bins. membery is filled // into the Y-bins. membery is filled into the Z-bins. For a more detailed // description see the class description above. // MH3::MH3(const char *memberx, const char *membery, const char *memberz) : fDimension(3) { fHist = new TH3F; fData[0] = new MDataChain(memberx); fData[1] = new MDataChain(membery); fData[2] = new MDataChain(memberz); fName = gsDefName; fTitle = Form(gsDefTitle.Data(), 3); fHist->SetDirectory(NULL); fScale[0] = 1; fScale[1] = 1; fScale[2] = 1; } // -------------------------------------------------------------------------- // // Deletes the histogram // MH3::~MH3() { delete fHist; for (int i=0; i<3; i++) if (fData[i]) delete fData[i]; } // -------------------------------------------------------------------------- // // Setup the Binning for the histograms automatically if the correct // instances of MBinning are found in the parameter list // For a more detailed description see class description above. // Bool_t MH3::SetupFill(const MParList *plist) { TString bname("Binning"); bname += fName; MBinning *binsx = NULL; MBinning *binsy = NULL; MBinning *binsz = NULL; switch (fDimension) { case 3: binsz = (MBinning*)plist->FindObject(bname+"Z"); if (!binsz) { *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl; return kFALSE; } fHist->SetZTitle(fData[2]->GetTitle()); if (!fData[2]->PreProcess(plist)) return kFALSE; case 2: binsy = (MBinning*)plist->FindObject(bname+"Y"); if (!binsy) { *fLog << err << dbginf << "MBinning '" << bname << "Y' not found... aborting." << endl; return kFALSE; } fHist->SetYTitle(fData[1]->GetTitle()); if (!fData[1]->PreProcess(plist)) return kFALSE; case 1: binsx = (MBinning*)plist->FindObject(bname+"X"); if (!binsx) { *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl; return kFALSE; } fHist->SetXTitle(fData[0]->GetTitle()); if (!fData[0]->PreProcess(plist)) return kFALSE; } fHist->SetName(fName); TString title("Histogram for "); title += fName; switch (fDimension) { case 1: fHist->SetTitle(title+" (1D)"); SetBinning(fHist, binsx); return kTRUE; case 2: fHist->SetTitle(title+" (2D)"); SetBinning((TH2*)fHist, binsx, binsy); return kTRUE; case 3: fHist->SetTitle(title+" (3D)"); SetBinning((TH3*)fHist, binsx, binsy, binsz); return kTRUE; } return kTRUE; } // -------------------------------------------------------------------------- // // Set the name of the histogram ant the MH3 container // void MH3::SetName(const char *name) { fHist->SetName(name); MParContainer::SetName(name); } // -------------------------------------------------------------------------- // // Set the title of the histogram ant the MH3 container // void MH3::SetTitle(const char *title) { fHist->SetTitle(title); MParContainer::SetTitle(title); } // -------------------------------------------------------------------------- // // Fills the one, two or three data members into our histogram // Bool_t MH3::Fill(const MParContainer *par) { Double_t x=0; Double_t y=0; Double_t z=0; switch (fDimension) { case 3: z = fData[2]->GetValue()*fScale[2]; case 2: y = fData[1]->GetValue()*fScale[1]; case 1: x = fData[0]->GetValue()*fScale[0]; } switch (fDimension) { case 3: ((TH3*)fHist)->Fill(x, y, z); return kTRUE; case 2: ((TH2*)fHist)->Fill(x, y); return kTRUE; case 1: fHist->Fill(x); return kTRUE; } return kFALSE; } // -------------------------------------------------------------------------- // // Draw clone of histogram. So that the object can be deleted // and the histogram is still visible in the canvas. // The cloned object are deleted together with the canvas if the canvas is // destroyed. If you want to handle destroying the canvas you can get a // pointer to it from this function // TObject *MH3::DrawClone(Option_t *opt) const { TCanvas &c = *MH::MakeDefCanvas(fHist); // // This is necessary to get the expected bahviour of DrawClone // gROOT->SetSelectedPad(NULL); fHist->DrawCopy(opt); TString str(opt); if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileX(); p->Draw("same"); p->SetBit(kCanDelete); p->SetDirectory(NULL); } if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileY(); p->Draw("same"); p->SetBit(kCanDelete); p->SetDirectory(NULL); } c.Modified(); c.Update(); return &c; } // -------------------------------------------------------------------------- // // Creates a new canvas and draws the histogram into it. // Be careful: The histogram belongs to this object and won't get deleted // together with the canvas. // void MH3::Draw(Option_t *opt) { if (!gPad) MH::MakeDefCanvas(fHist); fHist->Draw(opt); TString str(opt); if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileX(); p->Draw("same"); p->SetBit(kCanDelete); p->SetDirectory(NULL); } if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileY(); p->Draw("same"); p->SetBit(kCanDelete); p->SetDirectory(NULL); } gPad->Modified(); gPad->Update(); } // -------------------------------------------------------------------------- // // Implementation of SavePrimitive. Used to write the call to a constructor // to a macro. In the original root implementation it is used to write // gui elements to a macro-file. // void MH3::StreamPrimitive(ofstream &out) const { TString name = GetUniqueName(); out << " MH3 " << name << "(\""; out << fData[0]->GetRule() << "\""; if (fDimension>1) out << ", \"" << fData[1]->GetRule() << "\""; if (fDimension>2) out << ", \"" << fData[2]->GetRule() << "\""; out << ");" << endl; if (fName!=gsDefName) out << " " << name << ".SetName(\"" << fName << "\");" << endl; if (fTitle!=Form(gsDefTitle.Data(), fDimension)) out << " " << name << ".SetTitle(\"" << fTitle << "\");" << endl; switch (fDimension) { case 3: if (fScale[2]!=1) out << " " << name << ".SetScaleZ(" << fScale[2] << ");" << endl; case 2: if (fScale[1]!=1) out << " " << name << ".SetScaleY(" << fScale[1] << ");" << endl; case 1: if (fScale[0]!=1) out << " " << name << ".SetScaleX(" << fScale[0] << ");" << endl; } }