]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveCascadeEditor.cxx
MUON DAs
[u/mrichter/AliRoot.git] / EVE / EveBase / AliEveCascadeEditor.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 "AliEveCascadeEditor.h"
11 #include "AliEveCascade.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 AliEveCascade.
26 //
27
28 ClassImp(AliEveCascadeEditor)
29
30 //______________________________________________________________________________
31 AliEveCascadeEditor::AliEveCascadeEditor(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   fInfoLabelRadius(0),
36   fInfoLabelDCA(0),
37   fInfoLabelCharge(0),
38   fInfoLabelPhi(0),
39   fInfoLabelTheta(0),
40   fInfoLabelPtot(0),
41   fInfoLabelPt(0),
42   fInfoLabelEta(0),
43   fXButtonDetailedView(0),
44   fXButtonMassHyp(0)
45   // Initialize widget pointers to 0
46 {
47   // Constructor.
48
49   MakeTitle("AliEveCascade");
50
51   fInfoLabelRadius = new TGLabel(this);
52   fInfoLabelRadius->SetTextJustify(kTextLeft);
53   AddFrame(fInfoLabelRadius, new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 8, 0, 2, 0));
54   
55   fInfoLabelDCA = new TGLabel(this);
56   fInfoLabelDCA->SetTextJustify(kTextLeft);
57   AddFrame(fInfoLabelDCA,    new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 8, 0, 2, 0));
58
59   fInfoLabelCharge = new TGLabel(this);
60   fInfoLabelCharge->SetTextJustify(kTextLeft);
61   AddFrame(fInfoLabelCharge, new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 8, 0, 2, 0));
62
63   fInfoLabelPhi = new TGLabel(this);
64   fInfoLabelPhi->SetTextJustify(kTextLeft);
65   AddFrame(fInfoLabelPhi,    new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 8, 0, 2, 0));
66   
67   fInfoLabelTheta = new TGLabel(this);
68   fInfoLabelTheta->SetTextJustify(kTextLeft);
69   AddFrame(fInfoLabelTheta,  new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 8, 0, 2, 0));
70   
71   fInfoLabelPtot = new TGLabel(this);
72   fInfoLabelPtot->SetTextJustify(kTextLeft);
73   AddFrame(fInfoLabelPtot,   new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 8, 0, 2, 0));
74   
75   fInfoLabelPt = new TGLabel(this);
76   fInfoLabelPt->SetTextJustify(kTextLeft);
77   AddFrame(fInfoLabelPt,     new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 8, 0, 2, 0));
78   
79   fInfoLabelEta = new TGLabel(this);
80   fInfoLabelEta->SetTextJustify(kTextLeft);
81   AddFrame( fInfoLabelEta,   new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 8, 0, 2, 0));
82
83   
84   fXButtonDetailedView = new TGTextButton(this, "Detailed View");
85   AddFrame(fXButtonDetailedView, new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 1, 1, 0, 0));
86   fXButtonDetailedView->Connect("Clicked()", "AliEveCascadeEditor", this, "DisplayDetailed()");
87   
88   fXButtonMassHyp = new TGTextButton(this, "Mass Hypotheses according to charge");
89   AddFrame(fXButtonMassHyp, new TGLayoutHints(kLHintsLeft|kLHintsExpandX, 1, 1, 0, 0));
90   fXButtonMassHyp->Connect("Clicked()", "AliEveCascadeEditor", this, "DisplayMassHyp()");
91   
92 }
93
94 /******************************************************************************/
95
96 //______________________________________________________________________________
97 void AliEveCascadeEditor::SetModel(TObject* obj)
98 {
99   // Set model object.
100
101   fM = static_cast<AliEveCascade*>(obj);
102
103   // Set values of widgets
104   fInfoLabelRadius->SetText(Form("Radius = %f cm",     fM->GetRadius() ));
105   fInfoLabelDCA   ->SetText(Form("DCA (Xi dghters) = %f cm", fM->GetDaughterDCA()));
106   fInfoLabelCharge->SetText(Form("Charge = %d",        fM->GetCharge() ));
107   fInfoLabelPhi   ->SetText(Form("Phi     = %f deg",   fM->GetPhi()   * 180.0/TMath::Pi()    ));
108   fInfoLabelTheta ->SetText(Form("Theta  = %f deg",    fM->GetTheta() * 180.0/TMath::Pi()    ));
109   fInfoLabelPtot  ->SetText(Form("Ptot    = %f GeV/c", fM->GetPtot()   ));
110   fInfoLabelPt    ->SetText(Form("Pt      = %f GeV/c", fM->GetPt()     ));
111   fInfoLabelEta   ->SetText(Form("Eta    = %f",        fM->GetEta()    ));
112 }
113
114 /******************************************************************************/
115
116 // Implements callback/slot methods
117
118 //______________________________________________________________________________
119 // void AliEveCascadeEditor::DoXYZZ()
120 // {
121 //    // Slot for XYZZ.
122 //
123 //    fM->SetXYZZ(fXYZZ->GetValue());
124 //    Update();
125 // }
126
127
128 #include <TEveManager.h>
129 #include <TEveWindow.h>
130 #include <TEveViewer.h>
131 #include <TEveScene.h>
132 #include <TEveGeoNode.h>
133 #include <TEveProjectionManager.h>
134
135 #include <TGLCamera.h>
136 #include <TGLViewer.h>
137 #include "TGLCameraOverlay.h"
138
139 #include <TLatex.h>
140 #include <TRootEmbeddedCanvas.h>
141 #include <TInterpreter.h>
142
143
144 void AliEveCascadeEditor::DisplayDetailed()
145 {
146   // Display a detailed view (bending plane + transverse plane)
147         
148   printf("\n--> Detailed View :\n");
149   TEveWindowSlot *slot = TEveWindow::CreateWindowMainFrame();
150   TEveWindowPack *pack = slot->MakePack();
151   pack->SetShowTitleBar(kFALSE);
152   pack->SetHorizontal();
153
154   //----------
155   // Part 0 : get the different objects to display
156   //
157   
158   TEvePointSet *lPrimVtxLocation = new TEvePointSet("Prim Vtx location");
159    lPrimVtxLocation->SetNextPoint(fM->fRecBirthV.fX, fM->fRecBirthV.fY, fM->fRecBirthV.fZ);
160    lPrimVtxLocation->SetTitle("prim vtx location");
161    lPrimVtxLocation->SetMarkerStyle(4);
162    lPrimVtxLocation->SetMarkerSize(2.5);
163    lPrimVtxLocation->SetMarkerColor(kCyan);
164
165   TEvePointSet *lCascadeLocation = new TEvePointSet("Cascade decay location");
166    lCascadeLocation->SetNextPoint(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
167    lCascadeLocation->SetTitle("Cascade decay location");
168    lCascadeLocation->SetMarkerStyle(4);
169    lCascadeLocation->SetMarkerSize(2.5);
170    lCascadeLocation->SetMarkerColor(kMagenta-9);
171
172   TEveUtil::LoadMacro("clusters_from_index.C");
173   
174   AliEveTrack *bacTrack = fM->GetBacTrack();
175   AliEveTrack *negTrack = fM->GetNegTrack();
176   AliEveTrack *posTrack = fM->GetPosTrack();
177   
178   char macroWithIndex[100] = {0};
179   Int_t daughterIndex = 0;
180   
181   TEvePointSet *bacDaughterCluster = 0;
182   TEvePointSet *negDaughterCluster = 0;
183   TEvePointSet *posDaughterCluster = 0;
184   
185   // Clusters linked with the bachelor track
186   daughterIndex = bacTrack->GetIndex();
187   snprintf(macroWithIndex,100,"clusters_from_index(%d)",daughterIndex);
188   Long_t bacResult = gInterpreter->ProcessLine(macroWithIndex);
189   if (bacResult) {
190           bacDaughterCluster = reinterpret_cast<TEvePointSet*>(bacResult);
191           if (bacDaughterCluster){
192                   bacDaughterCluster->SetMarkerStyle(4);        
193                   bacDaughterCluster->SetMarkerSize(1.5);       
194                   bacDaughterCluster->SetMarkerColor(kMagenta);
195           }
196   }
197   else
198   {
199           Warning("DisplayDetailed", "Cascade : Import of bachelor clusters failed.");
200   }
201   
202   // Clusters linked with the negative daughter track (V0)
203   daughterIndex = negTrack->GetIndex();
204   snprintf(macroWithIndex,100,"clusters_from_index(%d)",daughterIndex);
205   Long_t negResult = gInterpreter->ProcessLine(macroWithIndex);
206   if (negResult) {
207           negDaughterCluster = reinterpret_cast<TEvePointSet*>(negResult);
208           if (negDaughterCluster){
209                   negDaughterCluster->SetMarkerStyle(4);        
210                   negDaughterCluster->SetMarkerSize(1.5);       
211                   negDaughterCluster->SetMarkerColor(kBlue+3);
212           }
213   }
214   else
215   {
216           Warning("DisplayDetailed", "Cascade : Import of negative daughter's clusters failed.");
217   }
218
219   // Clusters linked with the positive daughter track (V0)
220   daughterIndex = posTrack->GetIndex();
221   snprintf(macroWithIndex,100,"clusters_from_index(%d)",daughterIndex);
222   Long_t posResult = gInterpreter->ProcessLine(macroWithIndex);
223   if (posResult) {
224           posDaughterCluster = reinterpret_cast<TEvePointSet*>(posResult);
225           if (posDaughterCluster){
226                   posDaughterCluster->SetMarkerStyle(4);        
227                   posDaughterCluster->SetMarkerSize(1.5);       
228                   posDaughterCluster->SetMarkerColor(kRed+3);
229           }
230   }
231   else
232   {
233           Warning("DisplayDetailed", "Cascade : Import of positive daughter's clusters failed.");
234   }
235
236   
237   
238   
239   //----------
240   // Part 1 : bending plane view
241   //
242   pack->NewSlot()->MakeCurrent();
243   TEveViewer *bpViewer = gEve->SpawnNewViewer("Cascade bending plane");
244   TEveScene  *bpScene  = gEve->SpawnNewScene("Cascade bending plane Scene");
245
246   TEveProjectionManager *projMgr = new TEveProjectionManager(TEveProjection::kPT_RPhi);
247   bpScene->AddElement(projMgr);
248
249   TEveUtil::LoadMacro("geom_gentle.C");
250   Long_t result = gInterpreter->ProcessLine("geom_gentle_rphi()");
251   if (result)
252   {
253           TEveGeoShape *geomRPhi = reinterpret_cast<TEveGeoShape*>(result);
254           geomRPhi->IncDenyDestroy();
255           projMgr->SetCurrentDepth(-10);
256           projMgr->ImportElements(geomRPhi);
257           projMgr->SetCurrentDepth(0);
258   }
259   else
260   {
261           Warning("DisplayDetailed", "Cascade (bending plane view) : Import of R-Phi geometry failed.");
262   }
263
264   projMgr->ImportElements(fM);
265   projMgr->ImportElements(lPrimVtxLocation);
266   projMgr->ImportElements(lCascadeLocation);
267   bpViewer->AddScene(bpScene);
268
269   if (bacDaughterCluster) projMgr->ImportElements(bacDaughterCluster);
270   if (negDaughterCluster) projMgr->ImportElements(negDaughterCluster);
271   if (posDaughterCluster) projMgr->ImportElements(posDaughterCluster);
272   
273   // This is the to-do list for the bending plane:
274   // 1. show the V0 daughters track + corresponding clusters                        -> to do ...
275   // 2. show axis and tickles along X and Y                                         -> done!
276   //       -> note for the projection the cartesian scales are not very useful
277   //       -> propose a phi and R scale which rotates with a reset at 0;
278   //       -> propose a transformation for an eta scale (keep the z one);
279   // 3. show the center, the main vertex and the detectors for this view            -> done!
280   // 4. show V0 direction in the bending plane with arrow length proportional to pT -> done!
281   // 5. show angles with respect to axis (phi angle)                                -> almost.
282   // 6. show clusters in the ITS and in the TPC associated with the daughter tracks
283   //       -> include a radius cut for plotting only ITS and TPC                    -> done!
284   bpViewer->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
285   bpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
286   TGLViewer *lbpGLViewer = bpViewer->GetGLViewer();
287   TGLCameraOverlay* co = lbpGLViewer->GetCameraOverlay();
288   co->SetShowOrthographic(true); //(false);
289   co->SetOrthographicMode(TGLCameraOverlay::kAxis); // ::kPlaneIntersect or ::kBar
290   // lbpGLViewer->RequestDraw();
291   // end of the bending plane part
292
293   
294   
295   
296   //----------
297   // Part 2 : decay plane view
298   //
299   pack->NewSlot()->MakeCurrent();
300   TEveViewer *dpViewer = gEve->SpawnNewViewer("Cascade decay plane");
301   TEveScene  *dpScene  = gEve->SpawnNewScene("Cascade decay plane Scene");
302
303   
304
305   result = gInterpreter->ProcessLine("geom_gentle(kFALSE)");
306   if (result)
307   {
308           TEveGeoShape *geom = reinterpret_cast<TEveGeoShape*>(result);
309           geom->IncDenyDestroy();
310           geom->FindChild("TRD+TOF")->SetRnrState(kFALSE);
311           geom->FindChild("PHOS")   ->SetRnrState(kFALSE);
312           geom->FindChild("HMPID")  ->SetRnrState(kFALSE);
313           dpScene->AddElement(geom);
314   }
315   else
316   {
317           Warning("DisplayDetailed", "Cascade (decay plane view) : Import of 3D geometry  failed.");
318   }
319
320   dpViewer->AddScene(dpScene);
321   dpScene->AddElement(fM);
322   dpScene->AddElement(lPrimVtxLocation);
323   dpScene->AddElement(lCascadeLocation);
324   if (bacDaughterCluster) dpScene->AddElement(bacDaughterCluster);
325   if (negDaughterCluster) dpScene->AddElement(negDaughterCluster);
326   if (posDaughterCluster) dpScene->AddElement(posDaughterCluster);
327   
328   
329
330   // This is the to-do list for the decay plane:
331   // 1. fix the view to decay plane (no rotation allowed but moving the center ok)
332   // 3. show the center, the main vertex and the detectors for this view -> done!
333   // 4. show the x,y and z axis and the different angles
334   //       -> this needs a referential object that we can move around
335   //          or fix to a selected point (origin being the default)
336   // 5. draw the dca between Xi daughters and the extrapolation to the main vertex.
337   //       -> this is an issue since we only store the distance: check with J.Belikov
338   dpViewer->GetGLViewer()->ResetCamerasAfterNextUpdate();
339   TGLCamera& dpCam = dpViewer->GetGLViewer()->CurrentCamera();
340   dpCam.SetExternalCenter(kTRUE);
341   dpCam.SetCenterVec(fM->fRecDecayV.fX, fM->fRecDecayV.fY, fM->fRecDecayV.fZ);
342   dpCam.RotateRad(0,-TMath::Pi()/2.); // RotateRad rotates in radians (hRotate,vRotate)
343   //  Here rotate the _view_ (not the camera) by (fM->GetPhi() - TMath::Pi()/2.)
344
345   // In the end maybe truck and rotate properly...
346   //  dpCam.Truck(0,200);// Truck center wrt the view panel: (x=0 pixels, y pixels)
347   //  dpCam.Rotate(0,50,0,0); // Rotate in pixels (xDelta,yDelta)
348
349   // end of the decay plane part
350
351   
352   
353   
354   //----------
355   // Part 3 : displaying extra information
356   //
357   slot = pack->NewSlot();
358   
359   TEveWindowFrame *frame = slot->MakeFrame(new TRootEmbeddedCanvas());
360   frame->SetElementName("Details");
361
362   // Print and show detailed information about the V0
363   // Calculation of the invariant mass with the max prob PID hypothesis first
364   // pseudorapidity, phi angle, pt, radius, dcas
365   char info[100] = {0};
366   snprintf(info,100,"#phi = %.3f rad = %.1f deg",fM->GetPhi(),(180./TMath::Pi())*fM->GetPhi());
367   TLatex* ltx = new TLatex(0.05, 0.9, info);
368   ltx->SetTextSize(0.08);
369   ltx->Draw();
370
371   snprintf(info,100,"radius = %.3f [cm]",fM->GetRadius());
372   ltx->DrawLatex(0.05, 0.8, info);
373
374   snprintf(info,100,"p_{T} = %.3f [GeV/c]",fM->GetPt());
375   ltx->DrawLatex(0.05, 0.7, info);
376
377   snprintf(info,100,"Xi dghtrs dca = %.4f [cm]",fM->GetDaughterDCA());
378   ltx->DrawLatex(0.05, 0.6, info);
379
380   snprintf(info,100,"#eta = - ln( tan(#theta/2) ) = %.3f",fM->GetEta());
381   ltx->DrawLatex(0.05, 0.5, info);
382
383   
384   if(fM->GetCharge() < 0){
385           snprintf(info,100,"mass_{#Xi^{-}} : %.3f [GeV/c^{2}]", fM->GetXiMinusInvMass()    );
386           ltx->DrawLatex(0.05, 0.3, info);
387           snprintf(info,100,"mass_{#Omega^{-}} : %.3f [GeV/c^{2}]", fM->GetOmegaMinusInvMass() );
388           ltx->DrawLatex(0.05, 0.2, info);
389   }
390   else {
391           snprintf(info,100,"mass_{#Xi^{+}} : %.3f [GeV/c^{2}]", fM->GetXiPlusInvMass()     );
392           ltx->DrawLatex(0.05, 0.3, info);
393           snprintf(info,100,"mass_{#Omega^{+}} : %.3f [GeV/c^{2}]", fM->GetOmegaPlusInvMass()  );
394           ltx->DrawLatex(0.05, 0.2, info);
395   }
396   
397   gEve->Redraw3D();
398 }
399
400 void AliEveCascadeEditor::DisplayMassHyp()
401 {       
402    // Printf the proper mass hypotheses, according to the charge of the Bachelor
403         printf("\n--> Mass Hyp:");
404         
405         if(fM->GetCharge() < 0){
406                 printf("Xi-    mass hyp : %f GeV/c^2}\n", fM->GetXiMinusInvMass()    );
407                 printf("Omega- mass hyp : %f GeV/c^2\n\n", fM->GetOmegaMinusInvMass() );
408         }
409         else{
410                 printf("Xi+    mass hyp : %f GeV/c^2\n", fM->GetXiPlusInvMass()     );
411                 printf("Omega+ mass hyp : %f GeV/c^2\n\n", fM->GetOmegaPlusInvMass()  );
412         }
413 }