/* ======================================================================== *\ ! ! * ! * 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 // tolower #include #include #include #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MParList.h" #include "MBinning.h" #include "MDataChain.h" ClassImp(MH3); using namespace std; static const TString gsDefName = "MH3"; static const TString gsDefTitle = "Container for a %dD Mars Histogram"; // -------------------------------------------------------------------------- // // Default constructor. // MH3::MH3(const unsigned int dim) : fDimension(dim>3?3:dim), fHist(NULL) { switch (fDimension) { case 1: fHist = new TH1F; fHist->SetYTitle("Counts"); break; case 2: fHist = new TH2F; fHist->SetZTitle("Counts"); break; case 3: fHist = new TH3F; break; } fData[0] = NULL; fData[1] = NULL; fData[2] = NULL; fName = gsDefName; fTitle = Form(gsDefTitle.Data(), 1); if (fHist) { fHist->SetDirectory(NULL); fHist->UseCurrentStyle(); } fScale[0] = 1; fScale[1] = 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->UseCurrentStyle(); 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->UseCurrentStyle(); 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->UseCurrentStyle(); 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]; } // -------------------------------------------------------------------------- // // Return the data members used by the data chain to be used in // MTask::AddBranchToList // TString MH3::GetDataMember() const { TString str=fData[0]->GetDataMember(); if (fData[1]) { str += ";"; str += fData[1]->GetDataMember(); } if (fData[2]) { str += ";"; str += fData[2]->GetDataMember(); } return str; } // -------------------------------------------------------------------------- // // 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) { // reset histogram (necessary if the same eventloop is run more than once) fHist->Reset(); TString bname("Binning"); bname += fName; MBinning *binsx = NULL; MBinning *binsy = NULL; MBinning *binsz = NULL; switch (fDimension) { case 3: binsz = (MBinning*)plist->FindObject(bname+"Z", "MBinning"); if (!binsz) { *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl; return kFALSE; } if (binsz->IsLogarithmic()) fHist->SetBit(kIsLogz); if (fData[2]) fHist->SetZTitle(fData[2]->GetTitle()); if (fData[2] && !fData[2]->PreProcess(plist)) return kFALSE; case 2: binsy = (MBinning*)plist->FindObject(bname+"Y", "MBinning"); if (!binsy) { *fLog << err << dbginf << "MBinning '" << bname << "Y' not found... aborting." << endl; return kFALSE; } if (binsy->IsLogarithmic()) fHist->SetBit(kIsLogy); if (fData[1]) fHist->SetYTitle(fData[1]->GetTitle()); if (fData[1] && !fData[1]->PreProcess(plist)) return kFALSE; case 1: binsx = (MBinning*)plist->FindObject(bname+"X", "MBinning"); if (!binsx) { if (fDimension==1) binsx = (MBinning*)plist->FindObject(bname, "MBinning"); if (!binsx) { *fLog << err << dbginf << "Neither MBinning '" << bname << "X' nor '" << bname << "' found... aborting." << endl; return kFALSE; } } if (binsx->IsLogarithmic()) fHist->SetBit(kIsLogx); if (fData[0]!=NULL) fHist->SetXTitle(fData[0]->GetTitle()); if (fData[0] && !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; } cout << "Still alive...?" << endl; 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, const Stat_t w) { 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, w); return kTRUE; case 2: ((TH2*)fHist)->Fill(x, y, w); return kTRUE; case 1: fHist->Fill(x, w); return kTRUE; } return kFALSE; } /* // -------------------------------------------------------------------------- // // Set the palette you wanna use: // - you could set the root "Pretty Palette Violet->Red" by // gStyle->SetPalette(1, 0), but in some cases this may look // confusing // - The maximum colors root allowes us to set by ourself // is 50 (idx: 51-100). This colors are set to a grayscaled // palette // - the number of contours must be two less than the number // of palette entries // void MHStarMap::PrepareDrawing() const { const Int_t numg = 32; // number of gray scaled colors const Int_t numw = 32; // number of white Int_t palette[numg+numw]; // // The first half of the colors are white. // This is some kind of optical background supression // gROOT->GetColor(51)->SetRGB(1, 1, 1); Int_t i; for (i=0; iGetColor(52+i)->SetRGB(gray, gray, gray); palette[i] = 52+i; } // // Set the palette and the number of contour levels // gStyle->SetPalette(numg+numw, palette); fStarMap->SetContour(numg+numw-2); } */ // -------------------------------------------------------------------------- // // Setup a inversed deep blue sea palette for the fCenter histogram. // void MH3::SetColors() const { // FIXME: This must be redone each time the canvas is repainted.... gStyle->SetPalette(51, NULL); Int_t c[50]; for (int i=0; i<50; i++) c[49-i] = gStyle->GetColorPalette(i); gStyle->SetPalette(50, c); } // -------------------------------------------------------------------------- // // Draw clone of histogram. So that the object can be deleted // // Possible options are: // PROFX: Draw a x-profile into the histogram (for 2D histograms only) // PROFY: Draw a y-profile into the histogram (for 2D histograms only) // ONLY: Draw the profile histogram only (for 2D histograms only) // // If the kIsLog?-Bit is set the axis is displayed lkogarithmically. // eg this is set when applying a logarithmic MBinning // // 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 { TString str(opt); TVirtualPad *c = gPad; if (!str.Contains("nonew", TString::kIgnoreCase)) { c = MH::MakeDefCanvas(fHist); // // This is necessary to get the expected bahviour of DrawClone // gROOT->SetSelectedPad(NULL); } if (str.Contains("COL", TString::kIgnoreCase)) SetColors(); Bool_t only = str.Contains("ONLY", TString::kIgnoreCase) && fDimension==2; if (!only) fHist->DrawCopy(opt); if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileX("_pfx", -1, 9999, "s"); p->SetLineColor(kBlue); p->Draw(only?"":"same"); p->SetBit(kCanDelete); p->SetDirectory(NULL); } if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileY("_pfy", -1, 9999, "s"); p->SetLineColor(kBlue); p->Draw(only?"":"same"); p->SetBit(kCanDelete); p->SetDirectory(NULL); } if (fHist->TestBit(kIsLogx)) c->SetLogx(); if (fHist->TestBit(kIsLogy)) c->SetLogy(); if (fHist->TestBit(kIsLogz)) c->SetLogz(); c->Modified(); c->Update(); return c; } */ // -------------------------------------------------------------------------- // // Creates a new canvas and draws the histogram into it. // // Possible options are: // PROFX: Draw a x-profile into the histogram (for 2D histograms only) // PROFY: Draw a y-profile into the histogram (for 2D histograms only) // ONLY: Draw the profile histogram only (for 2D histograms only) // // If the kIsLog?-Bit is set the axis is displayed lkogarithmically. // eg this is set when applying a logarithmic MBinning // // Be careful: The histogram belongs to this object and won't get deleted // together with the canvas. // void MH3::Draw(Option_t *opt) { TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fHist); pad->SetBorderMode(0); AppendPad(""); TString str(opt); // FIXME: Do it in Paint() if (str.Contains("COL", TString::kIgnoreCase)) SetColors(); fHist->SetFillStyle(4000); Bool_t only = str.Contains("ONLY", TString::kIgnoreCase) && fDimension==2; if (!only) fHist->Draw(opt); if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileX("_pfx", -1, 9999, "s"); p->SetLineColor(kBlue); p->Draw(only?"":"same"); p->SetBit(kCanDelete); p->SetDirectory(NULL); } if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileY("_pfy", -1, 9999, "s"); p->SetLineColor(kBlue); p->Draw(only?"":"same"); p->SetBit(kCanDelete); p->SetDirectory(NULL); } if (fHist->TestBit(kIsLogx)) pad->SetLogx(); if (fHist->TestBit(kIsLogy)) pad->SetLogy(); if (fHist->TestBit(kIsLogz)) pad->SetLogz(); pad->Modified(); pad->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; } } // -------------------------------------------------------------------------- // // Used to rebuild a MH3 object of the same type (data members, // dimension, ...) // MParContainer *MH3::New() const { MH3 *h = NULL; if (fData[0] == NULL) h=new MH3(fDimension); else switch (fDimension) { case 1: h=new MH3(fData[0]->GetRule()); break; case 2: h=new MH3(fData[0]->GetRule(), fData[1]->GetRule()); break; case 3: h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule()); break; } switch (fDimension) { case 3: h->SetScaleZ(fScale[2]); case 2: h->SetScaleY(fScale[1]); case 1: h->SetScaleX(fScale[0]); } return h; } TString MH3::GetRule(const Char_t axis) const { switch (tolower(axis)) { case 'x': return fData[0] ? fData[0]->GetRule() : TString(""); case 'y': return fData[1] ? fData[1]->GetRule() : TString(""); case 'z': return fData[2] ? fData[2]->GetRule() : TString(""); default: return ""; } } // -------------------------------------------------------------------------- // // Returns the total number of bins in a histogram (excluding under- and // overflow bins) // Int_t MH3::GetNbins() const { Int_t num = 1; switch (fDimension) { case 3: num *= fHist->GetNbinsZ()+2; case 2: num *= fHist->GetNbinsY()+2; case 1: num *= fHist->GetNbinsX()+2; } return num; } // -------------------------------------------------------------------------- // // Returns the total number of bins in a histogram (excluding under- and // overflow bins) Return -1 if bin is underflow or overflow // Int_t MH3::FindFixBin(Double_t x, Double_t y, Double_t z) const { const TAxis &axex = *fHist->GetXaxis(); const TAxis &axey = *fHist->GetYaxis(); const TAxis &axez = *fHist->GetZaxis(); Int_t binz = 0; Int_t biny = 0; Int_t binx = 0; switch (fDimension) { case 3: binz = axez.FindFixBin(z); if (binz>axez.GetLast() || binzaxey.GetLast() || binyaxex.GetLast()) return -1; } const Int_t nx = fHist->GetNbinsX()+2; const Int_t ny = fHist->GetNbinsY()+2; return binx + nx*(biny +ny*binz); }