1 /**************************************************************************
2 * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /** @file AliFMDDisplay.cxx
17 @author Christian Holm Christensen <cholm@nbi.dk>
18 @date Mon Mar 27 12:39:09 2006
19 @brief FMD Event display
21 //___________________________________________________________________
23 // The classes defined here, are utility classes for reading in data
24 // for the FMD. They are put in a seperate library to not polute the
25 // normal libraries. The classes are intended to be used as base
26 // classes for customized class that do some sort of analysis on the
27 // various types of data produced by the FMD.
29 // Latest changes by Christian Holm Christensen
31 #include "AliFMDDisplay.h" // ALIFMDDISPLAY_H
32 #include "AliFMDHit.h" // ALIFMDHIT_H
33 #include "AliFMDDigit.h" // ALIFMDDIGIT_H
34 #include "AliFMDRecPoint.h" // ALIFMDRECPOINT_H
35 #include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
36 #include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
37 // #include <AliESDFMD.h> // ALIESDFMD_H
40 // #include <TArrayF.h>
41 #include <TMarker3DBox.h>
42 #include <TGeoManager.h>
44 #include <TApplication.h>
46 // #include <TParticle.h>
49 #include <TVirtualX.h>
51 //____________________________________________________________________
52 ClassImp(AliFMDDisplay)
54 ; // This is here to keep Emacs for indenting the next line
57 //____________________________________________________________________
58 AliFMDDisplay* AliFMDDisplay::fgInstance = 0;
60 //____________________________________________________________________
62 AliFMDDisplay::Instance()
64 // Return static instance
68 //____________________________________________________________________
69 AliFMDDisplay::AliFMDDisplay(const char* gAliceFile)
70 : AliFMDInput(gAliceFile),
79 // Constructor of an FMD display object.
81 fMarkers = new TObjArray;
82 fHits = new TObjArray;
83 fMarkers->SetOwner(kTRUE);
84 fHits->SetOwner(kFALSE);
88 //____________________________________________________________________
90 AliFMDDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
92 // AliInfo(Form("Event %d, at (%d,%d)", px, py));
93 if (px == 0 && py == 0) return;
94 if (!fZoomMode && fPad->GetView()) {
95 fPad->GetView()->ExecuteRotateView(event, px, py);
98 fPad->SetCursor(kCross);
101 fPad->TAttLine::Modify();
102 fX0 = fPad->AbsPixeltoX(px);
103 fY0 = fPad->AbsPixeltoY(py);
104 fXPixel = fOldXPixel = px;
105 fYPixel = fOldYPixel = py;
110 gVirtualX->DrawBox(fXPixel, fYPixel, fOldXPixel, fOldYPixel,
115 gVirtualX->DrawBox(fXPixel, fYPixel, fOldXPixel, fOldYPixel,
119 fPad->GetCanvas()->FeedbackMode(kFALSE);
120 if (px == fXPixel || py == fYPixel) return;
121 fX1 = fPad->AbsPixeltoX(px);
122 fY1 = fPad->AbsPixeltoY(py);
123 if (fX1 < fX0) std::swap(fX0, fX1);
124 if (fY1 < fY0) std::swap(fY0, fY1);
125 fPad->Range(fX0, fY0, fX1, fY1);
131 //____________________________________________________________________
133 AliFMDDisplay::DistancetoPrimitive(Int_t px, Int_t)
135 // AliInfo(Form("@ (%d,%d)", px, py));
136 fPad->SetCursor(kCross);
137 Float_t xmin = fPad->GetX1();
138 Float_t xmax = fPad->GetX2();
139 Float_t dx = .02 * (xmax - xmin);
140 Float_t x = fPad->AbsPixeltoX(px);
141 if (x < xmin + dx || x > xmax - dx) return 9999;
142 return (fZoomMode ? 0 : 7);
144 //____________________________________________________________________
146 AliFMDDisplay::Init()
148 // Initialize. GEt transforms and such,
149 if (!AliFMDInput::Init()) return kFALSE;
150 AliFMDGeometry* geom = AliFMDGeometry::Instance();
152 geom->InitTransformations();
153 // AliFMDParameters* parm = AliFMDParameters::Instance();
157 //____________________________________________________________________
159 AliFMDDisplay::Begin(Int_t event)
161 // Begin of event. Make canvas is not already done
163 gStyle->SetPalette(1);
164 fCanvas = new TCanvas("display", "Display", 700, 700);
165 fCanvas->SetFillColor(1);
166 fCanvas->ToggleEventStatus();
167 fPad = new TPad("view3D", "3DView", 0.0, 0.05, 1.0, 1.0, 1, 0, 0);
173 fButton = new TButton("Continue", "AliFMDDisplay::Instance()->Continue()",
176 fZoom = new TButton("Zoom", "AliFMDDisplay::Instance()->Zoom()",
179 fPick = new TButton("Pick", "AliFMDDisplay::Instance()->Pick()",
183 AliInfo("Clearing canvas");
186 Warning("End", "No geometry manager");
189 AliInfo("Drawing geometry");
191 fGeoManager->GetTopVolume()->Draw();
192 AliInfo("Adjusting view");
194 if (fPad->GetView()) {
195 fPad->GetView()->SetView(-200, -40, 80, irep);
196 fPad->GetView()->Zoom();
200 return AliFMDInput::Begin(event);
203 //____________________________________________________________________
207 // End of event. Draw everything
214 // fCanvas->Modified(kTRUE);
220 gApplication->StartIdleing();
221 gSystem->InnerLoop();
222 gApplication->StopIdleing();
224 AliInfo("After idle loop");
227 AliInfo("After clearing caches");
228 return AliFMDInput::End();
231 //____________________________________________________________________
233 AliFMDDisplay::LookupColor(Float_t x, Float_t max) const
236 Int_t idx = Int_t(x / max * gStyle->GetNumberOfColors());
237 return gStyle->GetColorPalette(idx);
241 //____________________________________________________________________
243 AliFMDDisplay::ProcessHit(AliFMDHit* hit, TParticle* p)
246 if (!hit) { AliError("No hit"); return kFALSE; }
247 if (!p) { AliError("No track"); return kFALSE; }
251 Float_t pt = TMath::Sqrt(hit->Py()*hit->Py()+hit->Px()*hit->Px());
252 Float_t theta = TMath::ATan2(pt, hit->Pz());
253 Float_t phi = TMath::ATan2(hit->Py(), hit->Px());
254 TMarker3DBox* marker = new TMarker3DBox(hit->X(), hit->Y(), hit->Z(),
255 size, size, size, theta, phi);
256 marker->SetLineColor(LookupColor(hit->Edep(), 1));
257 marker->SetRefObject(hit);
258 fMarkers->Add(marker);
262 //____________________________________________________________________
264 AliFMDDisplay::ProcessDigit(AliFMDDigit* digit)
267 if (!digit) { AliError("No digit"); return kFALSE; }
270 AliFMDGeometry* geom = AliFMDGeometry::Instance();
271 AliFMDParameters* parm = AliFMDParameters::Instance();
272 Double_t threshold = (parm->GetPedestal(digit->Detector(),
276 + 4 * parm->GetPedestalWidth(digit->Detector(),
280 if (digit->Counts() < threshold) return kTRUE;
282 geom->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(),
283 digit->Strip(), x, y, z);
285 Float_t r = TMath::Sqrt(x * x + y * y);
286 Float_t theta = TMath::ATan2(r, z);
287 Float_t phi = TMath::ATan2(y, x);
288 TMarker3DBox* marker = new TMarker3DBox(x,y,z,size,size,size,theta,phi);
289 marker->SetRefObject(digit);
290 marker->SetLineColor(LookupColor(digit->Counts(), 1024));
291 fMarkers->Add(marker);
295 //____________________________________________________________________
297 AliFMDDisplay::ProcessRaw(AliFMDDigit* digit)
300 return ProcessDigit(digit);
303 //____________________________________________________________________
305 AliFMDDisplay::ProcessRecPoint(AliFMDRecPoint* recpoint)
307 // Process reconstructed point
308 if (!recpoint) { AliError("No recpoint"); return kFALSE; }
309 if (recpoint->Particles() < .1) return kTRUE;
310 fHits->Add(recpoint);
312 AliFMDGeometry* geom = AliFMDGeometry::Instance();
313 geom->Detector2XYZ(recpoint->Detector(), recpoint->Ring(),
314 recpoint->Sector(), recpoint->Strip(), x, y, z);
317 Float_t r = TMath::Sqrt(x * x + y * y);
318 Float_t theta = TMath::ATan2(r, z);
319 Float_t phi = TMath::ATan2(y, x);
320 TMarker3DBox* marker = new TMarker3DBox(x,y,z,size,size,size,theta,phi);
321 marker->SetRefObject(recpoint);
322 marker->SetLineColor(LookupColor(recpoint->Particles(), 20));
323 fMarkers->Add(marker);
327 //____________________________________________________________________