]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveV0Editor.cxx
f2e02f9e0c91ba26ea47dfc7c46d8a219765a609
[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 the bending plane view
125   //
126   pack->NewSlot()->MakeCurrent();
127   TEveViewer *bpViewer = gEve->SpawnNewViewer("V0 bending plane View");
128   TEveScene  *bpScene  = gEve->SpawnNewScene("V0 bending plane Scene");
129
130   TEveUtil::LoadMacro("geom_gentle.C");
131   Long_t result = gInterpreter->ProcessLine("geom_gentle_rphi()");
132   if (result)
133   {
134     TEveGeoShape *geomRPhi = reinterpret_cast<TEveGeoShape*>(result);
135     geomRPhi->IncDenyDestroy();
136     TEveProjectionManager *projMgr = new TEveProjectionManager();
137     projMgr->ImportElements(geomRPhi, bpScene);
138   }
139   else
140   {
141     Warning("DisplayDetailed", "Import of R-Phi geometry failed.");
142   }
143
144   bpViewer->AddScene(bpScene);
145   bpScene->AddElement(fM);
146
147   char displayInfo[100] = {0};
148   sprintf(displayInfo,"pt = %.3f",fM->GetPt());
149   TEveLine *lv0TransverseMomentumDirection = new TEveLine(displayInfo);
150   lv0TransverseMomentumDirection->SetLineColor(kOrange+8);
151   lv0TransverseMomentumDirection->SetLineWidth(2);
152   lv0TransverseMomentumDirection->SetLineStyle(2);
153   lv0TransverseMomentumDirection->SetLineWidth(2);
154   Float_t scalePt = 100.;
155   lv0TransverseMomentumDirection->SetPoint(0,fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
156   lv0TransverseMomentumDirection->SetPoint(1,scalePt*fM->fRecDecayP.fX, scalePt*fM->fRecDecayP.fY,0);
157   
158   bpScene->AddElement(lv0TransverseMomentumDirection);
159
160   TEvePointSet *pvlocation = new TEvePointSet("PV location");
161   pvlocation->SetNextPoint(fM->fRecBirthV.fX, fM->fRecBirthV.fY, fM->fRecBirthV.fZ);
162   pvlocation->SetTitle("pv location");
163   pvlocation->SetMarkerStyle(4);
164   pvlocation->SetMarkerSize(2.5);
165   pvlocation->SetMarkerColor(7);
166   bpScene->AddElement(pvlocation);
167
168   TEvePointSet *v0location = new TEvePointSet("V0 location");
169   v0location->SetNextPoint(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
170   v0location->SetTitle("v0 location");
171   v0location->SetMarkerStyle(4);
172   v0location->SetMarkerSize(2.5);
173   v0location->SetMarkerColor(kOrange+8);
174   bpScene->AddElement(v0location);
175
176   // show V0 position with a marker in orange
177   // show the pv in light blue...
178
179   // This is the to-do list for the bending plane:
180   // 1. fix the view to orthographic XOY (no rotation allowed but moving the center ok) ->done! 
181   // 2. show axis and tickles along X and Y   ->done!
182   // 3. show the center, the main vertex and the detectors for this view ->done!
183   // 4. show V0 direction in the bending plane with arrow length proportional to pT ->done!
184   // 5. show angles with respect to axis (phi angle) ->almost.
185   bpViewer->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
186   bpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
187   TGLViewer *lbpGLViewer = bpViewer->GetGLViewer();
188   TGLCameraOverlay* co = lbpGLViewer->GetCameraOverlay();
189   co->SetShowOrthographic(true); //(false);
190   co->SetOrthographicMode(TGLCameraOverlay::kAxis); // ::kPlaneIntersect or ::kBar
191   // end of the bending plane part
192   
193
194
195   //
196   // This part is for the decay plane view
197   //
198   pack->NewSlot()->MakeCurrent();
199   TEveViewer *dpViewer = gEve->SpawnNewViewer("V0 decay plane View");
200   TEveScene  *dpScene  = gEve->SpawnNewScene("V0 decay plane Scene");
201   dpViewer->AddScene(dpScene);
202   dpScene->AddElement(fM);
203   // This is the to-do list for the decay plane:
204   // 1. fix the view to decay plane (no rotation allowed but moving the center ok)
205   // 2. show V0 direction with a vertical arrow length proportional to pT;
206   // 3. show the center, the main vertex and the detectors for this view;
207   // 4. show the x,y and z axis and the different angles;
208   // 5. draw the dca between daughters and the extrapolation to the main vertex.
209   dpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
210   TGLCamera& dpCam = dpViewer->GetGLViewer()->CurrentCamera();
211   dpCam.SetExternalCenter(kTRUE);
212   dpCam.SetCenterVec(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
213   // end of the decay plane part
214
215
216   // This part is for displaying the information
217   slot = pack->NewSlot();
218   
219   TEveWindowFrame *frame = slot->MakeFrame(new TRootEmbeddedCanvas());
220   frame->SetElementName("Details");
221
222   // Print and show detailed information about the V0
223   // Calculation of the invariant mass with the max prob PID hypothesis first
224   // pseudorapidity, phi angle, pt, radius, dcas
225   char info[100] = {0};
226   sprintf(info,"#phi = %.3frad = %.1f°",fM->GetPhi(),57.296*fM->GetPhi());
227   TLatex* ltx = new TLatex(0.05, 0.9, info);
228   ltx->SetTextSize(0.08);
229   ltx->Draw();
230
231   sprintf(info,"radius = %.3f [cm]",fM->GetRadius());
232   ltx->DrawLatex(0.05, 0.8, info);
233
234   sprintf(info,"p_{T} = %.3f [GeV/c]",fM->GetPt());
235   ltx->DrawLatex(0.05, 0.7, info);
236
237   sprintf(info,"daughters dca = %.3f [cm]",fM->GetDaughterDCA());
238   ltx->DrawLatex(0.05, 0.6, info);
239
240   sprintf(info,"#eta = - ln( tan(#theta/2) ) = %.3f",fM->GetEta());
241   ltx->DrawLatex(0.05, 0.5, info);
242
243   sprintf(info,"mass_{K^{0}_{s}} = %.3f [GeV/c^{2}]",fM->GetK0sInvMass());
244   ltx->DrawLatex(0.05, 0.3, info);
245
246   sprintf(info,"mass_{#Lambda} = %.3f [GeV/c^{2}]",fM->GetLambdaInvMass());
247   ltx->DrawLatex(0.05, 0.2, info);
248
249   sprintf(info,"mass_{#bar{#Lambda}} = %.3f [GeV/c^{2}]",fM->GetAntiLambdaInvMass());
250   ltx->DrawLatex(0.05, 0.1, info);
251
252   gEve->Redraw3D();
253 }