1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /* $Id: AliTRDgtuSim.cxx 28397 2008-09-02 09:33:00Z cblume $ */
18 ////////////////////////////////////////////////////////////////////////////
22 // Authors: J. Klein (Jochen.Klein@cern.ch) //
24 ////////////////////////////////////////////////////////////////////////////
32 #include "TClonesArray.h"
35 #include "AliRunLoader.h"
36 #include "AliLoader.h"
37 #include "AliTreeLoader.h"
39 #include "AliESDTrdTrack.h"
41 #include "AliTRDgtuSim.h"
42 #include "AliTRDgtuTMU.h"
43 #include "AliTRDtrackGTU.h"
44 #include "AliTRDtrackletWord.h"
45 #include "AliTRDtrackletMCM.h"
46 #include "AliESDEvent.h"
48 ClassImp(AliTRDgtuSim)
50 AliTRDgtuSim::AliTRDgtuSim(AliRunLoader *rl)
58 fTrackletTree = new TTree("gtutracklets", "Tree with GTU tracklets");
59 fTrackletTree->SetDirectory(0);
62 AliTRDgtuSim::~AliTRDgtuSim()
67 fTrackletArray->Delete();
68 delete fTrackletArray;
72 Bool_t AliTRDgtuSim::RunGTUFromTrackletFile(TString filename, Int_t event, Int_t noev)
74 // run the GTU from a file of tracklets
75 // used for comparison to VHDL simulation
77 AliInfo(Form("Running the GTU simulation on file: %s", filename.Data()));
78 ifstream input(filename.Data());
84 Int_t iEventPrev = -1;
85 Int_t iStackPrev = -1;
95 AliDebug(5,"--------- Reading from file ----------");
96 while (getline(input, str)) {
99 AliDebug(5,Form("Line %i : %s", lineno, string.Data()));
101 TObjArray *tokens = string.Tokenize(" ");
102 if (tokens->GetEntriesFast() < 7) {
103 AliWarning(Form("Invalid input in line %i, too few parameters", lineno));
107 if ( ((TObjString*) tokens->At(0))->GetString().Atoi() < event)
110 iEvent = ((TObjString*) tokens->At(0))->GetString().Atoi();
111 iSec = ((TObjString*) tokens->At(1))->GetString().Atoi();
112 iStack = ((TObjString*) tokens->At(2))->GetString().Atoi();
113 iLink = 2 * ((TObjString*) tokens->At(3))->GetString().Atoi() + ((TObjString*) tokens->At(4))->GetString().Atoi();
115 if (iEvent != iEventPrev || iStack != iStackPrev || iSec != iSecPrev) {
117 TList *listOfTracks = new TList();
118 fTMU->SetStack(iStackPrev);
119 fTMU->SetSector(iSecPrev);
120 fTMU->RunTMU(listOfTracks);
121 AliDebug(1,Form("--- There are %i tracks. Writing ...", listOfTracks->GetEntries()));
122 WriteTracksToTree(listOfTracks);
123 fTMU->WriteTrackletsToTree(fTrackletTree);
124 WriteTracksToDataFile(listOfTracks, iEventPrev);
125 if (listOfTracks->GetEntries() > 0)
126 AliDebug(2,Form(" %d GeV/c", ((AliTRDtrackGTU*) listOfTracks->At(0))->GetPt() ));
128 fTMU = new AliTRDgtuTMU();
132 fTMU = new AliTRDgtuTMU();
141 for (Int_t i = 5; i < tokens->GetEntriesFast(); i++) {
142 UInt_t trackletWord = 0;
143 sscanf(((TObjString*) tokens->At(i))->GetString().Data(), "%x", &trackletWord);
144 if (trackletWord == 0x10001000)
146 AliDebug(2,Form("%i. tracklet: %s -> 0x%08x", i-4, ((TObjString*) tokens->At(i))->GetString().Data(), trackletWord));
147 AliTRDtrackletWord *trkl = new AliTRDtrackletWord(trackletWord);
148 fTMU->AddTracklet(trkl, iLink);
152 if (fTMU && evcnt < noev) {
153 TList *listOfTracks = new TList();
154 fTMU->SetStack(iStackPrev);
155 fTMU->SetSector(iSecPrev);
156 fTMU->RunTMU(listOfTracks);
157 WriteTracksToTree(listOfTracks);
158 fTMU->WriteTrackletsToTree(fTrackletTree);
159 WriteTracksToDataFile(listOfTracks, iEventPrev);
165 AliInfo(Form("Analyzed %i events", evcnt));
169 Bool_t AliTRDgtuSim::RunGTU(AliLoader *loader, AliESDEvent *esd)
171 // run the GTU on tracklets taken from the loader
172 // if specified the GTU tracks are written to the ESD event
174 if (!LoadTracklets(loader)) {
175 AliError("Could not load the tracklets. Nothing done ...");
179 AliDebug(1, Form("running on %i tracklets", fTrackletArray->GetEntriesFast()));
181 Int_t iStackPrev = -1;
192 TList *listOfTracks = new TList();
194 TIter next(fTrackletArray);
195 AliTRDtrackletBase *trkl;
197 while ((trkl = (AliTRDtrackletBase*) next())) {
198 iSec = trkl->GetDetector() / 30;
199 iStack = (trkl->GetDetector() % 30) / 6;
200 iLink = 2 * (trkl->GetDetector() % 6) + (trkl->GetYbin() < 0 ? 0 : 1);
202 if (iStack != iStackPrev || iSec != iSecPrev) {
204 fTMU->SetStack(iStackPrev);
205 fTMU->SetSector(iSecPrev);
206 fTMU->RunTMU(listOfTracks);
207 WriteTracksToTree(listOfTracks);
208 fTMU->WriteTrackletsToTree(fTrackletTree);
209 WriteTracksToESD(listOfTracks, esd);
211 listOfTracks->Delete();
213 fTMU = new AliTRDgtuTMU();
218 AliDebug(1, Form("adding tracklet: 0x%08x", trkl->GetTrackletWord()));
219 fTMU->AddTracklet(trkl, iLink);
223 fTMU->SetStack(iStackPrev);
224 fTMU->SetSector(iSecPrev);
225 fTMU->RunTMU(listOfTracks);
226 WriteTracksToTree(listOfTracks);
227 fTMU->WriteTrackletsToTree(fTrackletTree);
228 WriteTracksToESD(listOfTracks, esd);
231 listOfTracks->Delete();
239 Bool_t AliTRDgtuSim::LoadTracklets(AliLoader *const loader)
241 // load the tracklets using the given loader
243 AliDebug(1,"Loading tracklets ...");
246 AliError("No loader given!");
250 AliDataLoader *trackletLoader = loader->GetDataLoader("tracklets");
251 if (!trackletLoader) {
252 AliError("No tracklet loader found!");
256 trackletLoader->Load();
257 TTree *trackletTree = 0x0;
259 // simulated tracklets
260 trackletTree = trackletLoader->Tree();
262 TBranch *trklbranch = trackletTree->GetBranch("mcmtrklbranch");
265 fTrackletArray = new TClonesArray("AliTRDtrackletMCM", 1000);
266 else if ((TClass::GetClass("AliTRDtrackletMCM"))->InheritsFrom(fTrackletArray->Class()))
267 fTrackletArray->Delete();
269 fTrackletArray->Delete();
270 delete fTrackletArray;
271 fTrackletArray = new TClonesArray("AliTRDtrackletMCM", 1000);
274 AliTRDtrackletMCM *trkl = new AliTRDtrackletMCM;
275 trklbranch->SetAddress(&trkl);
276 for (Int_t iTracklet = 0; iTracklet < trklbranch->GetEntries(); iTracklet++) {
277 trklbranch->GetEntry(iTracklet);
278 new ((*fTrackletArray)[fTrackletArray->GetEntries()]) AliTRDtrackletMCM(*trkl);
285 AliTreeLoader *tl = (AliTreeLoader*) trackletLoader->GetBaseLoader("tracklets-raw");
286 trackletTree = tl ? tl->Load(), tl->Tree() : 0x0;
290 fTrackletArray = new TClonesArray("AliTRDtrackletWord", 1000);
291 else if ((TClass::GetClass("AliTRDtrackletWord"))->InheritsFrom(fTrackletArray->Class()))
292 fTrackletArray->Delete();
294 fTrackletArray->Delete();
295 delete fTrackletArray;
296 fTrackletArray = new TClonesArray("AliTRDtrackletWord", 1000);
300 TClonesArray *ar = 0x0;
301 trackletTree->SetBranchAddress("hc", &hc);
302 trackletTree->SetBranchAddress("trkl", &ar);
304 for (Int_t iEntry = 0; iEntry < trackletTree->GetEntries(); iEntry++) {
305 trackletTree->GetEntry(iEntry);
306 printf("%i tracklets in HC %i\n", ar->GetEntriesFast(), hc);
307 for (Int_t iTracklet = 0; iTracklet < ar->GetEntriesFast(); iTracklet++) {
308 AliTRDtrackletWord *trklWord = (AliTRDtrackletWord*) (*ar)[iTracklet];
309 new((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(trklWord->GetTrackletWord(), hc);
315 AliError("No raw tracklet tree found\n");
320 Bool_t AliTRDgtuSim::WriteTracksToDataFile(TList *listOfTracks, Int_t event)
322 // write the found tracks to a data file
323 // used for comparison to VHDL simulation
329 out = fopen("test.data", "a");
331 AliDebug(1,Form("%i tracks found in event %i", listOfTracks->GetSize(), event));
332 fprintf(out, "0 %5i %2i %i 00000000\n", event, sm, stack);
333 for (Int_t i = 0; i < listOfTracks->GetSize(); i++) {
334 AliTRDtrackGTU *trk = (AliTRDtrackGTU*) listOfTracks->At(i);
335 sm = trk->GetSector();
336 stack = trk->GetStack();
337 fprintf(out, "1 %5i %2i %2i %3i %3i %3i %3i %3i %3i %3i %4i %f\n", event, sm, stack, trk->GetTrackletMask(),
338 trk->GetTrackletIndex(5),
339 trk->GetTrackletIndex(4),
340 trk->GetTrackletIndex(3),
341 trk->GetTrackletIndex(2),
342 trk->GetTrackletIndex(1),
343 trk->GetTrackletIndex(0),
351 Bool_t AliTRDgtuSim::WriteTracksToTree(TList *listOfTracks, Int_t /*event*/)
353 // write the tracks to the tree for intermediate storage
355 AliDebug(1,Form("Writing %i tracks to the tree...", listOfTracks->GetEntries()));
360 if (listOfTracks->GetEntries() <= 0)
364 fTrackTree = new TTree("gtutracks", "GTU tracks");
365 fTrackTree->SetDirectory(0);
368 AliTRDtrackGTU *trk = 0x0;
369 TBranch *branch = fTrackTree->GetBranch("TRDgtuTrack");
371 branch = fTrackTree->Branch("TRDgtuTrack", "AliTRDtrackGTU", &trk, 32000, 99);
374 TIter next(listOfTracks);
375 while ((trk = (AliTRDtrackGTU*) next())) {
377 branch->SetAddress(&trk);
380 fTrackTree->ResetBranchAddress(branch);
385 Bool_t AliTRDgtuSim::WriteTreesToFile() const {
386 // write the trees holding tracklets and tracks to file
388 TFile *f = TFile::Open("TRD.GtuTracking.root", "RECREATE");
391 f->WriteTObject(fTrackTree);
393 f->WriteTObject(fTrackletTree);
398 Bool_t AliTRDgtuSim::WriteTracksToESD(const TList * const listOfTracks, AliESDEvent *esd)
400 // fill the found tracks to the given ESD event
403 TIter next(listOfTracks);
404 while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) next()) {
405 AliESDTrdTrack *trdtrack = trk->CreateTrdTrack();
406 esd->AddTrdTrack(trdtrack);
413 Bool_t AliTRDgtuSim::WriteTracksToLoader()
415 // write the GTU tracks to the dedicated loader
416 // these tracks contain more information than the ones in the ESD
419 AliError("No track tree found!");
423 AliRunLoader *rl = AliRunLoader::Instance();
424 AliDataLoader *dl = 0x0;
426 dl = rl->GetLoader("TRDLoader")->GetDataLoader("gtutracks");
428 AliError("Could not get the GTU-track data loader!");
432 TTree *trackTree = dl->Tree();
435 trackTree = dl->Tree();
438 AliTRDtrackGTU *trk = 0x0;
439 if (!trackTree->GetBranch("TRDtrackGTU"))
440 trackTree->Branch("TRDtrackGTU", "AliTRDtrackGTU", &trk, 32000);
442 AliDebug(1,Form("Found %i tracks", fTrackTree->GetEntries()));
444 for (Int_t iTrack = 0; iTrack < fTrackTree->GetEntries(); iTrack++) {
445 fTrackTree->SetBranchAddress("TRDgtuTrack", &trk);
446 fTrackTree->GetEntry(iTrack);
447 trackTree->SetBranchAddress("TRDtrackGTU", &trk);
451 dl->WriteData("OVERWRITE");