/* ======================================================================== *\ ! ! * ! * 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 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MCamDisplay // // Camera Display. The Pixels are displayed in // contents/area [somthing/mm^2] // // To change the scale to a logarithmic scale SetLogz() of the Pad. // // //////////////////////////////////////////////////////////////////////////// #include "MCamDisplay.h" #include #include #include #include #include #include #include #include #include #include #include #include "MH.h" #include "MHexagon.h" #include "MGeomCam.h" #include "MRflEvtData.h" #include "MCerPhotPix.h" #include "MCerPhotEvt.h" #include "MPedestalPix.h" #include "MPedestalCam.h" #include "MCurrents.h" #include "MCamEvent.h" #include "MImgCleanStd.h" #define kItemsLegend 48 // see SetPalette(1,0) ClassImp(MCamDisplay); using namespace std; // ------------------------------------------------------------------------ // // Default Constructor. To be used by the root system ONLY. // MCamDisplay::MCamDisplay() : fGeomCam(NULL), fAutoScale(kTRUE), fColors(kItemsLegend) { fNumPixels = 0; fRange = 0; fPixels = NULL; fLegend = NULL; fLegText = NULL; fPhotons = NULL; fArrowX = NULL; fArrowY = NULL; fLegRadius = NULL; fLegDegree = NULL; fMinimum = 0; fMaximum = 1; fNotify = NULL; #if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06) SetPalette(1, 0); #else SetPalette(51, 0); #endif } // ------------------------------------------------------------------------ // // Constructor. Makes a clone of MGeomCam. // MCamDisplay::MCamDisplay(MGeomCam *geom) : fGeomCam(NULL), fAutoScale(kTRUE), fColors(kItemsLegend), fData(geom->GetNumPixels()), fMinimum(0), fMaximum(1) { fGeomCam = (MGeomCam*)geom->Clone(); fNotify = new TList; // // create the hexagons of the display // fNumPixels = fGeomCam->GetNumPixels(); fRange = fGeomCam->GetMaxRadius(); // 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. fPhotons = new TClonesArray("TMarker", 0); // // Construct all hexagons. Use new-operator with placement // fPixels = new TClonesArray("MHexagon", fNumPixels); for (UInt_t i=0; i ROOT_VERSION(3,01,06) pix.SetBit(/*kNoContextMenu|*/kCannotPick); #endif pix.SetFillColor(16); pix.ResetBit(kIsUsed); } // // set up the Legend // const Float_t H = 0.9*fRange; const Float_t h = 2./kItemsLegend; const Float_t w = fRange/sqrt((float)fNumPixels); fLegend = new TClonesArray("TBox", kItemsLegend); fLegText = new TClonesArray("TText", kItemsLegend+1); for (Int_t i=0; i ROOT_VERSION(3,01,06) newbox.SetBit(/*kNoContextMenu|*/kCannotPick); #endif } for (Int_t i=0; i ROOT_VERSION(3,01,06) newtxt.SetBit(/*kNoContextMenu|*/kCannotPick); #endif } fArrowX = new TArrow(-fRange*.9, -fRange*.9, -fRange*.6, -fRange*.9, 0.025); fArrowY = new TArrow(-fRange*.9, -fRange*.9, -fRange*.9, -fRange*.6, 0.025); TString text; text += (int)(fRange*.3); text += "mm"; fLegRadius = new TText(-fRange*.85, -fRange*.85, text); text = ""; text += (float)((int)(fRange*.3*fGeomCam->GetConvMm2Deg()*10))/10; text += "\\circ"; text = text.Strip(TString::kLeading); fLegDegree = new TLatex(-fRange*.85, -fRange*.75, text); fLegRadius->SetTextSize(0.04); fLegDegree->SetTextSize(0.04); #if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06) SetPalette(1, 0); #else SetPalette(51, 0); #endif } // ------------------------------------------------------------------------ // // Destructor. Deletes TClonesArrays for hexagons and legend elements. // MCamDisplay::~MCamDisplay() { if (fPixels) { fPixels->Delete(); delete fPixels; } if (fLegend) { fLegend->Delete(); delete fLegend; } if (fLegText) { fLegText->Delete(); delete fLegText; } if (fPhotons) { fPhotons->Delete(); delete fPhotons; } if (fArrowX) delete fArrowX; if (fArrowY) delete fArrowY; if (fLegRadius) delete fLegRadius; if (fLegDegree) delete fLegDegree; if (fGeomCam) delete fGeomCam; if (fNotify) delete fNotify; } // ------------------------------------------------------------------------ // // 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 draw a camera into its own pad do something like: // // TCanvas *c = new TCanvas; // c->Divide(2,1); // MGeomCamMagic m; // MCamDisplay *d=new MCamDisplay(&m); // d->FillRandom(); // c->cd(1); // gPad->SetBorderMode(0); // gPad->Divide(1,1); // gPad->cd(1); // d->Draw(); // d->SetBit(kCanDelete); // void MCamDisplay::Draw(Option_t *option) { // root 3.02: // gPad->SetFixedAspectRatio() TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas("CamDisplay", "Mars Camera Display", 656, 600); pad->SetBorderMode(0); pad->SetFillColor(16); AppendPad(""); } void MCamDisplay::SetRange() { // // Maintain aspect ratio // const float ratio = 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(-fRange, -fRange, (2*ratio-1)*fRange, fRange); // // 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); } void MCamDisplay::Update(Bool_t islog) { // FIXME: Don't do this if nothing changed! if (fAutoScale) { fMinimum = FLT_MAX; fMaximum = -FLT_MAX; for (UInt_t i=0; ifMaximum) fMaximum = fData[i]; } } if (fMinimum==FLT_MAX && fMaximum==-FLT_MAX) fMinimum = fMaximum = 0; if (fMinimum==fMaximum) fMaximum = fMinimum + 1; UpdateLegend(fMinimum, fMaximum, islog); for (UInt_t i=0; iGetLogz()); // Paint Legend fArrowX->Paint(">"); fArrowY->Paint(">"); fLegRadius->Paint(); fLegDegree->Paint(); // Paint primitives (pixels, color legend, photons, ...) { fPixels->ForEach( TObject, Paint)(); } { fLegend->ForEach( TObject, Paint)(); } { fLegText->ForEach(TObject, Paint)(); } { fPhotons->ForEach(TObject, Paint)(); } } // ------------------------------------------------------------------------ // // 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 MCamDisplay::SetPalette(Int_t ncolors, Int_t *colors) { // // If not enough colors are specified skip this. // if (ncolors>1 && ncolors<50) { cout << "MCamDisplay::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); if (!fPixels) { for (int i=0; iGetColorPalette(i); return; } // // Change the colors of the pixels // for (unsigned int i=0; iGetColorPalette(idx)); } // // Store the color palette used for a later reverse lookup // for (int i=0; iGetColorPalette(i); GetBox(i)->SetFillColor(fColors[i]); } } void MCamDisplay::SetPrettyPalette() { SetPalette(1, 0); } void MCamDisplay::SetDeepBlueSeaPalette() { SetPalette(51, 0); } void MCamDisplay::SetInvDeepBlueSeaPalette() { SetPalette(52, 0); } void MCamDisplay::SetPalette() { for (int i=0; iSetFillColor(fColors[i]); } void MCamDisplay::DrawPixelNumbers() { for (int i=0; iSetTextSize(0.3*h.GetD()/fGeomCam->GetMaxRadius()); } } // ------------------------------------------------------------------------ // // Call this function to fill the currents // void MCamDisplay::Fill(const MCamEvent &event, Int_t type) { Reset(); // FIXME: Security check missing! for (UInt_t idx=0; idxGetPixRatio(idx), type)) (*this)[idx].SetBit(kIsUsed); } } void MCamDisplay::FillRandom() { Reset(); // FIXME: Security check missing! for (UInt_t idx=0; idxUniform(); (*this)[idx].SetBit(kIsUsed); } } // ------------------------------------------------------------------------ // // Call this function to fill the currents // void MCamDisplay::Fill(const TArrayF &event, Bool_t ispos) { Reset(); fData = event; for (UInt_t idx=0; idx0) (*this)[idx].SetBit(kIsUsed); } // ------------------------------------------------------------------------ // // Fill the colors in respect to the cleaning levels // void MCamDisplay::FillLevels(const MCerPhotEvt &event, Float_t lvl1, Float_t lvl2) { Fill(event, 2); for (UInt_t i=0; ilvl1) fData[i] = 0; else if (fData[i]>lvl2) fData[i] = 1; else fData[i] = 2; } } // ------------------------------------------------------------------------ // // Fill the colors in respect to the cleaning levels // void MCamDisplay::FillLevels(const MCerPhotEvt &event, const MImgCleanStd &clean) { FillLevels(event, clean.GetCleanLvl1(), clean.GetCleanLvl2()); } // ------------------------------------------------------------------------ // // Show a reflector event. EMarkerStyle is defined in root/include/Gtypes.h // To remove the photons from the display call FillRflEvent(NULL) // void MCamDisplay::ShowRflEvent(const MRflEvtData *event, EMarkerStyle ms) { const Int_t num = event ? event->GetNumPhotons() : 0; fPhotons->ExpandCreate(num); if (num < 1) return; Int_t i=num-1; do { const MRflSinglePhoton &ph = event->GetPhoton(i); TMarker &m = (TMarker&)*fPhotons->UncheckedAt(i); m.SetX(ph.GetX()); m.SetY(ph.GetY()); m.SetMarkerStyle(ms); } while (i--); } // ------------------------------------------------------------------------ // // Fill a reflector event. Sums all pixels in each pixel as the // pixel contents. // // WARNING: Due to the estimation in DistanceToPrimitive and the // calculation in pixels instead of x, y this is only a // rough estimate. // void MCamDisplay::FillRflEvent(const MRflEvtData &event) { // // reset pixel colors to background color // Reset(); // // sum the photons content in each pixel // const Int_t entries = event.GetNumPhotons(); for (Int_t i=0; iGetPixRatio(id); } // // Set color of pixels // for (UInt_t id=0; id0) (*this)[id].SetBit(kIsUsed); } // ------------------------------------------------------------------------ // // Reset the all pixel colors to a default value // void MCamDisplay::Reset() { for (UInt_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 MCamDisplay::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog) { // // 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]; } // ------------------------------------------------------------------------ // // Change the text on the legend according to the range of the Display // void MCamDisplay::UpdateLegend(Float_t minphe, Float_t maxphe, Bool_t islog) { for (Int_t i=0; i0) val = pow(10, log10(maxphe/minphe)*pos) * minphe; else val = minphe + pos * (maxphe-minphe); TText &txt = *(TText*)fLegText->At(i); txt.SetText(txt.GetX(), txt.GetY(), Form(val<1e6?"%5.1f":"%5.1e", val)); } } // ------------------------------------------------------------------------ // // Save primitive as a C++ statement(s) on output stream out // void MCamDisplay::SavePrimitive(ofstream &out, Option_t *opt) { cout << "MCamDisplay::SavePrimitive: Must be rewritten!" << endl; /* if (!gROOT->ClassSaved(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 running in pixel coordinates // Int_t MCamDisplay::DistancetoPrimitive(Int_t px, Int_t py) { Int_t dist = 999999; for (UInt_t i=0; iDistancetoPrimitive(px, py); if (dDistancetoPrimitive(px, py)>0) 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 *MCamDisplay::GetObjectInfo(Int_t px, Int_t py) const { 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; } // ------------------------------------------------------------------------ // // Execute a mouse event on the camera // void MCamDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py) { //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; cout << "Software Pixel Index: " << idx << endl; cout << "Hardware Pixel Id: " << idx+1 << endl; cout << "Contents: " << fData[idx] << endl; //fNotify->Print(); if (fNotify->GetSize()>0) new TCanvas; fNotify->ForEach(MCamEvent, DrawPixelContent)(idx); }