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