/* ======================================================================== *\ ! ! * ! * 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-2002 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MHHillas // // This class contains histograms for every Hillas parameter // ///////////////////////////////////////////////////////////////////////////// #include "MHHillas.h" #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MGeomCam.h" #include "MHillas.h" #include "MParList.h" ClassImp(MHHillas); // -------------------------------------------------------------------------- // // Setup four histograms for Width, Length // MHHillas::MHHillas(const char *name, const char *title) : fMm2Deg(-1), fUseMmScale(kTRUE) { // // set the name and title of this object // fName = name ? name : "MHHillas"; 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 // fWidth = new TH1F("Width", "Width of Ellipse", 100, 0, 300); fLength = new TH1F("Length", "Length of Ellipse", 100, 0, 300); fLength->SetDirectory(NULL); fWidth->SetDirectory(NULL); fLength->SetXTitle("Length [mm]"); fWidth->SetXTitle("Width [mm]"); fLength->SetYTitle("Counts"); fWidth->SetYTitle("Counts"); } // -------------------------------------------------------------------------- // // Delete the four histograms // MHHillas::~MHHillas() { delete fWidth; delete fLength; } // -------------------------------------------------------------------------- // // Setup the Binning for the histograms automatically if the correct // instances of MBinning (with the names 'BinningWidth' and 'BinningLength') // 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 MHHillas::SetupFill(const MParList *plist) { const MBinning* binsw = (MBinning*)plist->FindObject("BinningWidth"); const MBinning* binsl = (MBinning*)plist->FindObject("BinningLength"); if (!binsw || !binsl) { *fLog << err << dbginf << "At least one MBinning not found... aborting." << endl; return kFALSE; } SetBinning(fWidth, binsw); SetBinning(fLength, binsl); const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam"); if (!geom) { *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl; return kTRUE; } fLength->SetXTitle("Length [\\circ]"); fWidth->SetXTitle("Width [\\circ]"); fMm2Deg = geom->GetConvMm2Deg(); fUseMmScale = kFALSE; 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 MHHillas::Fill(const MParContainer *par) { const MHillas &h = *(MHillas*)par; if (fUseMmScale) { fWidth ->Fill(h.GetWidth()); fLength->Fill(h.GetLength()); } else { fWidth ->Fill(fMm2Deg*h.GetWidth()); fLength->Fill(fMm2Deg*h.GetLength()); } 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 MHHillas::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 MHHillas::SetMmScale(Bool_t mmscale) { if (fUseMmScale == mmscale) return; if (fMm2Deg<0) { *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl; return; } if (fUseMmScale) { fLength->SetXTitle("Length [mm]"); fWidth->SetXTitle("Width [mm]"); fLength->Scale(1./fMm2Deg); fWidth->Scale(1./fMm2Deg); } else { fLength->SetXTitle("Length [\\circ]"); fWidth->SetXTitle("Width [\\circ]"); fLength->Scale(fMm2Deg); fWidth->Scale(fMm2Deg); } fUseMmScale = mmscale; } // -------------------------------------------------------------------------- // // Draw clones of all four histograms. So that the object can be deleted // and the histograms are still visible in the canvas. // The cloned object are deleted together with the canvas if the canvas is // destroyed. If you want to handle dostroying the canvas you can get a // pointer to it from this function // TObject *MHHillas::DrawClone(Option_t *opt) const { TCanvas *c = MakeDefCanvas("Hillas", "Histograms of Hillas Parameters", 350, 500); c->Divide(1, 2); gROOT->SetSelectedPad(NULL); // // This is necessary to get the expected bahviour of DrawClone // c->cd(1); fLength->DrawCopy(); c->cd(2); fWidth->DrawCopy(); c->Modified(); c->Update(); return c; } // -------------------------------------------------------------------------- // // Creates a new canvas and draws the four histograms into it. // Be careful: The histograms belongs to this object and won't get deleted // together with the canvas. // void MHHillas::Draw(Option_t *) { if (!gPad) MakeDefCanvas("Hillas", "Histograms of Hillas Parameters", 350, 500); gPad->Divide(1, 2); gPad->cd(1); fLength->Draw(); gPad->cd(2); fWidth->Draw(); gPad->Modified(); gPad->Update(); }