]>
Commit | Line | Data |
---|---|---|
1 | /************************************************************************** | |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | /* $Id: AliTRDgtuSim.cxx 28397 2008-09-02 09:33:00Z cblume $ */ | |
17 | ||
18 | //////////////////////////////////////////////////////////////////////////// | |
19 | // // | |
20 | // GTU simulation // | |
21 | // // | |
22 | // Authors: J. Klein (Jochen.Klein@cern.ch) // | |
23 | // // | |
24 | //////////////////////////////////////////////////////////////////////////// | |
25 | ||
26 | #include <stdio.h> | |
27 | #include <fstream> | |
28 | #include <string> | |
29 | ||
30 | #include "TFile.h" | |
31 | #include "TROOT.h" | |
32 | #include "TObject.h" | |
33 | #include "TClonesArray.h" | |
34 | ||
35 | #include "AliRun.h" | |
36 | #include "AliRunLoader.h" | |
37 | #include "AliLoader.h" | |
38 | #include "AliTreeLoader.h" | |
39 | #include "AliLog.h" | |
40 | #include "AliESDTrdTrack.h" | |
41 | #include "AliESDTrdTracklet.h" | |
42 | ||
43 | #include "AliTRDgtuSim.h" | |
44 | #include "AliTRDfeeParam.h" | |
45 | #include "AliTRDgtuTMU.h" | |
46 | #include "AliTRDtrackGTU.h" | |
47 | #include "AliTRDtrackletWord.h" | |
48 | #include "AliTRDtrackletMCM.h" | |
49 | #include "AliESDEvent.h" | |
50 | ||
51 | ClassImp(AliTRDgtuSim) | |
52 | ||
53 | AliTRDgtuSim::AliTRDgtuSim(AliRunLoader *rl) | |
54 | : TObject(), | |
55 | fRunLoader(rl), | |
56 | fFeeParam(AliTRDfeeParam::Instance()), | |
57 | fTMU(0x0), | |
58 | fTrackletArray(0x0) | |
59 | { | |
60 | ||
61 | } | |
62 | ||
63 | AliTRDgtuSim::~AliTRDgtuSim() | |
64 | { | |
65 | // destructor | |
66 | ||
67 | if (fTrackletArray) | |
68 | fTrackletArray->Clear(); | |
69 | delete fTrackletArray; | |
70 | } | |
71 | ||
72 | Bool_t AliTRDgtuSim::RunGTUFromTrackletFile(TString filename, Int_t event, Int_t noev) | |
73 | { | |
74 | // run the GTU from a file of tracklets | |
75 | // used for comparison to VHDL simulation | |
76 | ||
77 | ifstream input(filename.Data()); | |
78 | ||
79 | std::string str; | |
80 | TString string; | |
81 | int lineno = 0; | |
82 | ||
83 | Int_t iEventPrev = -1; | |
84 | Int_t iStackPrev = -1; | |
85 | Int_t iSecPrev = -1; | |
86 | Int_t iSec = -1; | |
87 | Int_t iStack = -1; | |
88 | Int_t iLink = -1; | |
89 | Int_t iEvent = -1; | |
90 | Int_t evcnt = -1; | |
91 | ||
92 | fTMU = 0x0; | |
93 | ||
94 | TClonesArray trklArray("AliTRDtrackletWord", 100); | |
95 | TClonesArray trklArrayGTU("AliTRDtrackletGTU", 100); | |
96 | ||
97 | AliDebug(1, Form("--------- Reading from %s ----------", filename.Data())); | |
98 | while (getline(input, str)) { | |
99 | lineno++; | |
100 | string = str; | |
101 | ||
102 | TObjArray *tokens = string.Tokenize(" "); | |
103 | if (tokens->GetEntriesFast() < 7) { | |
104 | AliWarning(Form("Invalid input in line %i, too few parameters", lineno)); | |
105 | continue; | |
106 | } | |
107 | ||
108 | if ( ((TObjString*) tokens->At(0))->GetString().Atoi() < event) | |
109 | continue; | |
110 | ||
111 | iEvent = ((TObjString*) tokens->At(0))->GetString().Atoi(); | |
112 | iSec = ((TObjString*) tokens->At(1))->GetString().Atoi(); | |
113 | iStack = ((TObjString*) tokens->At(2))->GetString().Atoi(); | |
114 | iLink = 2 * ((TObjString*) tokens->At(3))->GetString().Atoi() + ((TObjString*) tokens->At(4))->GetString().Atoi(); | |
115 | ||
116 | if ((iEvent != iEventPrev) || | |
117 | (iStack != iStackPrev) || | |
118 | (iSec != iSecPrev)) { | |
119 | if(fTMU) { | |
120 | TList *listOfTracks = new TList(); | |
121 | fTMU->SetStack(iStackPrev); | |
122 | fTMU->SetSector(iSecPrev); | |
123 | fTMU->RunTMU(listOfTracks); | |
124 | AliDebug(1,Form("--- There are %i tracks. Writing ...", listOfTracks->GetEntries())); | |
125 | WriteTracksToDataFile(listOfTracks, iEventPrev); | |
126 | if (listOfTracks->GetEntries() > 0) | |
127 | AliDebug(2,Form(" %4.1f GeV/c", ((AliTRDtrackGTU*) listOfTracks->At(0))->GetPt() )); | |
128 | delete fTMU; | |
129 | fTMU = new AliTRDgtuTMU(); | |
130 | delete listOfTracks; | |
131 | listOfTracks = 0x0; | |
132 | } else { | |
133 | fTMU = new AliTRDgtuTMU(); | |
134 | } | |
135 | iStackPrev = iStack; | |
136 | iSecPrev = iSec; | |
137 | iEventPrev = iEvent; | |
138 | evcnt++; | |
139 | if (evcnt == noev) | |
140 | break; | |
141 | } | |
142 | for (Int_t i = 5; i < tokens->GetEntriesFast(); i++) { | |
143 | UInt_t trackletWord = 0; | |
144 | sscanf(((TObjString*) tokens->At(i))->GetString().Data(), "%i", &trackletWord); | |
145 | if (trackletWord == 0x10001000) | |
146 | break; | |
147 | AliDebug(2, Form("link: %2i trkl: %2i - %s -> 0x%08x", | |
148 | iLink, i-4, ((TObjString*) tokens->At(i))->GetString().Data(), trackletWord)); | |
149 | AliTRDtrackletWord *tracklet = new (trklArray[trklArray.GetEntriesFast()]) AliTRDtrackletWord(trackletWord); | |
150 | AliTRDtrackletGTU *trkl = new (trklArrayGTU[trklArrayGTU.GetEntriesFast()]) AliTRDtrackletGTU(tracklet); | |
151 | if (fTMU) | |
152 | fTMU->AddTracklet(trkl, iLink); | |
153 | } | |
154 | } | |
155 | ||
156 | if (fTMU && evcnt < noev) { | |
157 | TList *listOfTracks = new TList(); | |
158 | fTMU->SetStack(iStackPrev); | |
159 | fTMU->SetSector(iSecPrev); | |
160 | fTMU->RunTMU(listOfTracks); | |
161 | WriteTracksToDataFile(listOfTracks, iEventPrev); | |
162 | delete fTMU; | |
163 | delete listOfTracks; | |
164 | fTMU = 0x0; | |
165 | } | |
166 | ||
167 | AliInfo(Form("Analyzed %i events", evcnt)); | |
168 | return kTRUE; | |
169 | } | |
170 | ||
171 | Bool_t AliTRDgtuSim::RunGTU(AliLoader *loader, AliESDEvent *esd) | |
172 | { | |
173 | // run the GTU on tracklets taken from the loader | |
174 | // if specified the GTU tracks are written to the ESD event | |
175 | ||
176 | if (!fFeeParam->GetTracklet()) | |
177 | return kFALSE; | |
178 | ||
179 | if (fTrackletArray) | |
180 | fTrackletArray->Clear(); | |
181 | ||
182 | if (loader) { | |
183 | if (!LoadTracklets(loader)) { | |
184 | AliError("Could not load the tracklets. Nothing done ..."); | |
185 | return kFALSE; | |
186 | } | |
187 | } | |
188 | else { | |
189 | LoadTracklets(esd); | |
190 | } | |
191 | ||
192 | AliDebug(1, Form("running on %i tracklets", fTrackletArray->GetEntriesFast())); | |
193 | ||
194 | Int_t iStackPrev = -1; | |
195 | Int_t iSecPrev = -1; | |
196 | Int_t iSec = -1; | |
197 | Int_t iStack = -1; | |
198 | Int_t iLink = -1; | |
199 | ||
200 | if (fTMU) { | |
201 | delete fTMU; | |
202 | fTMU = 0x0; | |
203 | } | |
204 | ||
205 | TList *listOfTracks = new TList(); | |
206 | ||
207 | TIter next(fTrackletArray); | |
208 | ||
209 | while (AliTRDtrackletGTU *trkl = (AliTRDtrackletGTU*) next()) { | |
210 | iSec = trkl->GetDetector() / 30; | |
211 | iStack = (trkl->GetDetector() % 30) / 6; | |
212 | iLink = trkl->GetHCId() % 12; | |
213 | ||
214 | if (iStack != iStackPrev || iSec != iSecPrev) { | |
215 | if(fTMU) { | |
216 | fTMU->SetStack(iStackPrev); | |
217 | fTMU->SetSector(iSecPrev); | |
218 | fTMU->RunTMU(listOfTracks); | |
219 | WriteTracksToLoader(listOfTracks); | |
220 | WriteTracksToESD(listOfTracks, esd); | |
221 | fTMU->Reset(); | |
222 | listOfTracks->Delete(); | |
223 | } else { | |
224 | fTMU = new AliTRDgtuTMU(); | |
225 | } | |
226 | iStackPrev = iStack; | |
227 | iSecPrev = iSec; | |
228 | AliDebug(1, Form("now in sec %i, stack %i", iSec, iStack)); | |
229 | } | |
230 | AliDebug(1, Form("adding tracklet: 0x%08x in sec %i stack %i link %i", | |
231 | trkl->GetTrackletWord(), trkl->GetDetector() / 30, (trkl->GetDetector() % 30) / 6, trkl->GetHCId() % 12)); | |
232 | if (fTMU) { | |
233 | fTMU->AddTracklet(trkl, iLink); | |
234 | } | |
235 | } | |
236 | ||
237 | if (fTMU) { | |
238 | fTMU->SetStack(iStackPrev); | |
239 | fTMU->SetSector(iSecPrev); | |
240 | fTMU->RunTMU(listOfTracks); | |
241 | WriteTracksToLoader(listOfTracks); | |
242 | WriteTracksToESD(listOfTracks, esd); | |
243 | delete fTMU; | |
244 | fTMU = 0x0; | |
245 | listOfTracks->Delete(); | |
246 | } | |
247 | ||
248 | delete listOfTracks; | |
249 | ||
250 | return kTRUE; | |
251 | } | |
252 | ||
253 | Bool_t AliTRDgtuSim::LoadTracklets(const AliESDEvent *const esd) | |
254 | { | |
255 | AliDebug(1,"Loading tracklets from ESD event ..."); | |
256 | ||
257 | if (!fTrackletArray) | |
258 | fTrackletArray = new TClonesArray("AliTRDtrackletGTU", 1000); | |
259 | ||
260 | for (Int_t iTracklet = 0; iTracklet < esd->GetNumberOfTrdTracklets(); iTracklet++) { | |
261 | AliESDTrdTracklet *trkl = esd->GetTrdTracklet(iTracklet); | |
262 | new ((*fTrackletArray)[fTrackletArray->GetEntries()]) AliTRDtrackletGTU(trkl); | |
263 | } | |
264 | ||
265 | return kTRUE; | |
266 | } | |
267 | ||
268 | Bool_t AliTRDgtuSim::LoadTracklets(AliLoader *const loader) | |
269 | { | |
270 | // load the tracklets using the given loader | |
271 | ||
272 | AliDebug(1,"Loading tracklets ..."); | |
273 | ||
274 | if (!fFeeParam->GetTracklet()) | |
275 | return kFALSE; | |
276 | ||
277 | if (!loader) { | |
278 | AliError("No loader given!"); | |
279 | return kFALSE; | |
280 | } | |
281 | ||
282 | AliDataLoader *trackletLoader = loader->GetDataLoader("tracklets"); | |
283 | if (!trackletLoader) { | |
284 | AliError("No tracklet loader found!"); | |
285 | return kFALSE; | |
286 | } | |
287 | ||
288 | trackletLoader->Load(); | |
289 | TTree *trackletTree = 0x0; | |
290 | ||
291 | // simulated tracklets | |
292 | trackletTree = trackletLoader->Tree(); | |
293 | if (trackletTree) { | |
294 | TBranch *trklbranch = trackletTree->GetBranch("mcmtrklbranch"); | |
295 | if (trklbranch) { | |
296 | if (!fTrackletArray) | |
297 | fTrackletArray = new TClonesArray("AliTRDtrackletGTU", 1000); | |
298 | ||
299 | AliTRDtrackletMCM *trkl = 0x0; | |
300 | trklbranch->SetAddress(&trkl); | |
301 | for (Int_t iTracklet = 0; iTracklet < trklbranch->GetEntries(); iTracklet++) { | |
302 | trklbranch->GetEntry(iTracklet); | |
303 | new ((*fTrackletArray)[fTrackletArray->GetEntries()]) AliTRDtrackletGTU(new AliTRDtrackletMCM(*trkl)); | |
304 | } | |
305 | return kTRUE; | |
306 | } | |
307 | } | |
308 | ||
309 | // raw tracklets | |
310 | AliTreeLoader *tl = (AliTreeLoader*) trackletLoader->GetBaseLoader("tracklets-raw"); | |
311 | trackletTree = tl ? tl->Load(), tl->Tree() : 0x0; | |
312 | ||
313 | if (trackletTree) { | |
314 | if (!fTrackletArray) | |
315 | fTrackletArray = new TClonesArray("AliTRDtrackletGTU", 1000); | |
316 | ||
317 | Int_t hc; | |
318 | TClonesArray *ar = 0x0; | |
319 | trackletTree->SetBranchAddress("hc", &hc); | |
320 | trackletTree->SetBranchAddress("trkl", &ar); | |
321 | ||
322 | for (Int_t iEntry = 0; iEntry < trackletTree->GetEntries(); iEntry++) { | |
323 | trackletTree->GetEntry(iEntry); | |
324 | AliDebug(2, Form("%i tracklets in HC %i", ar->GetEntriesFast(), hc)); | |
325 | for (Int_t iTracklet = 0; iTracklet < ar->GetEntriesFast(); iTracklet++) { | |
326 | AliTRDtrackletWord *trklWord = (AliTRDtrackletWord*) (*ar)[iTracklet]; | |
327 | new((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletGTU(new AliTRDtrackletWord(trklWord->GetTrackletWord(), hc)); | |
328 | } | |
329 | } | |
330 | return kTRUE; | |
331 | } | |
332 | ||
333 | AliError("No raw tracklet tree found\n"); | |
334 | ||
335 | return kFALSE; | |
336 | } | |
337 | ||
338 | Bool_t AliTRDgtuSim::WriteTracksToDataFile(TList *listOfTracks, Int_t event) | |
339 | { | |
340 | // write the found tracks to a data file | |
341 | // used for comparison to VHDL simulation | |
342 | ||
343 | Int_t sm = 0; | |
344 | Int_t stack = 0; | |
345 | ||
346 | FILE *out; | |
347 | out = fopen("test.data", "a"); | |
348 | ||
349 | AliDebug(1,Form("%i tracks found in event %i", listOfTracks->GetSize(), event)); | |
350 | // fprintf(out, "0 %5i %2i %i 00000000\n", event, sm, stack); | |
351 | for (Int_t i = 0; i < listOfTracks->GetSize(); i++) { | |
352 | AliTRDtrackGTU *trk = (AliTRDtrackGTU*) listOfTracks->At(i); | |
353 | sm = trk->GetSector(); | |
354 | stack = trk->GetStack(); | |
355 | ||
356 | ULong64_t trackWord = 1; | |
357 | AppendBits(trackWord, 1, 0); | |
358 | AppendBits(trackWord, 6, trk->GetTrackletMask()); | |
359 | AppendBits(trackWord, 18, (Int_t) trk->GetA()); | |
360 | AppendBits(trackWord, 18, (Int_t) trk->GetB()); | |
361 | AppendBits(trackWord, 12, (Int_t) trk->GetC()); | |
362 | AppendBits(trackWord, 8, trk->GetPID()); | |
363 | fprintf(out, "ev. %i sec. %i stack %i - track word: 0x%016llx, ", | |
364 | event, sm, stack, trackWord); | |
365 | ||
366 | trackWord = 0; | |
367 | AppendBits(trackWord, 11, 0); // flags | |
368 | AppendBits(trackWord, 3, 0); | |
369 | AppendBits(trackWord, 13, trk->GetYapprox()); | |
370 | AppendBits(trackWord, 6, trk->GetTrackletIndex(5)); | |
371 | AppendBits(trackWord, 6, trk->GetTrackletIndex(4)); | |
372 | AppendBits(trackWord, 6, trk->GetTrackletIndex(3)); | |
373 | AppendBits(trackWord, 6, trk->GetTrackletIndex(2)); | |
374 | AppendBits(trackWord, 6, trk->GetTrackletIndex(1)); | |
375 | AppendBits(trackWord, 6, trk->GetTrackletIndex(0)); | |
376 | fprintf(out, "extended track word: 0x%016llx\n", trackWord); | |
377 | ||
378 | fprintf(out, "1 %5i %2i %2i %3i %3i %3i %3i %3i %3i %3i %4i %f\n", event, sm, stack, trk->GetTrackletMask(), | |
379 | trk->GetTrackletIndex(5), | |
380 | trk->GetTrackletIndex(4), | |
381 | trk->GetTrackletIndex(3), | |
382 | trk->GetTrackletIndex(2), | |
383 | trk->GetTrackletIndex(1), | |
384 | trk->GetTrackletIndex(0), | |
385 | trk->GetPtInt(), | |
386 | trk->GetPt()); | |
387 | } | |
388 | fclose(out); | |
389 | return kTRUE; | |
390 | } | |
391 | ||
392 | Bool_t AliTRDgtuSim::WriteTracksToESD(const TList * const listOfTracks, AliESDEvent *esd) | |
393 | { | |
394 | // fill the found tracks to the given ESD event | |
395 | ||
396 | if (esd) { | |
397 | TIter next(listOfTracks); | |
398 | while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) next()) { | |
399 | AliESDTrdTrack *trdtrack = trk->CreateTrdTrack(); | |
400 | if (trdtrack->GetLabel() < 0) | |
401 | trdtrack->SetLabel(-2); | |
402 | esd->AddTrdTrack(trdtrack); | |
403 | delete trdtrack; | |
404 | } | |
405 | } | |
406 | return kTRUE; | |
407 | } | |
408 | ||
409 | Bool_t AliTRDgtuSim::WriteTracksToLoader(const TList * const listOfTracks) | |
410 | { | |
411 | // write the GTU tracks to the dedicated loader | |
412 | // these tracks contain more information than the ones in the ESD | |
413 | ||
414 | AliRunLoader *rl = AliRunLoader::Instance(); | |
415 | AliDataLoader *dl = 0x0; | |
416 | if (rl) | |
417 | dl = rl->GetLoader("TRDLoader")->GetDataLoader("gtutracks"); | |
418 | if (!dl) { | |
419 | AliError("Could not get the GTU-track data loader!"); | |
420 | return kFALSE; | |
421 | } | |
422 | ||
423 | TTree *trackTree = dl->Tree(); | |
424 | if (!trackTree) { | |
425 | dl->MakeTree(); | |
426 | trackTree = dl->Tree(); | |
427 | } | |
428 | ||
429 | AliTRDtrackGTU *trk = 0x0; | |
430 | ||
431 | if (!trackTree->GetBranch("TRDtrackGTU")) | |
432 | trackTree->Branch("TRDtrackGTU", "AliTRDtrackGTU", &trk, 32000); | |
433 | ||
434 | AliDebug(1, Form("Writing %i tracks to loader", listOfTracks->GetEntries())); | |
435 | TIter next(listOfTracks); | |
436 | while ((trk = (AliTRDtrackGTU*) next())) { | |
437 | trackTree->SetBranchAddress("TRDtrackGTU", &trk); | |
438 | trackTree->Fill(); | |
439 | } | |
440 | ||
441 | dl->WriteData("OVERWRITE"); | |
442 | ||
443 | return kTRUE; | |
444 | } |