]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRunDigitizer.cxx
8584fcd6999581cd2c0425cb3dd7b955a3a95f76
[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.4  2001/09/19 06:23:50  jchudoba
19 Move some tasks to AliStream and AliMergeCombi classes
20
21 Revision 1.3  2001/07/30 14:04:18  jchudoba
22 correct bug in the initialization
23
24 Revision 1.2  2001/07/28 10:44:32  hristov
25 Loop variable declared once; typos corrected
26
27 Revision 1.1  2001/07/27 12:59:00  jchudoba
28 Manager class for merging/digitization
29
30 */
31
32 ////////////////////////////////////////////////////////////////////////
33 //
34 // AliRunDigitizer.cxx
35 //
36 // Manager object for merging/digitization
37 //
38 // Instance of this class manages the digitization and/or merging of
39 // Sdigits into Digits. 
40 //
41 // Only one instance of this class is created in the macro:
42 //   AliRunDigitizer * manager = 
43 //      new AliRunDigitizer(nInputStreams,SPERB);
44 // where nInputStreams is number of input streams and SPERB is
45 // signals per background variable, which determines how combinations
46 // of signal and background events are generated.
47 // Then instances of specific detector digitizers are created:
48 //   AliMUONDigitizer *dMUON  = new AliMUONDigitizer(manager)
49 // and the I/O configured (you have to specify input files 
50 // and an output file). The manager connects appropriate trees from 
51 // the input files according a combination returned by AliMergeCombi 
52 // classcreates. It creates TreeD in the output and runs once per 
53 // event Digitize method of all existing AliDetDigitizers 
54 // (without any option). AliDetDigitizers ask manager
55 // for a TTree with input (manager->GetInputTreeS(Int_t i),
56 // merge all inputs, digitize it, and save it in the TreeD 
57 // obtained by manager->GetTreeD(). Output events are stored with 
58 // numbers from 0, this default can be changed by 
59 // manager->SetFirstOutputEventNr(Int_t) method. The particle numbers
60 // in the output are shifted by MASK, which is taken from manager.
61 //
62 // The default output is to the signal file (stream 0). This can be 
63 // changed with the SetOutputFile(TString fn)  method.
64 //
65 // Single input file is permitted. Maximum MAXSTREAMSTOMERGE can be merged.
66 // Input from the memory (on-the-fly merging) is not yet 
67 // supported, as well as access to the input data by invoking methods
68 // on the output data.
69 //
70 // Access to the some data is via gAlice for now (supposing the 
71 // same geometry in all input files), gAlice is taken from the first 
72 // input file on the first stream.
73 //
74 // Example with MUON digitizer:
75 //
76 //  AliRunDigitizer * manager = new AliRunDigitizer(2,1);
77 //  manager->SetInputStream(0,"1track_10events_phi45_60.root");
78 //  manager->SetInputStream(1,"1track_10events_phi120_135.root");
79 //  manager->SetOutputFile("digits.root");
80 //  AliMUONDigitizer *dMUON  = new AliMUONDigitizer(manager);
81 //  manager->SetNrOfEventsToWrite(1);
82 //  manager->Exec("");
83 //
84 //////////////////////////////////////////////////////////////////////// 
85
86 // system includes
87
88 #include <iostream.h>
89
90 // ROOT includes
91
92 #include "TFile.h"
93 #include "TTree.h"
94 #include "TList.h"
95
96 // AliROOT includes
97
98 #include "AliRunDigitizer.h"
99 #include "AliDigitizer.h"
100 #include "AliRun.h"
101 #include "AliHeader.h"
102 #include "TParticle.h"
103 #include "AliStream.h"
104 #include "AliMergeCombi.h"
105
106 ClassImp(AliRunDigitizer)
107
108 ////////////////////////////////////////////////////////////////////////
109 AliRunDigitizer::AliRunDigitizer(Int_t nInputStreams, Int_t sperb) : TTask("AliRunDigitizer","The manager for Merging")
110 {
111 // default ctor
112   if (nInputStreams == 0) {
113     Error("AliRunDigitizer","Specify nr of input streams");
114     return;
115   }
116   Int_t i;
117   fNinputs = nInputStreams;
118   fOutputFileName = "";
119   fOutputDirName = ".";
120   fCombination.Set(MAXSTREAMSTOMERGE);
121   for (i=0;i<MAXSTREAMSTOMERGE;i++) {
122     fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=NULL;
123     fCombination[i]=-1;
124   }
125   fkMASKSTEP = 10000000;
126   fkMASK[0] = 0;
127   for (i=1;i<MAXSTREAMSTOMERGE;i++) {
128     fkMASK[i] = fkMASK[i-1] + fkMASKSTEP;
129   }
130   fInputStreams = new TClonesArray("AliStream",nInputStreams);
131   TClonesArray &lInputStreams = *fInputStreams;
132 // the first Input is open RW to be output as well
133   new(lInputStreams[0]) AliStream("UPDATE");
134   for (i=1;i<nInputStreams;i++) {
135     new(lInputStreams[i]) AliStream();
136   }
137   fOutput = 0;
138   fEvent = 0;
139   fNrOfEventsToWrite = -1;
140   fNrOfEventsWritten = 0;
141   fCopyTreesFromInput = -1;
142   fCombi = new AliMergeCombi(nInputStreams,sperb);
143   fDebug = 0;
144   if (GetDebug()>2) 
145     cerr<<"AliRunDigitizer::AliRunDigitizer() called"<<endl;
146 }
147
148 ////////////////////////////////////////////////////////////////////////
149
150 AliRunDigitizer::~AliRunDigitizer() {
151 // dtor
152
153   if (fInputStreams) {
154     delete fInputStreams;
155     fInputStreams = 0;
156   }
157   if (fCombi) {
158     delete fCombi;
159     fCombi = 0;
160   }
161
162 }
163 ////////////////////////////////////////////////////////////////////////
164 void AliRunDigitizer::AddDigitizer(AliDigitizer *digitizer)
165 {
166 // add digitizer to the list of active digitizers
167   this->Add(digitizer);
168 }
169 ////////////////////////////////////////////////////////////////////////
170 void AliRunDigitizer::SetInputStream(Int_t i, char *inputFile)
171 {
172   if (i > fInputStreams->GetLast()) {
173     Error("SetInputStream","Input stream number too high");
174     return;
175   }
176   static_cast<AliStream*>(fInputStreams->At(i))->AddFile(inputFile);
177 }
178
179 ////////////////////////////////////////////////////////////////////////
180 void AliRunDigitizer::Digitize(Option_t* option)
181 {
182 // get a new combination of inputs, connect input trees and loop 
183 // over all digitizers
184
185   if (!InitGlobal()) {
186     cerr<<"False from InitGlobal"<<endl;
187     return;
188   }
189 // take gAlice from the first input file. It is needed to access
190 //  geometry data
191   if (!static_cast<AliStream*>(fInputStreams->At(0))->ImportgAlice()) {
192     cerr<<"gAlice object not found in the first file of "
193         <<"the 1st stream"<<endl;
194     return;
195   }
196   Int_t eventsCreated = 0;
197 // loop until there is anything on the input in case fNrOfEventsToWrite < 0
198   while ((eventsCreated++ < fNrOfEventsToWrite) || (fNrOfEventsToWrite < 0)) {
199     if (!ConnectInputTrees()) break;
200     InitEvent();
201 // loop over all registered digitizers and let them do the work
202     ExecuteTasks("");
203     CleanTasks();
204     FinishEvent();
205   }
206   FinishGlobal();
207 }
208
209 ////////////////////////////////////////////////////////////////////////
210 Bool_t AliRunDigitizer::ConnectInputTrees()
211 {
212 // fill arrays fArrayTreeS, fArrayTreeH and fArrayTreeTPCS with 
213 // pointers to the correct events according fCombination values
214 // null pointers can be in the output, AliDigitizer has to check it
215
216   TTree *tree;
217   char treeName[50];
218   Int_t serialNr;
219   Int_t eventNr[MAXSTREAMSTOMERGE], delta[MAXSTREAMSTOMERGE];
220   fCombi->Combination(eventNr, delta);
221   for (Int_t i=0;i<fNinputs;i++) {
222     if (delta[i] == 1) {
223       AliStream *iStream = static_cast<AliStream*>(fInputStreams->At(i));
224       if (!iStream->NextEventInStream(serialNr)) return kFALSE;
225       fInputFiles[i]=iStream->CurrentFile();
226       sprintf(treeName,"TreeS%d",serialNr);
227       tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
228       fArrayTreeS[i] = tree;
229       sprintf(treeName,"TreeH%d",serialNr);
230       tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
231       fArrayTreeH[i] = tree;
232       sprintf(treeName,"TreeS_75x40_100x60_%d",serialNr);
233       tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
234       fArrayTreeTPCS[i] = tree;
235     } else if (delta[i] != 0) {
236       Error("ConnectInputTrees","Only delta 0 or 1 is implemented");
237       return kFALSE;
238     }
239   }
240   return kTRUE;
241 }
242
243 ////////////////////////////////////////////////////////////////////////
244 Bool_t AliRunDigitizer::InitGlobal()
245 {
246 // called once before Digitize() is called, initialize digitizers and output
247
248   TList* subTasks = this->GetListOfTasks();
249   if (subTasks) {
250     subTasks->ForEach(AliDigitizer,Init)();
251   }  
252   return kTRUE;
253 }
254
255 ////////////////////////////////////////////////////////////////////////
256
257 void AliRunDigitizer::SetOutputFile(TString fn)
258 // the output will be to separate file, not to the signal file
259 {
260   fOutputFileName = fn;
261   (static_cast<AliStream*>(fInputStreams->At(0)))->ChangeMode("READ");
262   InitOutputGlobal();
263 }
264
265 ////////////////////////////////////////////////////////////////////////
266 Bool_t AliRunDigitizer::InitOutputGlobal()
267 {
268 // Creates the output file, called by InitEvent()
269
270   TString fn;
271   fn = fOutputDirName + '/' + fOutputFileName;
272   fOutput = new TFile(fn,"recreate");
273   if (GetDebug()>2) {
274     cerr<<"AliRunDigitizer::InitOutputGlobal(): file "<<fn.Data()<<" was opened"<<endl;
275   }
276   if (fOutput) return kTRUE;
277   Error("InitOutputGlobal","Could not create output file.");
278   return kFALSE;
279 }
280
281
282 ////////////////////////////////////////////////////////////////////////
283 void AliRunDigitizer::InitEvent()
284 {
285 // Creates TreeDxx in the output file, called from Digitize() once for 
286 //  each event. xx = fEvent
287
288   if (GetDebug()>2) 
289     cerr<<"AliRunDigitizer::InitEvent: fEvent = "<<fEvent<<endl;
290
291 // if fOutputFileName was not given, write output to signal file
292   if (fOutputFileName == "") {
293     fOutput = (static_cast<AliStream*>(fInputStreams->At(0)))->CurrentFile();
294   }
295   fOutput->cd();
296   char treeName[30];
297   sprintf(treeName,"TreeD%d",fEvent);
298   fTreeD = static_cast<TTree*>(fOutput->Get(treeName));
299   if (!fTreeD) {
300     fTreeD = new TTree(treeName,"Digits");
301     fTreeD->Write(0,TObject::kOverwrite);
302   }
303 }
304
305 ////////////////////////////////////////////////////////////////////////
306 void AliRunDigitizer::FinishEvent()
307 {
308 // called at the end of loop over digitizers
309
310   fOutput->cd();
311   if (fCopyTreesFromInput > -1) {
312     char treeName[20];
313     Int_t i = fCopyTreesFromInput; 
314     sprintf(treeName,"TreeK%d",fCombination[i]);
315     fInputFiles[i]->Get(treeName)->Clone()->Write();
316     sprintf(treeName,"TreeH%d",fCombination[i]);
317     fInputFiles[i]->Get(treeName)->Clone()->Write();
318   }
319   fEvent++;
320   fNrOfEventsWritten++;
321   if (fTreeD) {
322     delete fTreeD;
323     fTreeD = 0;
324   }
325 }
326 ////////////////////////////////////////////////////////////////////////
327 void AliRunDigitizer::FinishGlobal()
328 {
329 // called at the end of Exec
330 // save unique objects to the output file
331
332   fOutput->cd();
333   this->Write();
334   if (fCopyTreesFromInput > -1) {
335     fInputFiles[fCopyTreesFromInput]->Get("TE")->Clone()->Write();
336     gAlice->Write();
337   }
338   fOutput->Close();
339 }
340
341
342 ////////////////////////////////////////////////////////////////////////
343 Int_t  AliRunDigitizer::GetNParticles(Int_t event)
344 {
345 // return number of particles in all input files for a given
346 // event (as numbered in the output file)
347 // return -1 if some file cannot be accessed
348
349   Int_t sum = 0;
350   Int_t sumI;
351   for (Int_t i = 0; i < fNinputs; i++) {
352     sumI = GetNParticles(GetInputEventNumber(event,i), i);
353     if (sumI < 0) return -1;
354     sum += sumI;
355   }
356   return sum;
357 }
358
359 ////////////////////////////////////////////////////////////////////////
360 Int_t  AliRunDigitizer::GetNParticles(Int_t event, Int_t input)
361 {
362 // return number of particles in input file input for a given
363 // event (as numbered in this input file)
364 // return -1 if some error
365
366 // Must be revised in the version with AliStream
367
368   return -1;
369
370 /*
371   TFile *file = ConnectInputFile(input);
372   if (!file) {
373     Error("GetNParticles","Cannot open input file");
374     return -1;
375   }
376
377 // find the header and get Nprimaries and Nsecondaries
378   TTree* tE = (TTree *)file->Get("TE") ;
379   if (!tE) {
380     Error("GetNParticles","input file does not contain TE");
381     return -1;
382   }
383   AliHeader* header;
384   header = 0;
385   tE->SetBranchAddress("Header", &header);
386   if (!tE->GetEntry(event)) {
387     Error("GetNParticles","event %d not found",event);
388     return -1;
389   }
390   if (GetDebug()>2) {
391     cerr<<"Nprimary: "<< header->GetNprimary()<<endl;
392     cerr<<"Nsecondary: "<<header->GetNsecondary()<<endl;
393   }
394   return header->GetNprimary() + header->GetNsecondary();
395 */
396 }
397
398 ////////////////////////////////////////////////////////////////////////
399 Int_t* AliRunDigitizer::GetInputEventNumbers(Int_t event)
400 {
401 // return pointer to an int array with input event numbers which were
402 // merged in the output event event
403
404 // simplified for now, implement later
405   Int_t * a = new Int_t[MAXSTREAMSTOMERGE];
406   for (Int_t i = 0; i < fNinputs; i++) {
407     a[i] = event;
408   }
409   return a;
410 }
411 ////////////////////////////////////////////////////////////////////////
412 Int_t AliRunDigitizer::GetInputEventNumber(Int_t event, Int_t input)
413 {
414 // return an event number of an eventInput from input file input
415 // which was merged to create output event event
416
417 // simplified for now, implement later
418   return event;
419 }
420 ////////////////////////////////////////////////////////////////////////
421 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t event)
422 {
423 // return pointer to particle with index i (index with mask)
424
425 // decode the MASK
426   Int_t input = i/fkMASKSTEP;
427   return GetParticle(i,input,GetInputEventNumber(event,input));
428 }
429
430 ////////////////////////////////////////////////////////////////////////
431 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t input, Int_t event)
432 {
433 // return pointer to particle with index i in the input file input
434 // (index without mask)
435 // event is the event number in the file input
436 // return 0 i fit does not exist
437
438 // Must be revised in the version with AliStream
439
440   return 0;
441 /*
442   TFile *file = ConnectInputFile(input);
443   if (!file) {
444     Error("GetParticle","Cannot open input file");
445     return 0;
446   }
447
448 // find the header and get Nprimaries and Nsecondaries
449   TTree* tE = (TTree *)file->Get("TE") ;
450   if (!tE) {
451     Error("GetParticle","input file does not contain TE");
452     return 0;
453   }
454   AliHeader* header;
455   header = 0;
456   tE->SetBranchAddress("Header", &header);
457   if (!tE->GetEntry(event)) {
458     Error("GetParticle","event %d not found",event);
459     return 0;
460   }
461   
462 // connect TreeK  
463   char treeName[30];
464   sprintf(treeName,"TreeK%d",event);  
465   TTree* tK = static_cast<TTree*>(file->Get(treeName));
466   if (!tK) {
467     Error("GetParticle","input file does not contain TreeK%d",event);
468     return 0;
469   }
470   TParticle *particleBuffer;
471   particleBuffer = 0;
472   tK->SetBranchAddress("Particles", &particleBuffer);
473
474
475 // algorithmic way of getting entry index
476 // (primary particles are filled after secondaries)
477   Int_t entry;
478   if (i<header->GetNprimary())
479     entry = i+header->GetNsecondary();
480   else 
481     entry = i-header->GetNprimary();
482   Int_t bytesRead = tK->GetEntry(entry);
483 //  new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
484   if (bytesRead)
485     return particleBuffer;
486   return  0;
487 */
488 }
489
490 ////////////////////////////////////////////////////////////////////////
491 void AliRunDigitizer::ExecuteTask(Option_t* option)
492 {
493 // overwrite ExecuteTask to do Digitize only
494
495   if (!IsActive()) return;
496   Digitize(option);
497   fHasExecuted = kTRUE;
498   return;
499 }