/* ======================================================================== *\ ! ! * ! * 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 2001 ! ! Copyright: MAGIC Software Development, 2000-2001 ! ! \* ======================================================================== */ /////////////////////////////////////////////////////////////////////// // // MHHillasSrc // // This class contains histograms for every Hillas parameter // /////////////////////////////////////////////////////////////////////// #include "MHHillasSrc.h" #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MGeomCam.h" #include "MParList.h" #include "MHillas.h" #include "MHillasSrc.h" ClassImp(MHHillasSrc); using namespace std; // -------------------------------------------------------------------------- // // Setup four histograms for Alpha, and Dist // MHHillasSrc::MHHillasSrc(const char *name, const char *title) : fUseMmScale(kTRUE) { // // set the name and title of this object // fName = name ? name : "MHHillasSrc"; fTitle = title ? title : "Container for Hillas histograms"; // // loop over all Pixels and create two histograms // one for the Low and one for the High gain // connect all the histogram with the container fHist // fAlpha = new TH1F("Alpha", "Alpha of Ellipse", 90, -90, 90); fDist = new TH1F("Dist", "Dist of Ellipse", 70, 0, 623); fCosDA = new TH1F("CosDA", "cos(Delta,Alpha) of Ellipse", 101, -1, 1); fDCA = new TH1F("DCA", "Distance of closest aproach", 101, -500, 500); fDCADelta = new TH1F("DCADelta", "Angle between shower and x-axis", 80, 0, 360); fAlpha->SetDirectory(NULL); fDist->SetDirectory(NULL); fCosDA->SetDirectory(NULL); fDCA->SetDirectory(NULL); fDCADelta->SetDirectory(NULL); fAlpha->SetXTitle("\\alpha [\\circ]"); fDist->SetXTitle("Dist [mm]"); fCosDA->SetXTitle("cos(\\delta,\\alpha)"); fDCA->SetXTitle("DCA [\\circ]"); fDCADelta->SetXTitle("DCADelta [0, 2\\pi]"); fAlpha->SetYTitle("Counts"); fDist->SetYTitle("Counts"); fCosDA->SetYTitle("Counts"); fDCA->SetYTitle("Counts"); fDCADelta->SetYTitle("Counts"); } // -------------------------------------------------------------------------- // // Delete the four histograms // MHHillasSrc::~MHHillasSrc() { delete fAlpha; delete fDist; delete fCosDA; delete fDCA; delete fDCADelta; } // -------------------------------------------------------------------------- // // Setup the Binning for the histograms automatically if the correct // instances of MBinning (with the names 'BinningAlpha' and 'BinningDist') // are found in the parameter list // Use this function if you want to set the conversion factor which // is used to convert the mm-scale in the camera plain into the deg-scale // used for histogram presentations. The conversion factor is part of // the camera geometry. Please create a corresponding MGeomCam container. // Bool_t MHHillasSrc::SetupFill(const MParList *plist) { const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam"); if (!geom) *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl; else { fMm2Deg = geom->GetConvMm2Deg(); SetMmScale(kFALSE); } ApplyBinning(*plist, "Alpha", fAlpha); ApplyBinning(*plist, "Dist", fDist); ApplyBinning(*plist, "DCA", fDCA); ApplyBinning(*plist, "DCADelta", fDCADelta); return kTRUE; } // -------------------------------------------------------------------------- // // Fill the four histograms with data from a MHillas-Container. // Be careful: Only call this with an object of type MHillas // Bool_t MHHillasSrc::Fill(const MParContainer *par, const Stat_t w) { if (!par) { *fLog << err << "MHHillasSrc::Fill: Pointer (!=NULL) expected." << endl; return kFALSE; } const MHillasSrc &h = *(MHillasSrc*)par; fAlpha->Fill(h.GetAlpha(), w); fDist ->Fill(fUseMmScale ? h.GetDist() : fMm2Deg*h.GetDist(), w); fCosDA->Fill(h.GetCosDeltaAlpha(), w); fDCA ->Fill(fUseMmScale ? h.GetDCA() : fMm2Deg*h.GetDCA(), w); fDCADelta ->Fill(h.GetDCADelta(), w); return kTRUE; } // -------------------------------------------------------------------------- // // Use this function to setup your own conversion factor between degrees // and millimeters. The conversion factor should be the one calculated in // MGeomCam. Use this function with Caution: You could create wrong values // by setting up your own scale factor. // void MHHillasSrc::SetMm2Deg(Float_t mmdeg) { if (mmdeg<=0) { *fLog << warn << dbginf << "Warning - Conversion factor <= 0 - nonsense. Ignored." << endl; return; } if (fMm2Deg>0) *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl; fMm2Deg = mmdeg; } // -------------------------------------------------------------------------- // // With this function you can convert the histogram ('on the fly') between // degrees and millimeters. // void MHHillasSrc::SetMmScale(Bool_t mmscale) { if (fUseMmScale == mmscale) return; if (fMm2Deg<0) { *fLog << warn << GetDescriptor() << ": Warning - Sorry, no conversion factor for conversion available." << endl; return; } const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg; MH::ScaleAxis(fDist, scale); MH::ScaleAxis(fDCA, scale); fDist->SetXTitle(mmscale ? "Dist [mm]" : "Dist [\\circ]"); fDCA ->SetXTitle(mmscale ? "DCA [mm]" : "DCA [\\circ]"); fUseMmScale = mmscale; } // -------------------------------------------------------------------------- // // Creates a new canvas and draws the two histograms into it. // Be careful: The histograms belongs to this object and won't get deleted // together with the canvas. // void MHHillasSrc::Draw(Option_t *o) { TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this); pad->SetBorderMode(0); AppendPad(""); // FIXME: Display Source position // FIXME: If same-option given make two independant y-axis! const TString opt(o); const Bool_t same = opt.Contains("same"); if (!same) pad->Divide(2,2); else { fAlpha->SetLineColor(kGreen); fDist->SetLineColor(kGreen); fDCA->SetLineColor(kGreen); fCosDA->SetLineColor(kGreen); fDCADelta->SetLineColor(kGreen); } pad->cd(1); gPad->SetBorderMode(0); fAlpha->Draw(same?"same":""); pad->cd(2); gPad->SetBorderMode(0); fDist->Draw(same?"same":""); pad->cd(3); gPad->SetBorderMode(0); fDCA->Draw(same?"same":""); pad->cd(4); gPad->SetBorderMode(0); TVirtualPad *p = gPad; if (!same) p->Divide(1,2); p->cd(1); gPad->SetBorderMode(0); fCosDA->Draw(same?"same":""); p->cd(2); gPad->SetBorderMode(0); fDCADelta->Draw(same?"same":""); } void MHHillasSrc::Paint(Option_t *opt) { if (fCosDA->GetEntries()==0) return; TVirtualPad *savepad = gPad; savepad->cd(4); gPad->SetLogy(); gPad = savepad; } TH1 *MHHillasSrc::GetHistByName(const TString name) { if (name.Contains("Alpha", TString::kIgnoreCase)) return fAlpha; if (name.Contains("Dist", TString::kIgnoreCase)) return fDist; if (name.Contains("CosDA", TString::kIgnoreCase)) return fCosDA; return NULL; }