]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRunDigitizer.cxx
Rotation of modules in x-y plane corrected.
[u/mrichter/AliRoot.git] / STEER / AliRunDigitizer.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 /*
17 $Log$
18 Revision 1.20  2002/10/14 14:57:32  hristov
19 Merging the VirtualMC branch to the main development branch (HEAD)
20
21 Revision 1.13.6.2  2002/07/24 10:08:13  alibrary
22 Updating VirtualMC
23
24 Revision 1.19  2002/07/19 12:46:05  hristov
25 Write file instead of closing it
26
27 Revision 1.18  2002/07/17 08:59:39  jchudoba
28 Do not delete subtasks when AliRunDigitizer is deleted. Owner should delete them itself.
29
30 Revision 1.17  2002/07/16 13:47:53  jchudoba
31 Add methods to get access to names of files used in merging.
32
33 Revision 1.16  2002/06/07 09:18:47  jchudoba
34 Changes to enable merging of ITS fast rec points. Although this class should be responsible for a creation of digits only, other solutions would be more complicated.
35
36 Revision 1.15  2002/04/09 13:38:47  jchudoba
37 Add const to the filename argument
38
39 Revision 1.14  2002/04/04 09:28:04  jchudoba
40 Change default names of TPC trees. Use update instead of recreate for the output file. Overwrite the AliRunDigitizer object in the output if it exists.
41
42 Revision 1.13  2002/02/13 09:03:32  jchudoba
43 Pass option to subtasks. Delete input TTrees. Use gAlice from memory if it is present (user must delete the default one created by aliroot if he/she wants to use gAlice from the input file!). Add new data member to store name of the special TPC TTrees.
44
45 Revision 1.12  2001/12/10 16:40:52  jchudoba
46 Import gAlice from the signal file before InitGlobal() to allow detectors to use it during initialization
47
48 Revision 1.11  2001/12/03 07:10:13  jchudoba
49 Default ctor cannot create new objects, create dummy default ctor which leaves object in not well defined state - to be used only by root for I/O
50
51 Revision 1.10  2001/11/15 11:07:25  jchudoba
52 Set to zero new pointers to TPC and TRD special trees in the default ctor. Add const to all Get functions. Remove unused constant, rename constant according coding rules.
53
54 Revision 1.9  2001/11/15 09:00:11  jchudoba
55 Add special treatment for TPC and TRD, they use different trees than other detectors
56
57 Revision 1.8  2001/10/21 18:38:43  hristov
58 Several pointers were set to zero in the default constructors to avoid memory management problems
59
60 Revision 1.7  2001/10/04 15:56:07  jchudoba
61 TTask inheritance
62
63 Revision 1.4  2001/09/19 06:23:50  jchudoba
64 Move some tasks to AliStream and AliMergeCombi classes
65
66 Revision 1.3  2001/07/30 14:04:18  jchudoba
67 correct bug in the initialization
68
69 Revision 1.2  2001/07/28 10:44:32  hristov
70 Loop variable declared once; typos corrected
71
72 Revision 1.1  2001/07/27 12:59:00  jchudoba
73 Manager class for merging/digitization
74
75 */
76
77 ////////////////////////////////////////////////////////////////////////
78 //
79 // AliRunDigitizer.cxx
80 //
81 // Manager object for merging/digitization
82 //
83 // Instance of this class manages the digitization and/or merging of
84 // Sdigits into Digits. 
85 //
86 // Only one instance of this class is created in the macro:
87 //   AliRunDigitizer * manager = 
88 //      new AliRunDigitizer(nInputStreams,SPERB);
89 // where nInputStreams is number of input streams and SPERB is
90 // signals per background variable, which determines how combinations
91 // of signal and background events are generated.
92 // Then instances of specific detector digitizers are created:
93 //   AliMUONDigitizer *dMUON  = new AliMUONDigitizer(manager)
94 // and the I/O configured (you have to specify input files 
95 // and an output file). The manager connects appropriate trees from 
96 // the input files according a combination returned by AliMergeCombi 
97 // class. It creates TreeD in the output and runs once per 
98 // event Digitize method of all existing AliDetDigitizers 
99 // (without any option). AliDetDigitizers ask manager
100 // for a TTree with input (manager->GetInputTreeS(Int_t i),
101 // merge all inputs, digitize it, and save it in the TreeD 
102 // obtained by manager->GetTreeD(). Output events are stored with 
103 // numbers from 0, this default can be changed by 
104 // manager->SetFirstOutputEventNr(Int_t) method. The particle numbers
105 // in the output are shifted by MASK, which is taken from manager.
106 //
107 // The default output is to the signal file (stream 0). This can be 
108 // changed with the SetOutputFile(TString fn)  method.
109 //
110 // Single input file is permitted. Maximum kMaxStreamsToMerge can be merged.
111 // Input from the memory (on-the-fly merging) is not yet 
112 // supported, as well as access to the input data by invoking methods
113 // on the output data.
114 //
115 // Access to the some data is via gAlice for now (supposing the 
116 // same geometry in all input files), gAlice is taken from the first 
117 // input file on the first stream.
118 //
119 // Example with MUON digitizer, no merging, just digitization
120 //
121 //  AliRunDigitizer * manager = new AliRunDigitizer(1,1);
122 //  manager->SetInputStream(0,"galice.root");
123 //  AliMUONDigitizer *dMUON  = new AliMUONDigitizer(manager);
124 //  manager->Exec("");
125 //
126 // Example with MUON digitizer, merge all events from 
127 //   galice.root (signal) file with events from bgr.root 
128 //   (background) file. Number of merged events is
129 //   min(number of events in galice.root, number of events in bgr.root)
130 //
131 //  AliRunDigitizer * manager = new AliRunDigitizer(2,1);
132 //  manager->SetInputStream(0,"galice.root");
133 //  manager->SetInputStream(1,"bgr.root");
134 //  AliMUONDigitizer *dMUON  = new AliMUONDigitizer(manager);
135 //  manager->Exec("");
136 //
137 // Example with MUON digitizer, save digits in a new file digits.root,
138 //   process only 1 event
139 //
140 //  AliRunDigitizer * manager = new AliRunDigitizer(2,1);
141 //  manager->SetInputStream(0,"galice.root");
142 //  manager->SetInputStream(1,"bgr.root");
143 //  manager->SetOutputFile("digits.root");
144 //  AliMUONDigitizer *dMUON  = new AliMUONDigitizer(manager);
145 //  manager->SetNrOfEventsToWrite(1);
146 //  manager->Exec("");
147 //
148 //////////////////////////////////////////////////////////////////////// 
149
150 // system includes
151
152 #include <Riostream.h>
153
154 // ROOT includes
155
156 #include "TFile.h"
157 #include "TTree.h"
158 #include "TList.h"
159
160 // AliROOT includes
161
162 #include "AliRunDigitizer.h"
163 #include "AliDigitizer.h"
164 #include "AliRun.h"
165 #include "AliHeader.h"
166 #include "TParticle.h"
167 #include "AliStream.h"
168 #include "AliMergeCombi.h"
169
170 ClassImp(AliRunDigitizer)
171
172 ////////////////////////////////////////////////////////////////////////
173 AliRunDigitizer::AliRunDigitizer()
174 {
175 // root requires default ctor, where no new objects can be created
176 // do not use this ctor, it is supplied only for root needs
177   
178 // just set all pointers - data members to 0
179   fOutput = 0;
180   fTreeD = 0;
181   fTreeR = 0;
182   fTreeDTPC = 0;
183   fTreeDTRD = 0;
184   fInputStreams = 0;
185   for (Int_t i=0;i<kMaxStreamsToMerge;i++) {
186     fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
187     fInputFiles[i]=0;
188   }
189   fCombi = 0;
190
191 }
192
193 ////////////////////////////////////////////////////////////////////////
194 AliRunDigitizer::AliRunDigitizer(Int_t nInputStreams, Int_t sperb) : TTask("AliRunDigitizer","The manager for Merging")
195 {
196 // ctor which should be used to create a manager for merging/digitization
197   if (nInputStreams == 0) {
198     Error("AliRunDigitizer","Specify nr of input streams");
199     return;
200   }
201   Int_t i;
202   fNinputs = nInputStreams;
203   fOutputFileName = "";
204   fOutputDirName = ".";
205   fCombination.Set(kMaxStreamsToMerge);
206   for (i=0;i<kMaxStreamsToMerge;i++) {
207     fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
208     fCombination[i]=-1;
209   }
210   fkMASKSTEP = 10000000;
211   fkMASK[0] = 0;
212   for (i=1;i<kMaxStreamsToMerge;i++) {
213     fkMASK[i] = fkMASK[i-1] + fkMASKSTEP;
214   }
215   fInputStreams = new TClonesArray("AliStream",nInputStreams);
216   TClonesArray &lInputStreams = *fInputStreams;
217 // the first Input is open RW to be output as well
218   new(lInputStreams[0]) AliStream("UPDATE");
219   for (i=1;i<nInputStreams;i++) {
220     new(lInputStreams[i]) AliStream("READ");
221   }
222   fOutput = 0;
223   fEvent = 0;
224   fNrOfEventsToWrite = -1;
225   fNrOfEventsWritten = 0;
226   fCopyTreesFromInput = -1;
227   fCombi = new AliMergeCombi(nInputStreams,sperb);
228   fDebug = 0;
229   fTreeD = 0;
230   fTreeR = 0;
231   fTreeDTPC = 0;
232   fTreeDTRD = 0;
233   fTreeDTPCBaseName = "TreeD_75x40_100x60_150x60_";
234   fTreeTPCSBaseName = "TreeS_75x40_100x60_150x60_";
235
236   for (i=0; i<kMaxStreamsToMerge; i++) fInputFiles[i]=0;
237 }
238
239 ////////////////////////////////////////////////////////////////////////
240
241 AliRunDigitizer::~AliRunDigitizer() {
242 // dtor
243
244 // do not delete subtasks, let the creator delete them
245   if (GetListOfTasks()) 
246     GetListOfTasks()->Clear("nodelete");
247   
248   if (fInputStreams) {
249     delete fInputStreams;
250     fInputStreams = 0;
251   }
252   if (fCombi) {
253     delete fCombi;
254     fCombi = 0;
255   }
256
257 }
258 ////////////////////////////////////////////////////////////////////////
259 void AliRunDigitizer::AddDigitizer(AliDigitizer *digitizer)
260 {
261 // add digitizer to the list of active digitizers
262   this->Add(digitizer);
263 }
264 ////////////////////////////////////////////////////////////////////////
265 void AliRunDigitizer::SetInputStream(Int_t i, const char *inputFile)
266 {
267   if (i > fInputStreams->GetLast()) {
268     Error("SetInputStream","Input stream number too high");
269     return;
270   }
271   static_cast<AliStream*>(fInputStreams->At(i))->AddFile(inputFile);
272 }
273
274 ////////////////////////////////////////////////////////////////////////
275 void AliRunDigitizer::Digitize(Option_t* option)
276 {
277 // get a new combination of inputs, connect input trees and loop 
278 // over all digitizers
279
280 // take gAlice from the first input file. It is needed to access
281 //  geometry data
282 // If gAlice is already in memory, use it
283   if (!gAlice) {
284     if (!static_cast<AliStream*>(fInputStreams->At(0))->ImportgAlice()) {
285       cerr<<"gAlice object not found in the first file of "
286           <<"the 1st stream"<<endl;
287       return;
288     }
289   }
290   if (!InitGlobal()) {
291     cerr<<"False from InitGlobal"<<endl;
292     return;
293   }
294   Int_t eventsCreated = 0;
295 // loop until there is anything on the input in case fNrOfEventsToWrite < 0
296   while ((eventsCreated++ < fNrOfEventsToWrite) || (fNrOfEventsToWrite < 0)) {
297     if (!ConnectInputTrees()) break;
298     InitEvent();
299 // loop over all registered digitizers and let them do the work
300     ExecuteTasks(option);
301     CleanTasks();
302     FinishEvent();
303   }
304   FinishGlobal();
305 }
306
307 ////////////////////////////////////////////////////////////////////////
308 Bool_t AliRunDigitizer::ConnectInputTrees()
309 {
310 // fill arrays fArrayTreeS, fArrayTreeH and fArrayTreeTPCS with 
311 // pointers to the correct events according fCombination values
312 // null pointers can be in the output, AliDigitizer has to check it
313
314   TTree *tree;
315   char treeName[50];
316   Int_t serialNr;
317   Int_t eventNr[kMaxStreamsToMerge], delta[kMaxStreamsToMerge];
318   fCombi->Combination(eventNr, delta);
319   for (Int_t i=0;i<fNinputs;i++) {
320     if (delta[i] == 1) {
321       AliStream *iStream = static_cast<AliStream*>(fInputStreams->At(i));
322       if (!iStream->NextEventInStream(serialNr)) return kFALSE;
323       fInputFiles[i]=iStream->CurrentFile();
324       sprintf(treeName,"TreeS%d",serialNr);
325       tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
326       if (fArrayTreeS[i]) {
327         delete fArrayTreeS[i];
328         fArrayTreeS[i] = 0;
329       }
330       fArrayTreeS[i] = tree;
331       sprintf(treeName,"TreeH%d",serialNr);
332       tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
333       if (fArrayTreeH[i]) {
334         delete fArrayTreeH[i];
335         fArrayTreeH[i] = 0;
336       }
337       fArrayTreeH[i] = tree;
338       sprintf(treeName,"%s%d",fTreeTPCSBaseName,serialNr);
339       tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
340       if (fArrayTreeTPCS[i]) {
341         delete fArrayTreeTPCS[i];
342         fArrayTreeTPCS[i] = 0;
343       }
344       fArrayTreeTPCS[i] = tree;
345       sprintf(treeName,"TreeS%d_TRD",serialNr);
346       tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
347       if (fArrayTreeTRDS[i]) {
348         delete fArrayTreeTRDS[i];
349         fArrayTreeTRDS[i] = 0;
350       }
351       fArrayTreeTRDS[i] = tree;
352     } else if (delta[i] != 0) {
353       Error("ConnectInputTrees","Only delta 0 or 1 is implemented");
354       return kFALSE;
355     }
356   }
357   return kTRUE;
358 }
359
360 ////////////////////////////////////////////////////////////////////////
361 Bool_t AliRunDigitizer::InitGlobal()
362 {
363 // called once before Digitize() is called, initialize digitizers and output
364
365   TList* subTasks = this->GetListOfTasks();
366   if (subTasks) {
367     subTasks->ForEach(AliDigitizer,Init)();
368   }  
369   return kTRUE;
370 }
371
372 ////////////////////////////////////////////////////////////////////////
373
374 void AliRunDigitizer::SetOutputFile(TString fn)
375 // the output will be to separate file, not to the signal file
376 {
377   fOutputFileName = fn;
378   (static_cast<AliStream*>(fInputStreams->At(0)))->ChangeMode("READ");
379   InitOutputGlobal();
380 }
381
382 ////////////////////////////////////////////////////////////////////////
383 Bool_t AliRunDigitizer::InitOutputGlobal()
384 {
385 // Creates the output file, called by InitEvent()
386
387   TString fn;
388   fn = fOutputDirName + '/' + fOutputFileName;
389   fOutput = new TFile(fn,"update");
390   if (GetDebug()>2) {
391     cerr<<"AliRunDigitizer::InitOutputGlobal(): file "<<fn.Data()<<" was opened"<<endl;
392   }
393   if (fOutput) return kTRUE;
394   Error("InitOutputGlobal","Could not create output file.");
395   return kFALSE;
396 }
397
398
399 ////////////////////////////////////////////////////////////////////////
400 void AliRunDigitizer::InitEvent()
401 {
402 // Creates TreeDxx in the output file, called from Digitize() once for 
403 //  each event. xx = fEvent
404
405   if (GetDebug()>2) 
406     cerr<<"AliRunDigitizer::InitEvent: fEvent = "<<fEvent<<endl;
407
408 // if fOutputFileName was not given, write output to signal file
409   if (fOutputFileName == "") {
410     fOutput = (static_cast<AliStream*>(fInputStreams->At(0)))->CurrentFile();
411   }
412   fOutput->cd();
413   char treeName[30];
414   sprintf(treeName,"TreeD%d",fEvent);
415   fTreeD = static_cast<TTree*>(fOutput->Get(treeName));
416   if (!fTreeD) {
417     fTreeD = new TTree(treeName,"Digits");
418     fTreeD->Write(0,TObject::kOverwrite);
419   }
420
421 // tree for ITS fast points
422   sprintf(treeName,"TreeR%d",fEvent);
423   fTreeR = static_cast<TTree*>(fOutput->Get(treeName));
424   if (!fTreeR) {
425     fTreeR = new TTree(treeName,"Reconstruction");
426     fTreeR->Write(0,TObject::kOverwrite);
427   }
428
429 // special tree for TPC
430   sprintf(treeName,"%s%d",fTreeDTPCBaseName,fEvent);
431   fTreeDTPC = static_cast<TTree*>(fOutput->Get(treeName));
432   if (!fTreeDTPC) {
433     fTreeDTPC = new TTree(treeName,"TPC_Digits");
434     fTreeDTPC->Write(0,TObject::kOverwrite);
435   }
436
437 // special tree for TRD
438   sprintf(treeName,"TreeD%d_TRD",fEvent);
439   fTreeDTRD = static_cast<TTree*>(fOutput->Get(treeName));
440   if (!fTreeDTRD) {
441     fTreeDTRD = new TTree(treeName,"TRD_Digits");
442     fTreeDTRD->Write(0,TObject::kOverwrite);
443   }
444
445 }
446
447 ////////////////////////////////////////////////////////////////////////
448 void AliRunDigitizer::FinishEvent()
449 {
450 // called at the end of loop over digitizers
451
452   Int_t i;
453   fOutput->cd();
454   if (fCopyTreesFromInput > -1) {
455     char treeName[20];
456     i = fCopyTreesFromInput; 
457     sprintf(treeName,"TreeK%d",fCombination[i]);
458     fInputFiles[i]->Get(treeName)->Clone()->Write();
459     sprintf(treeName,"TreeH%d",fCombination[i]);
460     fInputFiles[i]->Get(treeName)->Clone()->Write();
461   }
462   fEvent++;
463   fNrOfEventsWritten++;
464   if (fTreeD) {
465     delete fTreeD;
466     fTreeD = 0;
467   }
468   if (fTreeR) {
469     delete fTreeR;
470     fTreeR = 0;
471   }
472   if (fTreeDTPC) {
473     delete fTreeDTPC;
474     fTreeDTPC = 0;
475   }
476   if (fTreeDTRD) {
477     delete fTreeDTRD;
478     fTreeDTRD = 0;
479   }
480 }
481 ////////////////////////////////////////////////////////////////////////
482 void AliRunDigitizer::FinishGlobal()
483 {
484 // called at the end of Exec
485 // save unique objects to the output file
486
487   fOutput->cd();
488   this->Write(0,TObject::kOverwrite);
489   if (fCopyTreesFromInput > -1) {
490     fInputFiles[fCopyTreesFromInput]->Get("TE")->Clone()->Write();
491     gAlice->Write();
492   }
493   fOutput->Write();
494 }
495
496
497 ////////////////////////////////////////////////////////////////////////
498 Int_t  AliRunDigitizer::GetNParticles(Int_t event) const
499 {
500 // return number of particles in all input files for a given
501 // event (as numbered in the output file)
502 // return -1 if some file cannot be accessed
503
504   Int_t sum = 0;
505   Int_t sumI;
506   for (Int_t i = 0; i < fNinputs; i++) {
507     sumI = GetNParticles(GetInputEventNumber(event,i), i);
508     if (sumI < 0) return -1;
509     sum += sumI;
510   }
511   return sum;
512 }
513
514 ////////////////////////////////////////////////////////////////////////
515 Int_t  AliRunDigitizer::GetNParticles(Int_t event, Int_t input) const
516 {
517 // return number of particles in input file input for a given
518 // event (as numbered in this input file)
519 // return -1 if some error
520
521 // Must be revised in the version with AliStream
522
523   return -1;
524
525 /*
526   TFile *file = ConnectInputFile(input);
527   if (!file) {
528     Error("GetNParticles","Cannot open input file");
529     return -1;
530   }
531
532 // find the header and get Nprimaries and Nsecondaries
533   TTree* tE = (TTree *)file->Get("TE") ;
534   if (!tE) {
535     Error("GetNParticles","input file does not contain TE");
536     return -1;
537   }
538   AliHeader* header;
539   header = 0;
540   tE->SetBranchAddress("Header", &header);
541   if (!tE->GetEntry(event)) {
542     Error("GetNParticles","event %d not found",event);
543     return -1;
544   }
545   if (GetDebug()>2) {
546     cerr<<"Nprimary: "<< header->GetNprimary()<<endl;
547     cerr<<"Nsecondary: "<<header->GetNsecondary()<<endl;
548   }
549   return header->GetNprimary() + header->GetNsecondary();
550 */
551 }
552
553 ////////////////////////////////////////////////////////////////////////
554 Int_t* AliRunDigitizer::GetInputEventNumbers(Int_t event) const
555 {
556 // return pointer to an int array with input event numbers which were
557 // merged in the output event event
558
559 // simplified for now, implement later
560   Int_t * a = new Int_t[kMaxStreamsToMerge];
561   for (Int_t i = 0; i < fNinputs; i++) {
562     a[i] = event;
563   }
564   return a;
565 }
566 ////////////////////////////////////////////////////////////////////////
567 Int_t AliRunDigitizer::GetInputEventNumber(Int_t event, Int_t input) const
568 {
569 // return an event number of an eventInput from input file input
570 // which was merged to create output event event
571
572 // simplified for now, implement later
573   return event;
574 }
575 ////////////////////////////////////////////////////////////////////////
576 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t event) const
577 {
578 // return pointer to particle with index i (index with mask)
579
580 // decode the MASK
581   Int_t input = i/fkMASKSTEP;
582   return GetParticle(i,input,GetInputEventNumber(event,input));
583 }
584
585 ////////////////////////////////////////////////////////////////////////
586 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t input, Int_t event) const
587 {
588 // return pointer to particle with index i in the input file input
589 // (index without mask)
590 // event is the event number in the file input
591 // return 0 i fit does not exist
592
593 // Must be revised in the version with AliStream
594
595   return 0;
596 /*
597   TFile *file = ConnectInputFile(input);
598   if (!file) {
599     Error("GetParticle","Cannot open input file");
600     return 0;
601   }
602
603 // find the header and get Nprimaries and Nsecondaries
604   TTree* tE = (TTree *)file->Get("TE") ;
605   if (!tE) {
606     Error("GetParticle","input file does not contain TE");
607     return 0;
608   }
609   AliHeader* header;
610   header = 0;
611   tE->SetBranchAddress("Header", &header);
612   if (!tE->GetEntry(event)) {
613     Error("GetParticle","event %d not found",event);
614     return 0;
615   }
616   
617 // connect TreeK  
618   char treeName[30];
619   sprintf(treeName,"TreeK%d",event);  
620   TTree* tK = static_cast<TTree*>(file->Get(treeName));
621   if (!tK) {
622     Error("GetParticle","input file does not contain TreeK%d",event);
623     return 0;
624   }
625   TParticle *particleBuffer;
626   particleBuffer = 0;
627   tK->SetBranchAddress("Particles", &particleBuffer);
628
629
630 // algorithmic way of getting entry index
631 // (primary particles are filled after secondaries)
632   Int_t entry;
633   if (i<header->GetNprimary())
634     entry = i+header->GetNsecondary();
635   else 
636     entry = i-header->GetNprimary();
637   Int_t bytesRead = tK->GetEntry(entry);
638 //  new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
639   if (bytesRead)
640     return particleBuffer;
641   return  0;
642 */
643 }
644
645 ////////////////////////////////////////////////////////////////////////
646 void AliRunDigitizer::ExecuteTask(Option_t* option)
647 {
648 // overwrite ExecuteTask to do Digitize only
649
650   if (!IsActive()) return;
651   Digitize(option);
652   fHasExecuted = kTRUE;
653   return;
654 }
655 ////////////////////////////////////////////////////////////////////////
656 TString AliRunDigitizer::GetInputFileName(const Int_t input, const Int_t order) const 
657 {
658 // returns file name of the order-th file in the input stream input
659 // returns empty string if such file does not exist
660 // first input stream is 0
661 // first file in the input stream is 0
662   TString fileName("");
663   if (input >= fNinputs) return fileName;
664   AliStream * stream = static_cast<AliStream*>(fInputStreams->At(input));
665   if (order > stream->GetNInputFiles()) return fileName;
666   fileName = stream->GetFileName(order);
667   return fileName;
668 }
669 ////////////////////////////////////////////////////////////////////////