/* ======================================================================== *\ ! ! * ! * 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 ! Author(s): Wolfgang Wittek 2002 ! ! Copyright: MAGIC Software Development, 2000-2004 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MHHillas // // This class contains histograms for the source independent image parameters // ///////////////////////////////////////////////////////////////////////////// #include "MHHillas.h" #include #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MParList.h" #include "MHillas.h" #include "MGeomCam.h" #include "MBinning.h" #include "MHCamera.h" ClassImp(MHHillas); using namespace std; // -------------------------------------------------------------------------- // // Setup four histograms for Width, Length // MHHillas::MHHillas(const char *name, const char *title) : fGeomCam(0), fMm2Deg(1), fUseMmScale(kTRUE) { // // set the name and title of this object // fName = name ? name : "MHHillas"; fTitle = title ? title : "Source independent image parameters"; fLength = new TH1F("Length", "Length of Ellipse", 100, 0, 296.7); fWidth = new TH1F("Width", "Width of Ellipse", 100, 0, 296.7); fDistC = new TH1F("DistC", "Distance from center of camera", 100, 0, 445); fDelta = new TH1F("Delta", "Angle (Main axis - x-axis)", 101, -90, 90); fLength->SetLineColor(kBlue); fLength->SetDirectory(NULL); fWidth->SetDirectory(NULL); fDistC->SetDirectory(NULL); fDelta->SetDirectory(NULL); fLength->SetXTitle("Length [mm]"); fWidth->SetXTitle("Width [mm]"); fDistC->SetXTitle("Distance [mm]"); fDelta->SetXTitle("Delta [\\circ]"); fLength->SetYTitle("Counts"); fWidth->SetYTitle("Counts"); fDistC->SetYTitle("Counts"); fDelta->SetYTitle("Counts"); MBinning bins; bins.SetEdgesLog(50, 1, 1e7); fSize = new TH1F; fSize->SetName("Size"); fSize->SetTitle("Number of Photons"); fSize->SetDirectory(NULL); fSize->SetXTitle("Size"); fSize->SetYTitle("Counts"); fSize->GetXaxis()->SetTitleOffset(1.2); fSize->GetXaxis()->SetLabelOffset(-0.015); fSize->SetFillStyle(4000); fSize->UseCurrentStyle(); bins.Apply(*fSize); fCenter = new TH2F("Center", "Center of Ellipse", 51, -445, 445, 51, -445, 445); fCenter->SetDirectory(NULL); fCenter->SetXTitle("x [mm]"); fCenter->SetYTitle("y [mm]"); fCenter->SetZTitle("Counts"); } // -------------------------------------------------------------------------- // // Delete the histograms // MHHillas::~MHHillas() { delete fLength; delete fWidth; delete fDistC; delete fDelta; delete fSize; delete fCenter; } // -------------------------------------------------------------------------- // // 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) { fGeomCam = (MGeomCam*)plist->FindObject("MGeomCam"); if (!fGeomCam) *fLog << warn << GetDescriptor() << ": No Camera Geometry available. Using mm-scale for histograms." << endl; else { fMm2Deg = fGeomCam->GetConvMm2Deg(); SetMmScale(kFALSE); } ApplyBinning(*plist, "Width", fWidth); ApplyBinning(*plist, "Length", fLength); ApplyBinning(*plist, "Dist", fDistC); ApplyBinning(*plist, "Delta", fDelta); ApplyBinning(*plist, "Size", fSize); const MBinning *bins = (MBinning*)plist->FindObject("BinningCamera"); if (!bins) { float r = fGeomCam ? fGeomCam->GetMaxRadius() : 600; r *= 0.9; if (!fUseMmScale) r *= fMm2Deg; MBinning b; b.SetEdges(61, -r, r); SetBinning(fCenter, &b, &b); } else SetBinning(fCenter, bins, bins); 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; } const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg; MH::ScaleAxis(fLength, scale); MH::ScaleAxis(fWidth, scale); MH::ScaleAxis(fDistC, scale); MH::ScaleAxis(fCenter, scale, scale); if (mmscale) { fLength->SetXTitle("Length [mm]"); fWidth->SetXTitle("Width [mm]"); fDistC->SetXTitle("Distance [mm]"); fCenter->SetXTitle("x [mm]"); fCenter->SetYTitle("y [mm]"); } else { fLength->SetXTitle("Length [\\circ]"); fWidth->SetXTitle("Width [\\circ]"); fDistC->SetXTitle("Distance [\\circ]"); fCenter->SetXTitle("x [\\circ]"); fCenter->SetYTitle("y [\\circ]"); } fUseMmScale = mmscale; } // -------------------------------------------------------------------------- // // Fill the 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 Stat_t w) { if (!par) { *fLog << err << "MHHillas::Fill: Pointer (!=NULL) expected." << endl; return kFALSE; } const MHillas &h = *(MHillas*)par; const Double_t d = sqrt(h.GetMeanX()*h.GetMeanX() + h.GetMeanY()*h.GetMeanY()); const Double_t scale = fUseMmScale ? 1 : fMm2Deg; fLength->Fill(scale*h.GetLength(), w); fWidth ->Fill(scale*h.GetWidth(), w); fDistC ->Fill(scale*d, w); fCenter->Fill(scale*h.GetMeanX(), scale*h.GetMeanY(), w); fDelta ->Fill(kRad2Deg*h.GetDelta(), w); fSize ->Fill(h.GetSize(), w); return kTRUE; } // -------------------------------------------------------------------------- // // Setup a inversed deep blue sea palette for the fCenter histogram. // void MHHillas::SetColors() const { 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); } // -------------------------------------------------------------------------- // // 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 *o) { TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this); pad->SetBorderMode(0); AppendPad(""); TString opt(o); opt.ToLower(); // FIXME: If same-option given make two independant y-axis! const Bool_t same = opt.Contains("same"); if (!same) pad->Divide(2,3); else { fDistC->SetLineColor(kGreen); fSize->SetLineColor(kGreen); fDelta->SetLineColor(kGreen); fWidth->SetLineColor(kMagenta); fLength->SetLineColor(kCyan); } pad->cd(1); gPad->SetBorderMode(0); MH::DrawSame(*fWidth, *fLength, "Width'n'Length", same); pad->cd(2); gPad->SetBorderMode(0); fDistC->Draw(same?"same":""); pad->cd(3); gPad->SetBorderMode(0); gPad->SetLogx(); gPad->SetLogy(); fSize->Draw(same?"same":""); if (!same) { pad->cd(4); gPad->SetBorderMode(0); gPad->SetPad(0.51, 0.01, 0.99, 0.65); SetColors(); fCenter->Draw("colz"); if (fGeomCam) { MHCamera *cam = new MHCamera(*fGeomCam); cam->Draw("same"); cam->SetBit(kCanDelete); } } pad->cd(5); gPad->SetBorderMode(0); fDelta->Draw(same?"same":""); pad->cd(6); if (gPad && !same) delete gPad; pad->Modified(); pad->Update(); } TH1 *MHHillas::GetHistByName(const TString name) { if (name.Contains("Width", TString::kIgnoreCase)) return fWidth; if (name.Contains("Length", TString::kIgnoreCase)) return fLength; if (name.Contains("Size", TString::kIgnoreCase)) return fSize; if (name.Contains("Delta", TString::kIgnoreCase)) return fDelta; if (name.Contains("DistC", TString::kIgnoreCase)) return fDistC; if (name.Contains("Center", TString::kIgnoreCase)) return fCenter; return NULL; } void MHHillas::Paint(Option_t *opt) { SetColors(); MH::Paint(); }