3 /**************************************************************************
4 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
5 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for *
6 * full copyright notice. *
7 **************************************************************************/
10 /// \file muon_digits.C
11 /// \brief Macro to visualise digits from MUON spectrometer
12 /// (both tracker and trigger).
14 /// Use muon_digits() in order to run it
16 /// Needs that alieve_init() is already called
18 /// \author P. Pillot, L. Aphecetche; Subatech
20 #if !defined(__CINT__) || defined(__MAKECINT__)
21 #include <Riostream.h>
24 #include <TEveManager.h>
25 #include <TEveQuadSet.h>
27 #include <MUON/AliMUONGeometryTransformer.h>
28 #include <MUON/AliMUONVDigit.h>
29 #include <MUON/AliMUONVDigitStore.h>
30 #include <MUON/mapping/AliMpPad.h>
31 #include <MUON/mapping/AliMpSegmentation.h>
32 #include <MUON/mapping/AliMpVSegmentation.h>
33 #include <MUON/mapping/AliMpCDB.h>
34 #include <STEER/STEER/AliRunLoader.h>
35 #include <EVE/EveBase/AliEveEventManager.h>
38 //______________________________________________________________________________
39 void add_muon_digits(TIter* next, TEveQuadSet* bending, TEveQuadSet* nonBending, Bool_t fromRaw)
42 AliMpCDB::LoadAll(kFALSE);
45 static AliMUONGeometryTransformer* gMUONGeometryTransformer = 0x0;
46 if (!gMUONGeometryTransformer)
48 AliEveEventManager::AssertGeometry();
49 gMUONGeometryTransformer = new AliMUONGeometryTransformer();
50 gMUONGeometryTransformer->LoadGeometryData();
53 // loop over digits and produce corresponding graphic objects
55 while ( ( digit = static_cast<AliMUONVDigit*>((*next)() ) ) )
57 if (!digit->IsTrigger() && !fromRaw && digit->Charge() < 1.e-3) continue;
59 Int_t detElemId = digit->DetElemId();
60 Int_t manuId = digit->ManuId();
62 const AliMpVSegmentation* vseg =
63 AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::GetCathodType(digit->Cathode()));
66 cout << Form("Could not get segmentation for DE %4d MANU %4d",detElemId,manuId) << endl;
67 continue; // should not happen, unless we got a readout error and thus a bad de,manu pair
70 AliMpPad pad = vseg->PadByLocation(manuId,digit->ManuChannel());
72 Double_t local[] = { pad.GetPositionX(), pad.GetPositionY(), 0.0 };
73 Double_t global[] = { 0.0, 0.0, 0.0 };
75 gMUONGeometryTransformer->Local2Global(detElemId,
76 local[0], local[1], local[2],
77 global[0], global[1], global[2]);
79 TEveQuadSet* pads = bending;
80 if (vseg->PlaneType()==AliMp::kNonBendingPlane) pads = nonBending;
82 pads->AddQuad(global[0]-pad.GetDimensionX(),global[1]-pad.GetDimensionY(),global[2],
83 2.*pad.GetDimensionX(),2.*pad.GetDimensionY());
85 if (fromRaw && !digit->IsTrigger()) pads->QuadValue(digit->ADC());
86 else pads->QuadValue((Int_t) digit->Charge());
91 //______________________________________________________________________________
95 AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
96 rl->LoadDigits("MUON");
97 TTree* dt = rl->GetTreeD("MUON", kFALSE);
99 AliMUONVDigitStore *digitStore = AliMUONVDigitStore::Create(*dt);
101 digitStore->Connect(*dt,0);
103 rl->UnloadDigits("MUON");
105 if (digitStore->GetSize() == 0 && !gEve->GetKeepEmptyCont()) {
110 // container for graphic representation of digits
111 TEveElementList* cont = new TEveElementList("MUON Digits");
113 TEveQuadSet* bending = new TEveQuadSet(TEveQuadSet::kQT_RectangleXY, kFALSE, 32);
114 bending->SetName("Bending");
115 bending->SetRenderMode(TEveDigitSet::kRM_Fill);
116 bending->SetPickable(kFALSE);
117 cont->AddElement(bending);
119 TEveQuadSet* nonBending = new TEveQuadSet(TEveQuadSet::kQT_RectangleXY, kFALSE, 32);
120 nonBending->SetName("Non bending");
121 nonBending->SetRenderMode(TEveDigitSet::kRM_Line);
122 nonBending->SetPickable(kFALSE);
123 cont->AddElement(nonBending);
125 // add digits to the containers
126 TIter next(digitStore->CreateIterator());
127 add_muon_digits(&next, bending, nonBending, kFALSE);
130 // set containers' title
131 Int_t nDigitB = bending->GetPlex()->Size();
132 Int_t nDigitNB = nonBending->GetPlex()->Size();
133 cont->SetTitle(Form("N=%d",nDigitB+nDigitNB));
134 bending->SetTitle(Form("N=%d",nDigitB));
135 nonBending->SetTitle(Form("N=%d",nDigitNB));
138 gStyle->SetPalette(1);
139 bending->AssertPalette();
140 nonBending->AssertPalette();
142 // add graphic containers
143 gEve->DisableRedraw();
144 gEve->AddElement(cont);
145 gEve->EnableRedraw();