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