/* ======================================================================== *\ ! ! * ! * 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, 05/2002 ! Author(s): Harald Kornmayer, 1/2001 ! ! Copyright: MAGIC Software Development, 2000-2003 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MHCamera // // Camera Display, based on a TH1D. Pleas be carefull using the // underlaying TH1D. // // To change the scale to a logarithmic scale SetLogy() of the Pad. // // Be carefull: Entries in this context means Entries/bin or Events // // FIXME? Maybe MHCamera can take the fLog object from MGeomCam? // //////////////////////////////////////////////////////////////////////////// #include "MHCamera.h" #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MH.h" #include "MHexagon.h" #include "MGeomPix.h" #include "MGeomCam.h" #include "MCamEvent.h" #define kItemsLegend 48 // see SetPalette(1,0) ClassImp(MHCamera); using namespace std; void MHCamera::Init() { UseCurrentStyle(); SetDirectory(NULL); SetLineColor(kGreen); SetMarkerStyle(kFullDotMedium); SetXTitle("Pixel Index"); fNotify = new TList; fNotify->SetBit(kMustCleanup); gROOT->GetListOfCleanups()->Add(fNotify); #if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06) SetPalette(1, 0); #else SetInvDeepBlueSeaPalette(); #endif } // ------------------------------------------------------------------------ // // Default Constructor. To be used by the root system ONLY. // MHCamera::MHCamera() : TH1D(), fGeomCam(NULL), fColors(kItemsLegend) { Init(); } // ------------------------------------------------------------------------ // // Constructor. Makes a clone of MGeomCam. Removed the TH1D from the // current directory. Calls Sumw2(). Set the histogram line color // (for error bars) to Green and the marker style to kFullDotMedium. // MHCamera::MHCamera(const MGeomCam &geom, const char *name, const char *title) : fGeomCam(NULL), fColors(kItemsLegend) { //fGeomCam = (MGeomCam*)geom.Clone(); SetGeometry(geom, name, title); Init(); // // root 3.02 // * base/inc/TObject.h: // register BIT(8) as kNoContextMenu. If an object has this bit set it will // not get an automatic context menu when clicked with the right mouse button. } void MHCamera::SetGeometry(const MGeomCam &geom, const char *name, const char *title) { SetNameTitle(name, title); TAxis &x = *GetXaxis(); SetBins(geom.GetNumPixels(), 0, 1); x.Set(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5); //SetBins(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5); //Rebuild(); Sumw2(); // necessary? if (fGeomCam) delete fGeomCam; fGeomCam = (MGeomCam*)geom.Clone(); fUsed.Set(geom.GetNumPixels()); for (Int_t i=0; i ROOT_VERSION(3,05,00) if (fBuffer) return BufferFill(x,1); #endif const Int_t bin = (Int_t)x+1; AddBinContent(bin); if (fSumw2.fN) fSumw2.fArray[bin]++; if (bin<=0 || bin>fNcells-2) return -1; fTsumw++; fTsumw2++; fTsumwx += x; fTsumwx2 += x*x; return bin; } // ------------------------------------------------------------------------ // // Taken from TH1D::Fill(). Uses the argument directly as bin index. // Doesn't increment the number of entries. // // -*-*-*-*-*-*Increment bin with abscissa X with a weight w*-*-*-*-*-*-*-* // ============================================= // // if x is less than the low-edge of the first bin, the Underflow bin is incremented // if x is greater than the upper edge of last bin, the Overflow bin is incremented // // If the storage of the sum of squares of weights has been triggered, // via the function Sumw2, then the sum of the squares of weights is incremented // by w^2 in the bin corresponding to x. // // -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* Int_t MHCamera::Fill(Axis_t x, Stat_t w) { #if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00) if (fBuffer) return BufferFill(x,w); #endif const Int_t bin = (Int_t)x+1; AddBinContent(bin, w); if (fSumw2.fN) fSumw2.fArray[bin] += w*w; if (bin<=0 || bin>fNcells-2) return -1; const Stat_t z = (w > 0 ? w : -w); fTsumw += z; fTsumw2 += z*z; fTsumwx += z*x; fTsumwx2 += z*x*x; return bin; } // ------------------------------------------------------------------------ // // Use x and y in millimeters // Int_t MHCamera::Fill(Axis_t x, Axis_t y, Stat_t w) { if (fNcells<=1 || IsFreezed()) return -1; for (Int_t idx=0; idx0) continue; SetUsed(idx); return Fill(idx, w); } return -1; } // ------------------------------------------------------------------------ // // Return the mean value of all entries which are used if all=kFALSE and // of all entries if all=kTRUE if sector<0. If sector>=0 only // entries with match the given sector are taken into account. // Stat_t MHCamera::GetMeanSector(Int_t sector, Bool_t all) const { if (fNcells<=1) return 0; Int_t n=0; Stat_t mean = 0; for (int i=0; i=0 only // entries with match the given sector are taken into account. // Stat_t MHCamera::GetRmsSector(Int_t sector, Bool_t all) const { if (fNcells<=1) return -1; Int_t n=0; Stat_t sum = 0; Stat_t sq = 0; for (int i=0; i=0 // only pixels with matching sector number are taken into account. // Double_t MHCamera::GetMinimumSector(Int_t sector, Bool_t all) const { if (fMinimum != -1111) return fMinimum; if (fNcells<=1) return 0; Double_t minimum=FLT_MAX; for (Int_t idx=0; idx=0 // only pixels with matching sector number are taken into account. // Double_t MHCamera::GetMaximumSector(Int_t sector, Bool_t all) const { if (fMaximum!=-1111) return fMaximum; if (fNcells<=1) return 1; Double_t maximum=-FLT_MAX; for (Int_t idx=0; idxmaximum) maximum = fArray[idx+1]; return Profile(maximum); } // ------------------------------------------------------------------------ // // Call this function to draw the camera layout into your canvas. // Setup a drawing canvas. Add this object and all child objects // (hexagons, etc) to the current pad. If no pad exists a new one is // created. (To access the 'real' pad containing the camera you have // to do a cd(1) in the current layer. // // To draw a camera into its own pad do something like: // // MGeomCamMagic m; // MHCamera *d=new MHCamera(m); // // TCanvas *c = new TCanvas; // c->Divide(2,1); // c->cd(1); // // d->FillRandom(); // d->Draw(); // d->SetBit(kCanDelete); // // There are several drawing options: // 'hist' Draw as a standard TH1 histogram (value vs. pixel index) // 'box' Draw hexagons which size is in respect to its contents // 'nocol' Leave the 'boxed' hexagons empty // 'pixelindex' Display the pixel index in each pixel // 'sectorindex' Display the sector index in each pixel // 'content' Display the relative content aligned to GetMaximum() and // GeMinimum() ((val-min)/(max-min)) // 'proj' Display the y-projection of the histogram // void MHCamera::Draw(Option_t *option) { // root 3.02: // gPad->SetFixedAspectRatio() const Color_t col = gPad ? gPad->GetFillColor() : 16; TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas("CamDisplay", "Mars Camera Display", 656, 600); pad->SetBorderMode(0); pad->SetFillColor(col); // // Create an own pad for the MHCamera-Object which can be // resized in paint to keep the correct aspect ratio // pad->Divide(1, 1, 0, 0, col); pad->cd(1); gPad->SetBorderMode(0); TString opt(option); opt.ToLower(); AppendPad(option); // // Do not change gPad. The user should not see, that Draw // changes gPad... // pad->cd(); } // ------------------------------------------------------------------------ // // This is TObject::DrawClone but completely ignores // gROOT->GetSelectedPad(). tbretz had trouble with this in the past. // If this makes trouble please write a bug report. // TObject *MHCamera::DrawClone(Option_t *option) const { // Draw a clone of this object in the current pad //TVirtualPad *pad = gROOT->GetSelectedPad(); TVirtualPad *padsav = gPad; //if (pad) pad->cd(); TObject *newobj = Clone(); if (!newobj) return 0; /* if (pad) { if (strlen(option)) pad->GetListOfPrimitives()->Add(newobj,option); else pad->GetListOfPrimitives()->Add(newobj,GetDrawOption()); pad->Modified(kTRUE); pad->Update(); if (padsav) padsav->cd(); return newobj; } */ TString opt(option); opt.ToLower(); newobj->Draw(opt.IsNull() ? GetDrawOption() : option); if (padsav) padsav->cd(); return newobj; } // ------------------------------------------------------------------------ // // Creates a TH1D which contains the projection of the contents of the // MHCamera onto the y-axis. The maximum and minimum are calculated // such that a slighly wider range than (GetMinimum(), GetMaximum()) is // displayed using THLimitsFinder::OptimizeLimits. // // If no name is given the newly allocated histogram is removed from // the current directory calling SetDirectory(0) in any other case // the newly created histogram is removed from the current directory // and added to gROOT such the gROOT->FindObject can find the histogram. // // If the standard name "_py" is given "_py" is appended to the name // of the MHCamera and the corresponding histogram is searched using // gROOT->FindObject and updated with the present projection. // // It is the responsibility of the user to make sure, that the newly // created histogram is freed correctly. // // Currently the new histogram is restrictred to 50 bins. // Maybe a optimal number can be calulated from the number of // bins on the x-axis of the MHCamera? // // The code was taken mainly from TH2::ProjectX such the interface // is more or less the same than to TH2-projections. // // If sector>=0 only entries with matching sector index are taken // into account. // TH1D *MHCamera::ProjectionS(Int_t sector, const char *name) const { Int_t nbins = 50; // Create the projection histogram TString pname(name); if (name=="_py") { pname.Prepend(GetName()); if (sector>=0) pname += sector; } TH1D *h1=0; //check if histogram with identical name exist TObject *h1obj = gROOT->FindObject(pname); if (h1obj && h1obj->InheritsFrom("TH1D")) { h1 = (TH1D*)h1obj; h1->Reset(); } if (!h1) { Double_t min = GetMinimumSector(sector); Double_t max = GetMaximumSector(sector); Int_t newbins=0; THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE); h1 = new TH1D(pname, GetTitle(), nbins, min, max); h1->SetDirectory(pname.IsNull() ? NULL : gROOT); h1->SetXTitle(GetYaxis()->GetTitle()); h1->SetYTitle("Counts"); //h1->Sumw2(); } // Fill the projected histogram for (Int_t idx=0; idxFill(GetBinContent(idx+1)); return h1; } // ------------------------------------------------------------------------ // // Resizes the current pad so that the camera is displayed in its // correct aspect ratio // void MHCamera::SetRange() { const Float_t range = fGeomCam->GetMaxRadius()*1.05; // // Maintain aspect ratio // const float ratio = TestBit(kNoLegend) ? 1 : 1.15; // // Calculate width and height of the current pad in pixels // Float_t w = gPad->GetWw(); Float_t h = gPad->GetWh()*ratio; // // This prevents the pad from resizing itself wrongly // if (gPad->GetMother() != gPad) { w *= gPad->GetMother()->GetAbsWNDC(); h *= gPad->GetMother()->GetAbsHNDC(); } // // Set Range (coordinate system) of pad // gPad->Range(-range, -range, (2*ratio-1)*range, range); // // Resize Pad to given ratio // if (hSetPad((1.-h/w)/2, 0, (h/w+1.)/2, 1); else gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1.)/2); } // ------------------------------------------------------------------------ // // Updates the pixel colors and paints the pixels // void MHCamera::Update(Bool_t islog, Bool_t isbox, Bool_t iscol) { Double_t min = GetMinimum(kFALSE); Double_t max = GetMaximum(kFALSE); if (min==FLT_MAX) { min = 0; max = 1; } if (min==max) max += 1; UpdateLegend(min, max, islog); MHexagon hex; for (Int_t i=0; i - Pixel Index #" << i << " contents is NaN (Not a Number)..." << endl; hex.SetFillColor(GetColor(GetBinContent(i+1), min, max, islog)); } else hex.SetFillColor(10); MGeomPix &pix = (*fGeomCam)[i]; if (!isbox) hex.PaintHexagon(pix.GetX(), pix.GetY(), pix.GetD()); else if (IsUsed(i) && !TMath::IsNaN(fArray[i+1])) { Float_t size = pix.GetD()*(GetBinContent(i+1)-min)/(max-min); if (size>pix.GetD()) size=pix.GetD(); hex.PaintHexagon(pix.GetX(), pix.GetY(), size); } } } // ------------------------------------------------------------------------ // // Print minimum and maximum // void MHCamera::Print(Option_t *) const { gLog << all << "Minimum: " << GetMinimum(); if (fMinimum==-1111) gLog << " "; gLog << endl; gLog << "Maximum: " << GetMaximum(); if (fMaximum==-1111) gLog << " "; gLog << endl; } // ------------------------------------------------------------------------ // // Paint the y-axis title // void MHCamera::PaintAxisTitle() { const Float_t range = fGeomCam->GetMaxRadius()*1.05; const Float_t w = (1 + 1.5/sqrt((float)(fNcells-2)))*range; TLatex *ptitle = new TLatex(w, -.90*range, GetYaxis()->GetTitle()); ptitle->SetTextSize(0.05); ptitle->SetTextAlign(21); // box with the histogram title ptitle->SetTextColor(gStyle->GetTitleTextColor()); #if ROOT_VERSION_CODE > ROOT_VERSION(3,05,01) ptitle->SetTextFont(gStyle->GetTitleFont("")); #endif ptitle->Paint(); } // ------------------------------------------------------------------------ // // Paints the camera. // void MHCamera::Paint(Option_t *o) { if (fNcells<=1) return; TString opt(o); opt.ToLower(); if (opt.Contains("hist")) { opt.ReplaceAll("hist", ""); TH1D::Paint(opt); return; } if (opt.Contains("proj")) { opt.ReplaceAll("proj", ""); Projection(GetName())->Paint(opt); return; } gPad->Clear(); // Maintain aspect ratio SetRange(); Bool_t isbox = opt.Contains("box"); Bool_t iscol = isbox ? !opt.Contains("nocol") : 1; if (GetPainter()) { // Paint statistics if (!TestBit(TH1::kNoStats)) fPainter->PaintStat(gStyle->GetOptStat(), NULL); // Paint primitives (pixels, color legend, photons, ...) if (fPainter->InheritsFrom(THistPainter::Class())) { static_cast(fPainter)->MakeChopt(""); static_cast(fPainter)->PaintTitle(); } } // Update Contents of the pixels and paint legend Update(gPad->GetLogy(), isbox, iscol); PaintAxisTitle(); if (opt.Contains("pixelindex")) PaintIndices(0); if (opt.Contains("sectorindex")) PaintIndices(1); if (opt.Contains("content")) PaintIndices(2); } // ------------------------------------------------------------------------ // // With this function you can change the color palette. For more // information see TStyle::SetPalette. Only palettes with 50 colors // are allowed. // In addition you can use SetPalette(52, 0) to create an inverse // deep blue sea palette // void MHCamera::SetPalette(Int_t ncolors, Int_t *colors) { // // If not enough colors are specified skip this. // if (ncolors>1 && ncolors<50) { gLog << err << "MHCamera::SetPalette: Only default palettes with 50 colors are allowed... ignored." << endl; return; } // // If ncolors==52 create a reversed deep blue sea palette // if (ncolors==52) { gStyle->SetPalette(51, NULL); TArrayI c(kItemsLegend); for (int i=0; iGetColorPalette(i); gStyle->SetPalette(kItemsLegend, c.GetArray()); } else gStyle->SetPalette(ncolors, colors); fColors.Set(kItemsLegend); for (int i=0; iGetColorPalette(i); } void MHCamera::SetPrettyPalette() { if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) SetPalette(1, 0); } void MHCamera::SetDeepBlueSeaPalette() { if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) SetPalette(51, 0); } void MHCamera::SetInvDeepBlueSeaPalette() { if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) SetPalette(52, 0); } void MHCamera::PaintIndices(Int_t type) { if (fNcells<=1) return; const Double_t min = GetMinimum(); const Double_t max = GetMaximum(); if (type==2 && max==min) return; TText txt; txt.SetTextFont(122); txt.SetTextAlign(22); // centered/centered for (Int_t i=0; iGetMaxRadius()/1.05); txt.PaintText(h.GetX(), h.GetY(), num); } } // ------------------------------------------------------------------------ // // Call this function to add a MCamEvent on top of the present contents. // void MHCamera::AddCamContent(const MCamEvent &event, Int_t type) { if (fNcells<=1 || IsFreezed()) return; // FIXME: Security check missing! for (Int_t idx=0; idx0) // error on the mean { const Double_t error = fSumw2.fArray[bin]/GetEntries(); const Double_t val = fArray[bin]/GetEntries(); rc = TMath::Sqrt(error - val*val); } else { rc = TH1D::GetBinError(bin); rc = Profile(rc); } return rc; } // ------------------------------------------------------------------------ // // Call this function to add a MHCamera on top of the present contents. // Type: // 0) bin content // 1) errors // 2) rel. errors // void MHCamera::AddCamContent(const MHCamera &d, Int_t type) { if (fNcells!=d.fNcells || IsFreezed()) return; // FIXME: Security check missing! for (Int_t idx=0; idxGetSize()!=fNcells-2) return; for (Int_t idx=0; idx(event)[idx]); // FIXME: Slow! if (!used || (*const_cast(used))[idx]) SetUsed(idx); } fEntries++; } // ------------------------------------------------------------------------ // // Call this function to add a MCamEvent on top of the present contents. // 1 is added to each pixel if the contents of MCamEvent>threshold // void MHCamera::CntCamContent(const MCamEvent &event, Double_t threshold, Int_t type) { if (fNcells<=1 || IsFreezed()) return; // FIXME: Security check missing! for (Int_t idx=0; idxthreshold) Fill(idx); } fEntries++; } // ------------------------------------------------------------------------ // // Call this function to add a TArrayD on top of the present contents. // 1 is added to each pixel if the contents of MCamEvent>threshold // void MHCamera::CntCamContent(const TArrayD &event, Double_t threshold, Bool_t ispos) { if (event.GetSize()!=fNcells-2 || IsFreezed()) return; for (Int_t idx=0; idx(event)[idx]>threshold) Fill(idx); if (!ispos || fArray[idx+1]>0) SetUsed(idx); } fEntries++; } // ------------------------------------------------------------------------ // // Fill the pixels with random contents. // void MHCamera::FillRandom() { if (fNcells<=1 || IsFreezed()) return; Reset(); // FIXME: Security check missing! for (Int_t idx=0; idxUniform()*fGeomCam->GetPixRatio(idx)); SetUsed(idx); } fEntries=1; } // ------------------------------------------------------------------------ // // The array must be in increasing order, eg: 2.5, 3.7, 4.9 // The values in each bin are replaced by the interval in which the value // fits. In the example we have four intervals // (<2.5, 2.5-3.7, 3.7-4.9, >4.9). Maximum and minimum are set // accordingly. // void MHCamera::SetLevels(const TArrayF &arr) { if (fNcells<=1) return; for (Int_t i=0; iSetPalette(1,0) // for the display. So we have to convert the value "wert" into // a color index that fits the color palette. // The range of the color palette is defined by the values fMinPhe // and fMaxRange. Between this values we have 50 color index, starting // with 0 up to 49. // Int_t MHCamera::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog) { if (TMath::IsNaN(val)) // FIXME: gLog! return 10; // // first treat the over- and under-flows // const Int_t maxcolidx = kItemsLegend-1; if (val >= max) return fColors[maxcolidx]; if (val <= min) return fColors[0]; // // calculate the color index // Float_t ratio; if (islog && min>0) ratio = log10(val/min) / log10(max/min); else ratio = (val-min) / (max-min); const Int_t colidx = (Int_t)(ratio*maxcolidx + .5); return fColors[colidx]; } TPaveStats *MHCamera::GetStatisticBox() { TObject *obj = 0; TIter Next(fFunctions); while ((obj = Next())) if (obj->InheritsFrom(TPaveStats::Class())) return static_cast(obj); return NULL; } // ------------------------------------------------------------------------ // // Change the text on the legend according to the range of the Display // void MHCamera::UpdateLegend(Float_t min, Float_t max, Bool_t islog) { const Float_t range = fGeomCam->GetMaxRadius()*1.05; TArrow arr; arr.PaintArrow(-range*.9, -range*.9, -range*.6, -range*.9, 0.025); arr.PaintArrow(-range*.9, -range*.9, -range*.9, -range*.6, 0.025); TString text; text += (int)(range*.3); text += "mm"; TText newtxt2; newtxt2.SetTextSize(0.04); newtxt2.PaintText(-range*.85, -range*.85, text); text = ""; text += (float)((int)(range*.3*fGeomCam->GetConvMm2Deg()*10))/10; text += "\\circ"; text = text.Strip(TString::kLeading); TLatex latex; latex.PaintLatex(-range*.85, -range*.75, 0, 0.04, text); if (TestBit(kNoLegend)) return; TPaveStats *stats = GetStatisticBox(); const Float_t hndc = 0.92 - (stats ? stats->GetY1NDC() : 1); const Float_t H = (0.75-hndc)*range; const Float_t offset = hndc*range; const Float_t h = 2./kItemsLegend; const Float_t w = range/sqrt((float)(fNcells-2)); TBox newbox; TText newtxt; newtxt.SetTextSize(0.03); newtxt.SetTextAlign(12); #if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06) newtxt.SetBit(/*kNoContextMenu|*/kCannotPick); newbox.SetBit(/*kNoContextMenu|*/kCannotPick); #endif const Float_t step = (islog && min>0 ? log10(max/min) : max-min) / kItemsLegend; const Int_t firsts = step*3 < 1e-8 ? 8 : (Int_t)floor(log10(step*3)); const TString opt = Form("%%.%if", firsts>0 ? 0 : TMath::Abs(firsts)); for (Int_t i=0; i0) val = pow(10, step*i) * min; else val = min + step*i; //const bool dispexp = max-min>1.5 && fabs(val)>0.1 && fabs(val)<1e6; newtxt.PaintText(range+1.5*w, H*(i*h-1)-offset, Form(opt, val)); } for (Int_t i=0; iClassSaved(TCanvas::Class())) fDrawingPad->SavePrimitive(out, opt); out << " " << fDrawingPad->GetName() << "->SetWindowSize("; out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl; */ } // ------------------------------------------------------------------------ // // compute the distance of a point (px,py) to the Camera // this functions needed for graphical primitives, that // means without this function you are not able to interact // with the graphical primitive with the mouse!!! // // All calcutations are done in pixel coordinates // Int_t MHCamera::DistancetoPrimitive(Int_t px, Int_t py) { if (fNcells<=1) return 999999; const Int_t kMaxDiff = 7; TPaveStats *box = (TPaveStats*)gPad->GetPrimitive("stats"); if (box) { const Double_t w = box->GetY2NDC()-box->GetY1NDC(); box->SetX1NDC(gStyle->GetStatX()-gStyle->GetStatW()); box->SetY1NDC(gStyle->GetStatY()-w); box->SetX2NDC(gStyle->GetStatX()); box->SetY2NDC(gStyle->GetStatY()); } if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) return TH1D::DistancetoPrimitive(px, py); for (Int_t i=0; iDistancetoPrimitive(px, py); if (dist > kMaxDiff) return 999999; gPad->SetSelected(box); return dist; } // ------------------------------------------------------------------------ // // Int_t MHCamera::GetPixelIndex(Int_t px, Int_t py) const { if (fNcells<=1) return -1; Int_t i; for (i=0; i0) continue; return i; } return -1; } // ------------------------------------------------------------------------ // // Returns string containing info about the object at position (px,py). // Returned string will be re-used (lock in MT environment). // char *MHCamera::GetObjectInfo(Int_t px, Int_t py) const { if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) return TH1D::GetObjectInfo(px, py); static char info[128]; const Int_t idx=GetPixelIndex(px, py); if (idx<0) return TObject::GetObjectInfo(px, py); sprintf(info, "Software Pixel Index: %d (Hardware Id=%d)", idx, idx+1); return info; } // ------------------------------------------------------------------------ // // Add a MCamEvent which should be displayed when the user clicks on a // pixel. // Warning: The object MUST inherit from TObject AND MCamEvent // void MHCamera::AddNotify(TObject *obj) { // Make sure, that the object derives from MCamEvent! MCamEvent *evt = dynamic_cast(obj); if (!evt) { gLog << err << "ERROR: MHCamera::AddNotify - TObject doesn't inherit from MCamEvent... ignored." << endl; return; } // Make sure, that it is deleted from the list too, if the obj is deleted obj->SetBit(kMustCleanup); // Add object to list fNotify->Add(obj); } // ------------------------------------------------------------------------ // // Execute a mouse event on the camera // void MHCamera::ExecuteEvent(Int_t event, Int_t px, Int_t py) { if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) { TH1D::ExecuteEvent(event, px, py); return; } //if (event==kMouseMotion && fStatusBar) // fStatusBar->SetText(GetObjectInfo(px, py), 0); if (event!=kButton1Down) return; const Int_t idx = GetPixelIndex(px, py); if (idx<0) return; gLog << all << GetTitle() << " <" << GetName() << ">" << endl; gLog << "Software Pixel Idx: " << idx << endl; gLog << "Hardware Pixel Id: " << idx+1 << endl; gLog << "Contents: " << GetBinContent(idx+1); if (GetBinError(idx+1)>0) gLog << " +/- " << GetBinError(idx+1); gLog << " <" << (IsUsed(idx)?"on":"off") << ">" << endl; if (fNotify && fNotify->GetSize()>0) { // FIXME: Is there a simpler and more convinient way? // The name which is created here depends on the instance of // MHCamera and on the pad on which it is drawn --> The name // is unique. For ExecuteEvent gPad is always correctly set. const TString name = Form("%p;%p;PixelContent", this, gPad); TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name); if (old) old->cd(); else new TCanvas(name); /* TIter Next(gPad->GetListOfPrimitives()); TObject *o; while (o=Next()) cout << o << ": " << o->GetName() << " " << o->IsA()->GetName() << endl; */ // FIXME: Make sure, that the old histograms are really deleted. // Are they already deleted? // The dynamic_cast is necessary here: We cannot use ForEach TIter Next(fNotify); MCamEvent *evt; while ((evt=dynamic_cast(Next()))) evt->DrawPixelContent(idx); gPad->Modified(); gPad->Update(); } } UInt_t MHCamera::GetNumPixels() const { return fGeomCam->GetNumPixels(); } TH1 *MHCamera::DrawCopy() const { gPad=NULL; return TH1D::DrawCopy(fName+";cpy"); }