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