]>
Commit | Line | Data |
---|---|---|
20dae051 | 1 | // $Id$ |
6226bf2b | 2 | // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 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 | ||
ba978640 | 10 | #if !defined(__CINT__) || defined(__MAKECINT__) |
f944f125 | 11 | #include <iostream> |
12 | ||
13 | #include <TGeoGlobalMagField.h> | |
14 | #include <TGeoManager.h> | |
15 | #include <TMap.h> | |
16 | #include <TString.h> | |
17 | #include <TSystem.h> | |
ba978640 | 18 | #include <TStyle.h> |
19 | #include <TTree.h> | |
f944f125 | 20 | #include <TRegexp.h> |
ba978640 | 21 | #include <TEveManager.h> |
22 | #include <TEveElement.h> | |
23 | #include <TEvePointSet.h> | |
24 | #include <TEveTreeTools.h> | |
25 | ||
6c49a8e1 | 26 | #include <AliAltroRawStreamV3.h> |
6c49a8e1 | 27 | #include <AliEveEventManager.h> |
28 | #include <AliEveTPCData.h> | |
29 | #include <AliEveTPCSector2D.h> | |
30 | #include <AliEveTPCSector3D.h> | |
f944f125 | 31 | #include <AliGeomManager.h> |
32 | #include <AliMagF.h> | |
33 | #include <AliRawReader.h> | |
34 | ||
35 | #include <AliGRPObject.h> | |
36 | #include <AliCDBManager.h> | |
37 | #include <AliCDBPath.h> | |
38 | #include <AliCDBEntry.h> | |
39 | ||
40 | #include <STEER/STEER/AliGRPManager.h> | |
41 | ||
42 | #include <TPC/AliTPCcalibDB.h> | |
43 | #include <TPC/AliTPCTransform.h> | |
44 | #include <TPC/AliTPCclusterMI.h> | |
45 | #include <TPC/AliTPCClustersRow.h> | |
46 | ||
47 | #include <TPC/AliTPCParam.h> | |
48 | #include <TPC/AliTPCRawStreamV3.h> | |
49 | #include <TPC/AliTPCRecoParam.h> | |
50 | #include <TPC/AliTPCReconstructor.h> | |
51 | ||
52 | #include <HLT/BASE/AliHLTOUT.h> | |
53 | #include <HLT/BASE/AliHLTSystem.h> | |
54 | #include <HLT/BASE/AliHLTPluginBase.h> | |
55 | ||
56 | #include <RAW/AliRawHLTManager.h> | |
ba978640 | 57 | #endif |
58 | ||
f944f125 | 59 | // read TPC clusters compressed by HLT from RawReader and return output in OutputTree |
60 | TTree* readHLTClusters(AliRawReader *RawReader); | |
61 | ||
62 | // add the HLT clusters from clustersTree to visualization | |
63 | void renderHLTClusters(TTree* clustersTree); | |
64 | ||
65 | ||
6226bf2b | 66 | // Macro to visualise rootified raw-data from TPC. |
67 | // | |
68 | // Use tpc_raw(Int_t mode) in order to run it | |
69 | // Needs that alieve_init() is already called | |
70 | // mode = 1 - show only 2D sectors | |
71 | // mode = 2 - show only 3D sectors | |
72 | // mode = 3 - show both 2D and 3D sectors | |
6226bf2b | 73 | void tpc_raw(Int_t mode = 3) |
74 | { | |
75 | gStyle->SetPalette(1, 0); | |
76 | ||
f944f125 | 77 | AliEveEventManager::AssertGeometry(); |
78 | AliEveEventManager::AssertMagField(); | |
79 | ||
6226bf2b | 80 | AliRawReader *reader = AliEveEventManager::AssertRawReader(); |
81 | reader->Reset(); | |
f944f125 | 82 | |
375fd576 | 83 | AliTPCRawStreamV3 input(reader); |
6226bf2b | 84 | reader->Select("TPC"); // ("TPC", firstRCU, lastRCU); |
85 | ||
4d62585e | 86 | AliEveTPCData *x = new AliEveTPCData; |
6226bf2b | 87 | // x->SetLoadPedestal(5); |
88 | x->SetLoadThreshold(5); | |
89 | x->SetAutoPedestal(kTRUE); | |
90 | ||
6226bf2b | 91 | x->LoadRaw(input, kTRUE, kTRUE); |
92 | ||
da94dab4 | 93 | gEve->DisableRedraw(); |
94 | ||
de33999e | 95 | TEveElementList* sec2d = new TEveElementList("TPC 2D"); |
96 | gEve->AddElement(sec2d); | |
97 | ||
98 | TEveElementList* sec3d = new TEveElementList("TPC 3D"); | |
99 | gEve->AddElement(sec3d); | |
100 | ||
ba978640 | 101 | AliEveTPCSector2D *s; |
102 | AliEveTPCSector3D *t; | |
103 | ||
6226bf2b | 104 | for (Int_t i=0; i<=35; ++i) { |
105 | if (mode & 1) { | |
de33999e | 106 | s = new AliEveTPCSector2D(Form("2D sector %d",i)); |
6226bf2b | 107 | s->SetSectorID(i); |
108 | s->SetAutoTrans(kTRUE); // place on proper 3D coordinates | |
109 | s->SetDataSource(x); | |
110 | s->SetFrameColor(36); | |
de33999e | 111 | sec2d->AddElement(s); |
6226bf2b | 112 | s->IncRTS(); |
113 | } | |
114 | if (mode & 2) { | |
de33999e | 115 | t = new AliEveTPCSector3D(Form("3D sector %d",i)); |
6226bf2b | 116 | t->SetSectorID(i); |
117 | t->SetAutoTrans(kTRUE); | |
118 | t->SetDataSource(x); | |
de33999e | 119 | sec3d->AddElement(t); |
6226bf2b | 120 | t->IncRTS(); |
121 | } | |
122 | } | |
f944f125 | 123 | |
124 | // Display TPC clusters compressed by HLT | |
125 | TTree* hltClustersTree = readHLTClusters(reader); // read HLT compressed clusters from TPC from raw reader and output them in hltClustersTree | |
126 | if(hltClustersTree) renderHLTClusters(hltClustersTree); | |
127 | ||
6226bf2b | 128 | gEve->EnableRedraw(); |
f944f125 | 129 | gEve->Redraw3D(); |
130 | ||
131 | } | |
132 | ||
133 | // read TPC Clusters compressed by HLT from a raw filename | |
134 | TTree* readHLTClusters(AliRawReader *fRawReader) | |
135 | { | |
136 | /************************* | |
137 | * Raw IO Stuff | |
138 | **************************/ | |
139 | AliRawReader* rawReader = AliRawHLTManager::CreateRawReaderHLT(fRawReader, "TPC"); | |
140 | rawReader->Select("TPC"); | |
141 | ||
142 | /************************* | |
143 | * HLT IO Stuff | |
144 | **************************/ | |
145 | AliHLTSystem* hltSystem=AliHLTPluginBase::GetInstance(); | |
146 | ||
147 | AliHLTOUT* hltOut = AliHLTOUT::New(rawReader); | |
148 | hltOut->Init(); | |
149 | hltSystem->InitHLTOUT(hltOut); | |
150 | ||
151 | /************************* | |
152 | * GRP Stuff | |
153 | **************************/ | |
154 | AliGRPManager manGRP; | |
155 | AliRecoParam* recoParam = AliEveEventManager::AssertRecoParams(); | |
156 | ||
157 | AliEveEventManager* curEventMan = (AliEveEventManager*)gEve->GetCurrentEvent(); | |
158 | ||
159 | const AliRunInfo* runInfo = manGRP.GetRunInfo(); | |
160 | const THashTable* cosmicTriggers = manGRP.GetCosmicTriggers(); | |
161 | const AliEventInfo* eventinfo = curEventMan->GetEventInfo(); | |
162 | ||
163 | recoParam->SetEventSpecie(runInfo,*eventinfo,cosmicTriggers); | |
164 | ||
165 | AliTPCRecoParam *tpcRecoParam = (AliTPCRecoParam*)recoParam->GetDetRecoParam(1); // TPC has index 1 | |
166 | tpcRecoParam->SetUseHLTClusters(3); // reconstruct only HLT clusters | |
167 | ||
168 | /************************** | |
169 | * Reconstruction of Clusters | |
170 | **************************/ | |
171 | TTree* outputClustersTree = new TTree; | |
172 | ||
173 | AliTPCReconstructor *reconstructor = new AliTPCReconstructor(); | |
174 | reconstructor->SetOption("useHLT"); | |
175 | reconstructor->CreateTracker(); // this will set the option to reconstruct for only HLT clusters | |
176 | ||
177 | reconstructor->SetRecoParam(tpcRecoParam); | |
178 | reconstructor->SetRunInfo((AliRunInfo*)runInfo); | |
179 | reconstructor->SetEventInfo((AliEventInfo*)eventinfo); | |
180 | ||
181 | reconstructor->Reconstruct(rawReader, outputClustersTree); | |
182 | ||
183 | delete reconstructor; | |
184 | ||
185 | hltSystem->ReleaseHLTOUT(hltOut); | |
186 | ||
187 | return outputClustersTree; | |
188 | } | |
189 | ||
190 | void renderHLTClusters(TTree* clustersTree) | |
191 | { | |
192 | ||
193 | /************************** | |
194 | * Visualization of Clusters | |
195 | **************************/ | |
196 | const Int_t kMaxCl=100*160; | |
197 | Int_t fNColorBins = 5; | |
198 | ||
199 | TEvePointSet* clusters = new TEvePointSet(kMaxCl); | |
200 | clusters->SetOwnIds(kTRUE); | |
201 | ||
202 | TEvePointSetArray * cc = new TEvePointSetArray("TPC Clusters Colorized"); | |
203 | cc->SetMainColor(kRed); | |
204 | cc->SetMarkerStyle(4); | |
205 | cc->SetMarkerSize(0.4); | |
206 | cc->InitBins("Cluster Charge", fNColorBins, 0., fNColorBins*60.); | |
207 | ||
208 | cc->GetBin(0)->SetMainColor(kGray); | |
209 | cc->GetBin(0)->SetMarkerSize(0.4); | |
210 | cc->GetBin(1)->SetMainColor(kBlue); | |
211 | cc->GetBin(1)->SetMarkerSize(0.42); | |
212 | cc->GetBin(2)->SetMainColor(kCyan); | |
213 | cc->GetBin(2)->SetMarkerSize(0.44); | |
214 | cc->GetBin(3)->SetMainColor(kGreen); | |
215 | cc->GetBin(3)->SetMarkerSize(0.46); | |
216 | cc->GetBin(4)->SetMainColor(kYellow); | |
217 | cc->GetBin(4)->SetMarkerSize(0.48); | |
218 | cc->GetBin(5)->SetMainColor(kRed); | |
219 | cc->GetBin(5)->SetMarkerSize(0.50); | |
220 | cc->GetBin(6)->SetMainColor(kMagenta); | |
221 | cc->GetBin(6)->SetMarkerSize(0.52); | |
222 | ||
223 | ||
224 | // Loop over clusters | |
225 | Int_t nentries = clustersTree->GetEntriesFast(); | |
226 | ||
227 | AliTPCClustersRow *clrow = new AliTPCClustersRow(); | |
228 | clrow->SetClass("AliTPCclusterMI"); | |
229 | //clrow->SetArray(kMaxCl); | |
230 | clustersTree->SetBranchAddress("Segment", &clrow); | |
231 | ||
232 | for (Int_t i=0; i<nentries; i++) { | |
233 | if (!clustersTree->GetEvent(i)) continue; | |
234 | ||
235 | TClonesArray *cl = clrow->GetArray(); | |
236 | Int_t ncl = cl->GetEntriesFast(); | |
237 | ||
238 | while (ncl--) | |
239 | { | |
240 | AliTPCclusterMI* clusterMI = (AliTPCclusterMI*) cl->At(ncl); | |
241 | ||
242 | AliCluster *c = (AliCluster*) cl->UncheckedAt(ncl); | |
243 | Float_t g[3]; //global coordinates | |
244 | c->GetGlobalXYZ(g); | |
245 | ||
246 | cc->Fill(g[0], g[1], g[2], clusterMI->GetQ()); | |
247 | clusters->SetNextPoint(g[0], g[1], g[2]); | |
248 | AliCluster *atp = new AliCluster(*clusterMI); | |
249 | clusters->SetPointId(atp); | |
250 | ||
251 | } | |
252 | cl->Clear(); | |
253 | } | |
254 | ||
255 | delete clrow; | |
256 | ||
257 | clusters->SetName("TPC Clusters"); | |
258 | clusters->SetTitle(Form("N=%d", clusters->Size())); | |
259 | ||
260 | const TString viz_tag("REC Clusters TPC"); // to be changed | |
261 | clusters->ApplyVizTag(viz_tag, "Clusters"); | |
262 | ||
263 | cc->SetRnrSelf(kTRUE); | |
264 | ||
265 | gEve->AddElement(cc); | |
266 | ||
267 | return; | |
6226bf2b | 268 | } |