]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveV0Editor.cxx
e23eb2857f130ce46bbb6ab6913a0fe4702b85d2
[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   //
152   // This part is for the bending plane view
153   //
154   pack->NewSlot()->MakeCurrent();
155   TEveViewer *bpViewer = gEve->SpawnNewViewer("V0 bending plane View");
156   TEveScene  *bpScene  = gEve->SpawnNewScene("V0 bending plane Scene");
157
158   TEveUtil::LoadMacro("geom_gentle.C");
159   Long_t result = gInterpreter->ProcessLine("geom_gentle_rphi()");
160   if (result)
161   {
162     TEveGeoShape *geomRPhi = reinterpret_cast<TEveGeoShape*>(result);
163     geomRPhi->IncDenyDestroy();
164     TEveProjectionManager *projMgr = new TEveProjectionManager();
165     projMgr->ImportElements(geomRPhi, bpScene);
166   }
167   else
168   {
169     Warning("DisplayDetailed", "Import of R-Phi geometry failed.");
170   }
171
172   bpViewer->AddScene(bpScene);
173   bpScene->AddElement(fM);
174   bpScene->AddElement(lv0TransverseMomentumDirection);
175   bpScene->AddElement(pvlocation);
176   bpScene->AddElement(v0location);
177
178   // This is the to-do list for the bending plane:
179   // 1. fix the view to orthographic XOY (no rotation allowed but moving the center ok) ->done! 
180   // 2. show axis and tickles along X and Y ->done!
181   //       -> note for the projection the cartesian scales are not very useful
182   //       -> propose a phi and R scale which rotates with a reset at 0;
183   //       -> propose a transformation for an eta scale (keep the z one);
184   // 3. show the center, the main vertex and the detectors for this view ->done!
185   // 4. show V0 direction in the bending plane with arrow length proportional to pT ->done!
186   // 5. show angles with respect to axis (phi angle) ->almost.
187   // 6. show clusters in the ITS and in the TPC associated with the daughter tracks
188   //       -> include a radius cut for plotting only ITS and TPC
189   bpViewer->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
190   bpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
191   TGLViewer *lbpGLViewer = bpViewer->GetGLViewer();
192   TGLCameraOverlay* co = lbpGLViewer->GetCameraOverlay();
193   co->SetShowOrthographic(true); //(false);
194   co->SetOrthographicMode(TGLCameraOverlay::kAxis); // ::kPlaneIntersect or ::kBar
195   // end of the bending plane part
196
197   //
198   // This part is for the decay plane view
199   //
200   pack->NewSlot()->MakeCurrent();
201   TEveViewer *dpViewer = gEve->SpawnNewViewer("V0 decay plane View");
202   TEveScene  *dpScene  = gEve->SpawnNewScene("V0 decay plane Scene");
203
204   dpViewer->AddScene(dpScene);
205
206   result = gInterpreter->ProcessLine("geom_gentle(kFALSE)");
207   if (result)
208   {
209     TEveGeoShape *geom = reinterpret_cast<TEveGeoShape*>(result);
210     geom->IncDenyDestroy();
211     geom->FindChild("TRD+TOF")->SetRnrState(kFALSE);
212     geom->FindChild("PHOS")   ->SetRnrState(kFALSE);
213     geom->FindChild("HMPID")  ->SetRnrState(kFALSE);
214     dpScene->AddElement(geom);
215   }
216   else
217   {
218     Warning("DisplayDetailed", "Import of 3D geometry failed.");
219   }
220
221   dpScene->AddElement(fM);
222   dpScene->AddElement(lv0TransverseMomentumDirection);
223   dpScene->AddElement(pvlocation);
224   dpScene->AddElement(v0location);
225
226   // This is the to-do list for the decay plane:
227   // 1. fix the view to decay plane (no rotation allowed but moving the center ok)
228   // 2. show V0 direction with a vertical arrow length proportional to pT -> done!
229   // 3. show the center, the main vertex and the detectors for this view -> done!
230   // 4. show the x,y and z axis and the different angles
231   //       -> this needs a referential object that we can move around
232   //          or fix to a selected point (origin being the default)
233   // 5. draw the dca between daughters and the extrapolation to the main vertex.
234   //       -> this is an issue since we only store the distance: check with J.Belikov
235   // 6. show clusters in the ITS and in the TPC associated with the daughter tracks
236   //       -> include a radius cut for plotting only ITS and TPC
237   dpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
238   TGLCamera& dpCam = dpViewer->GetGLViewer()->CurrentCamera();
239   dpCam.SetExternalCenter(kTRUE);
240   dpCam.SetCenterVec(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
241   // end of the decay plane part
242
243   //
244   // This part is for displaying the information
245   //
246   slot = pack->NewSlot();
247   
248   TEveWindowFrame *frame = slot->MakeFrame(new TRootEmbeddedCanvas());
249   frame->SetElementName("Details");
250
251   // Print and show detailed information about the V0
252   // Calculation of the invariant mass with the max prob PID hypothesis first
253   // pseudorapidity, phi angle, pt, radius, dcas
254   char info[100] = {0};
255   sprintf(info,"#phi = %.3frad = %.1f°",fM->GetPhi(),57.296*fM->GetPhi());
256   TLatex* ltx = new TLatex(0.05, 0.9, info);
257   ltx->SetTextSize(0.08);
258   ltx->Draw();
259
260   sprintf(info,"radius = %.3f [cm]",fM->GetRadius());
261   ltx->DrawLatex(0.05, 0.8, info);
262
263   sprintf(info,"p_{T} = %.3f [GeV/c]",fM->GetPt());
264   ltx->DrawLatex(0.05, 0.7, info);
265
266   sprintf(info,"daughters dca = %.3f [cm]",fM->GetDaughterDCA());
267   ltx->DrawLatex(0.05, 0.6, info);
268
269   sprintf(info,"#eta = - ln( tan(#theta/2) ) = %.3f",fM->GetEta());
270   ltx->DrawLatex(0.05, 0.5, info);
271
272   sprintf(info,"mass_{K^{0}_{s}} = %.3f [GeV/c^{2}]",fM->GetK0sInvMass());
273   ltx->DrawLatex(0.05, 0.3, info);
274
275   sprintf(info,"mass_{#Lambda} = %.3f [GeV/c^{2}]",fM->GetLambdaInvMass());
276   ltx->DrawLatex(0.05, 0.2, info);
277
278   sprintf(info,"mass_{#bar{#Lambda}} = %.3f [GeV/c^{2}]",fM->GetAntiLambdaInvMass());
279   ltx->DrawLatex(0.05, 0.1, info);
280
281   gEve->Redraw3D();
282 }