]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveV0Editor.cxx
From Boris: Proper detailed-view for V0s with two GL views and details in a canvas.
[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
106 #include <TGLCamera.h>
107 #include <TGLViewer.h>
108
109 #include <TLatex.h>
110 #include <TRootEmbeddedCanvas.h>
111
112 void AliEveV0Editor::DisplayDetailed()
113 {
114   TEveWindowSlot *slot = TEveWindow::CreateWindowMainFrame();
115   TEveWindowPack *pack = slot->MakePack();
116   pack->SetShowTitleBar(kFALSE);
117   pack->SetHorizontal();
118
119   //
120   // This part is for the bending plane view
121   //
122   pack->NewSlot()->MakeCurrent();
123   TEveViewer *bpViewer = gEve->SpawnNewViewer("V0 bending plane View");
124   TEveScene  *bpScene  = gEve->SpawnNewScene("V0 bending plane Scene");
125   bpViewer->AddScene(bpScene);
126   bpScene->AddElement(fM);
127   // This is the to-do list for the bending plane:
128   // 1. fix the view to orthographic XOY (no rotation allowed but moving the center ok) --tbu 
129   // 2. show axis and tickles along X and Y
130   // 3. show the center, the main vertex and the detectors for this view;
131   // 4. show V0 direction in the bending plane with arrow length proportional to pT;
132   // 5. show angles with respect to axis (phi angle).
133   bpViewer->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
134   //  bpViewer->GetGLViewer()->GetOrthoXOYCamera()->SetEnableRotate(kFALSE);
135   bpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
136   //  bpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
137   TGLCamera& bpCam = bpViewer->GetGLViewer()->CurrentCamera();
138   bpCam.SetExternalCenter(kTRUE);
139   //  bpCam.SetEnableRotate(kFALSE);
140   bpCam.SetCenterVec(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
141   // end of the bending plane part
142
143   //
144   // This part is for the decay plane view
145   //
146   pack->NewSlot()->MakeCurrent();
147   TEveViewer *dpViewer = gEve->SpawnNewViewer("V0 decay plane View");
148   TEveScene  *dpScene  = gEve->SpawnNewScene("V0 decay plane Scene");
149   dpViewer->AddScene(dpScene);
150   dpScene->AddElement(fM);
151   // This is the to-do list for the decay plane:
152   // 1. fix the view to decay plane (no rotation allowed but moving the center ok)
153   // 2. show V0 direction with a vertical arrow length proportional to pT;
154   // 3. show the center, the main vertex and the detectors for this view;
155   // 4. show the x,y and z axis and the different angles;
156   // 5. draw the dca between daughters and the extrapolation to the main vertex.
157   dpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
158   TGLCamera& dpCam = dpViewer->GetGLViewer()->CurrentCamera();
159   dpCam.SetExternalCenter(kTRUE);
160   dpCam.SetCenterVec(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
161   // end of the decay plane part
162
163
164   // This part is for displaying the information
165   slot = pack->NewSlot();
166   
167   TEveWindowFrame *frame = slot->MakeFrame(new TRootEmbeddedCanvas());
168   frame->SetElementName("Details");
169
170   // Print and show detailed information about the V0
171   // Calculation of the invariant mass with the max prob PID hypothesis first
172   // pseudorapidity, phi angle, pt, radius, dcas
173   char info[100] = {0};
174   sprintf(info,"#phi = %.3f",fM->GetPhi());
175   TLatex* ltx = new TLatex(0.1, 0.9, info);
176   ltx->SetTextSize(0.08);
177   ltx->Draw();
178
179   sprintf(info,"radius = %.3f [cm]",fM->GetRadius());
180   ltx->DrawLatex(0.1, 0.8, info);
181
182   sprintf(info,"daughters dca = %.3f [cm]",fM->GetDaughterDCA());
183   ltx->DrawLatex(0.1, 0.7, info);
184
185   sprintf(info,"#eta = #frac{1}{2} #times Ln(#frac{E+p_{z}}{E-p_{z}} ) = %.3f",fM->GetEta());
186   ltx->DrawLatex(0.1, 0.6, info);
187
188   sprintf(info,"mass_{K^{0}_{s}} = %.3f [GeV/c^{2}]",fM->GetK0sInvMass());
189   ltx->DrawLatex(0.1, 0.3, info);
190
191   sprintf(info,"mass_{#Lambda} = %.3f [GeV/c^{2}]",fM->GetLambdaInvMass());
192   ltx->DrawLatex(0.1, 0.2, info);
193
194   sprintf(info,"mass_{#bar{#Lambda}} = %.3f [GeV/c^{2}]",fM->GetAntiLambdaInvMass());
195   ltx->DrawLatex(0.1, 0.1, info);
196
197   gEve->Redraw3D();
198 }