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