]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveV0Editor.cxx
-Fixed bug in offline buffer that deleted event object.
[u/mrichter/AliRoot.git] / EVE / EveBase / AliEveV0Editor.cxx
1 // @(#)root/eve:$Id$
2 // Author: Matevz Tadel 2007
3
4 /**************************************************************************
5  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
7  * full copyright notice.                                                 *
8  **************************************************************************/
9
10 #include "AliEveV0Editor.h"
11 #include "AliEveV0.h"
12
13 #include "TVirtualPad.h"
14 #include "TColor.h"
15
16 // Cleanup these includes:
17 #include "TGLabel.h"
18 #include "TGButton.h"
19 #include "TGNumberEntry.h"
20 #include "TGColorSelect.h"
21 #include "TGDoubleSlider.h"
22
23
24 //______________________________________________________________________________
25 // GUI editor for AliEveV0.
26 //
27
28 ClassImp(AliEveV0Editor)
29
30 //______________________________________________________________________________
31 AliEveV0Editor::AliEveV0Editor(const TGWindow *p, Int_t width, Int_t height,
32                                UInt_t options, Pixel_t back) :
33   TGedFrame(p, width, height, options | kVerticalFrame, back),
34   fM(0),
35   fInfoLabel0(0),
36   fInfoLabel1(0),
37   fInfoLabelNegDaughter(0),
38   fInfoLabelPosDaughter(0),
39   fXButton(0)
40   // Initialize widget pointers to 0
41 {
42   // Constructor.
43
44   MakeTitle("AliEveV0");
45
46   fInfoLabel0 = new TGLabel(this);
47   fInfoLabel0->SetTextJustify(kTextLeft);
48   AddFrame(fInfoLabel0, new TGLayoutHints(kLHintsLeft|kLHintsExpandX,
49                                           8, 0, 2, 0));
50
51   fInfoLabel1 = new TGLabel(this);
52   fInfoLabel1->SetTextJustify(kTextLeft);
53   AddFrame(fInfoLabel1, new TGLayoutHints(kLHintsLeft|kLHintsExpandX,
54                                           8, 0, 2, 0));
55
56   fInfoLabelNegDaughter = new TGLabel(this);
57   fInfoLabelNegDaughter->SetTextJustify(kTextLeft);
58   AddFrame(fInfoLabelNegDaughter, new TGLayoutHints(kLHintsLeft|kLHintsExpandX,
59                                                     8, 0, 2, 0));
60
61   fInfoLabelPosDaughter = new TGLabel(this);
62   fInfoLabelPosDaughter->SetTextJustify(kTextLeft);
63   AddFrame(fInfoLabelPosDaughter, new TGLayoutHints(kLHintsLeft|kLHintsExpandX,
64                                                     8, 0, 2, 0));
65
66   fXButton = new TGTextButton(this, "Detailed View");
67   AddFrame(fXButton, new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 1, 1, 0, 0));
68   fXButton->Connect("Clicked()", "AliEveV0Editor", this, "DisplayDetailed()");
69 }
70
71 /******************************************************************************/
72
73 //______________________________________________________________________________
74 void AliEveV0Editor::SetModel(TObject* obj)
75 {
76   // Set model object.
77
78   fM = dynamic_cast<AliEveV0*>(obj);
79
80   // Set values of widgets
81   fInfoLabel0->SetText(Form("Radius = %f, DCA = %f", fM->GetRadius(), fM->GetDaughterDCA()));
82   fInfoLabel1->SetText(Form("Pt = %f", fM->GetPt()));
83   fInfoLabelNegDaughter->SetText(Form("Neg. Daughter Prob= %.2f for Pdg= %d", fM->GetNegMaxProbPid(), fM->GetNegMaxProbPdg()));
84   fInfoLabelPosDaughter->SetText(Form("Pos. Daughter Prob= %.2f for Pdg= %d", fM->GetPosMaxProbPid(), fM->GetPosMaxProbPdg()));
85
86 }
87
88 /******************************************************************************/
89
90 // Implements callback/slot methods
91
92 //______________________________________________________________________________
93 // void AliEveV0Editor::DoXYZZ()
94 // {
95 //    // Slot for XYZZ.
96 //
97 //    fM->SetXYZZ(fXYZZ->GetValue());
98 //    Update();
99 // }
100
101 #include <TEveManager.h>
102 #include <TEveWindow.h>
103 #include <TEveViewer.h>
104 #include <TEveScene.h>
105 #include <TEveGeoNode.h>
106 #include <TEveProjectionManager.h>
107
108 #include <TGLCamera.h>
109 #include <TGLViewer.h>
110 #include "TGLCameraOverlay.h"
111
112 #include <TLatex.h>
113 #include <TRootEmbeddedCanvas.h>
114 #include <TInterpreter.h>
115
116 void AliEveV0Editor::DisplayDetailed()
117 {
118   TEveWindowSlot *slot = TEveWindow::CreateWindowMainFrame();
119   TEveWindowPack *pack = slot->MakePack();
120   pack->SetShowTitleBar(kFALSE);
121   pack->SetHorizontal();
122
123   //
124   // This part is for getting the different objects to display
125   //
126   char displayInfo[100] = {0};
127   sprintf(displayInfo,"pt = %.3f",fM->GetPt());
128   TEveLine *lv0TransverseMomentumDirection = new TEveLine(displayInfo);
129   lv0TransverseMomentumDirection->SetLineColor(kOrange+8);
130   lv0TransverseMomentumDirection->SetLineWidth(2);
131   lv0TransverseMomentumDirection->SetLineStyle(2);
132   lv0TransverseMomentumDirection->SetLineWidth(2);
133   Float_t scalePt = 100.; // this needs to be available as a ruler
134   lv0TransverseMomentumDirection->SetPoint(0,fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
135   lv0TransverseMomentumDirection->SetPoint(1,scalePt*fM->fRecDecayP.fX, scalePt*fM->fRecDecayP.fY,0);
136
137   TEvePointSet *pvlocation = new TEvePointSet("PV location");
138   pvlocation->SetNextPoint(fM->fRecBirthV.fX, fM->fRecBirthV.fY, fM->fRecBirthV.fZ);
139   pvlocation->SetTitle("pv location");
140   pvlocation->SetMarkerStyle(4);
141   pvlocation->SetMarkerSize(2.5);
142   pvlocation->SetMarkerColor(7);
143
144   TEvePointSet *v0location = new TEvePointSet("V0 location");
145   v0location->SetNextPoint(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
146   v0location->SetTitle("v0 location");
147   v0location->SetMarkerStyle(4);
148   v0location->SetMarkerSize(2.5);
149   v0location->SetMarkerColor(kOrange+8);
150
151   TEveUtil::LoadMacro("clusters_from_index.C");
152
153   AliEveTrack *negTrack = fM->GetNegTrack();
154   AliEveTrack *posTrack = fM->GetPosTrack();
155
156
157   char macroWithIndex[100] = {0};
158   Int_t daughterIndex = 0;
159   TEvePointSet *negDaughterCluster = 0;
160   TEvePointSet *posDaughterCluster = 0;
161   
162   daughterIndex = negTrack->GetIndex();
163   sprintf(macroWithIndex,"clusters_from_index(%d)",daughterIndex);
164   Long_t negResult = gInterpreter->ProcessLine(macroWithIndex);
165   if (negResult) {
166     negDaughterCluster = reinterpret_cast<TEvePointSet*>(negResult);
167     if (negDaughterCluster){
168       negDaughterCluster->SetMarkerStyle(4);    
169       negDaughterCluster->SetMarkerSize(1.5);   
170       negDaughterCluster->SetMarkerColor(kBlue+3);
171     }
172   }
173   else
174   {
175     Warning("DisplayDetailed", "Import of negative daughter's clusters failed.");
176   }
177
178   daughterIndex = posTrack->GetIndex();
179   sprintf(macroWithIndex,"clusters_from_index(%d)",daughterIndex);
180   Long_t posResult = gInterpreter->ProcessLine(macroWithIndex);
181   if (posResult) {
182     posDaughterCluster = reinterpret_cast<TEvePointSet*>(posResult);
183     if (posDaughterCluster){
184       posDaughterCluster->SetMarkerStyle(4);    
185       posDaughterCluster->SetMarkerSize(1.5);   
186       posDaughterCluster->SetMarkerColor(kRed+3);
187     }
188   }
189   else
190   {
191     Warning("DisplayDetailed", "Import of positive daughter's clusters failed.");
192   }
193
194   //
195   // This part is for the bending plane view
196   //
197   pack->NewSlot()->MakeCurrent();
198   TEveViewer *bpViewer = gEve->SpawnNewViewer("V0 bending plane View");
199   TEveScene  *bpScene  = gEve->SpawnNewScene("V0 bending plane Scene");
200   
201   TEveProjectionManager *projMgr = new TEveProjectionManager(TEveProjection::kPT_RPhi);
202   bpScene->AddElement(projMgr);
203   bpViewer->AddScene(bpScene);
204
205   TEveUtil::LoadMacro("geom_gentle.C");
206   Long_t result = gInterpreter->ProcessLine("geom_gentle_rphi()");
207   if (result)
208   {
209     TEveGeoShape *geomRPhi = reinterpret_cast<TEveGeoShape*>(result);
210     geomRPhi->IncDenyDestroy();
211     projMgr->SetCurrentDepth(-10); // to put the geometry behind the projection of the V0 -> clearer
212     projMgr->ImportElements(geomRPhi);
213     projMgr->SetCurrentDepth(0);
214   }
215   else
216   {
217     Warning("DisplayDetailed", "Import of R-Phi geometry failed.");
218   }
219   
220   // Projection of the different elements onto the 2D view
221   projMgr->ImportElements(fM);
222   projMgr->ImportElements(lv0TransverseMomentumDirection);
223   projMgr->ImportElements(pvlocation);
224   projMgr->ImportElements(v0location);
225
226   if (negDaughterCluster) projMgr->ImportElements(negDaughterCluster);
227   if (posDaughterCluster) projMgr->ImportElements(posDaughterCluster);
228
229   // This is the to-do list for the bending plane:
230   // 1. fix the view to orthographic XOY (no rotation allowed but moving the center ok) ->done! 
231   // 2. show axis and tickles along X and Y ->done!
232   //       -> note for the projection the cartesian scales are not very useful
233   //       -> propose a phi and R scale which rotates with a reset at 0;
234   //       -> propose a transformation for an eta scale (keep the z one);
235   // 3. show the center, the main vertex and the detectors for this view ->done!
236   // 4. show V0 direction in the bending plane with arrow length proportional to pT ->done!
237   // 5. show angles with respect to axis (phi angle) ->almost.
238   // 6. show clusters in the ITS and in the TPC associated with the daughter tracks
239   //       -> include a radius cut for plotting only ITS and TPC ->done!
240   bpViewer->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
241   bpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
242   TGLViewer *lbpGLViewer = bpViewer->GetGLViewer();
243   TGLCameraOverlay* co = lbpGLViewer->GetCameraOverlay();
244   co->SetShowOrthographic(true); //(false);
245   co->SetOrthographicMode(TGLCameraOverlay::kAxis); // ::kPlaneIntersect or ::kBar
246   // end of the bending plane part
247
248   //
249   // This part is for the decay plane view
250   //
251   pack->NewSlot()->MakeCurrent();
252   TEveViewer *dpViewer = gEve->SpawnNewViewer("V0 decay plane View");
253   TEveScene  *dpScene  = gEve->SpawnNewScene("V0 decay plane Scene");
254
255   dpViewer->AddScene(dpScene);
256
257   result = gInterpreter->ProcessLine("geom_gentle(kFALSE)");
258   if (result)
259   {
260     TEveGeoShape *geom = reinterpret_cast<TEveGeoShape*>(result);
261     geom->IncDenyDestroy();
262     geom->FindChild("TRD+TOF")->SetRnrState(kFALSE);
263     geom->FindChild("PHOS")   ->SetRnrState(kFALSE);
264     geom->FindChild("HMPID")  ->SetRnrState(kFALSE);
265     dpScene->AddElement(geom);
266   }
267   else
268   {
269     Warning("DisplayDetailed", "Import of 3D geometry failed.");
270   }
271
272   dpScene->AddElement(fM);
273   dpScene->AddElement(lv0TransverseMomentumDirection);
274   dpScene->AddElement(pvlocation);
275   dpScene->AddElement(v0location);
276   if (negDaughterCluster) dpScene->AddElement(negDaughterCluster);
277   if (posDaughterCluster) dpScene->AddElement(posDaughterCluster);
278
279   // This is the to-do list for the decay plane:
280   // 1. fix the view to decay plane (no rotation allowed but moving the center ok)
281   // 2. show V0 direction with a vertical arrow length proportional to pT -> done!
282   // 3. show the center, the main vertex and the detectors for this view -> done!
283   // 4. show the x,y and z axis and the different angles
284   //       -> this needs a referential object that we can move around
285   //          or fix to a selected point (origin being the default)
286   // 5. draw the dca between daughters and the extrapolation to the main vertex.
287   //       -> this is an issue since we only store the distance: check with J.Belikov
288   // 6. show clusters in the ITS and in the TPC associated with the daughter tracks
289   //       -> include a radius cut for plotting only ITS and TPC ->done!
290   dpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
291   TGLCamera& dpCam = dpViewer->GetGLViewer()->CurrentCamera();
292   dpCam.SetExternalCenter(kTRUE);
293   dpCam.SetCenterVec(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
294   dpCam.RotateRad(0,-TMath::Pi()/2.); // RotateRad rotates in radians (hRotate,vRotate)
295   //  Here rotate the _view_ (not the camera) by (fM->GetPhi() - TMath::Pi()/2.)
296
297   // In the end maybe truck and rotate properly...
298   //  dpCam.Truck(0,200);// Truck center wrt the view panel: (x=0 pixels, y pixels)
299   //  dpCam.Rotate(0,50,0,0); // Rotate in pixels (xDelta,yDelta)
300
301   // end of the decay plane part
302
303   //
304   // This part is for displaying the information
305   //
306   slot = pack->NewSlot();
307   
308   TEveWindowFrame *frame = slot->MakeFrame(new TRootEmbeddedCanvas());
309   frame->SetElementName("Details");
310
311   // Print and show detailed information about the V0
312   // Calculation of the invariant mass with the max prob PID hypothesis first
313   // pseudorapidity, phi angle, pt, radius, dcas
314   char info[100] = {0};
315   sprintf(info,"#phi = %.3f rad = %.1f deg",fM->GetPhi(),(180./TMath::Pi())*fM->GetPhi());
316   TLatex* ltx = new TLatex(0.05, 0.9, info);
317   ltx->SetTextSize(0.08);
318   ltx->Draw();
319
320   sprintf(info,"radius = %.3f [cm]",fM->GetRadius());
321   ltx->DrawLatex(0.05, 0.8, info);
322
323   sprintf(info,"p_{T} = %.3f [GeV/c]",fM->GetPt());
324   ltx->DrawLatex(0.05, 0.7, info);
325
326   sprintf(info,"daughters dca = %.3f [cm]",fM->GetDaughterDCA());
327   ltx->DrawLatex(0.05, 0.6, info);
328
329   sprintf(info,"#eta = - ln( tan(#theta/2) ) = %.3f",fM->GetEta());
330   ltx->DrawLatex(0.05, 0.5, info);
331
332   sprintf(info,"mass_{K^{0}_{s}} = %.3f [GeV/c^{2}]",fM->GetK0sInvMass());
333   ltx->DrawLatex(0.05, 0.3, info);
334
335   sprintf(info,"mass_{#Lambda} = %.3f [GeV/c^{2}]",fM->GetLambdaInvMass());
336   ltx->DrawLatex(0.05, 0.2, info);
337
338   sprintf(info,"mass_{#bar{#Lambda}} = %.3f [GeV/c^{2}]",fM->GetAntiLambdaInvMass());
339   ltx->DrawLatex(0.05, 0.1, info);
340
341   gEve->Redraw3D();
342 }