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