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 **************************************************************************/
18 ////////////////////////////////////////////////////////////////////////
20 // AliRunDigitizer.cxx
22 // Manager object for merging/digitization
24 // Instance of this class manages the digitization and/or merging of
25 // Sdigits into Digits.
27 // Only one instance of this class is created in the macro:
28 // AliRunDigitizer * manager =
29 // new AliRunDigitizer(nInputStreams,SPERB);
30 // where nInputStreams is number of input streams and SPERB is
31 // signals per background variable, which determines how combinations
32 // of signal and background events are generated.
33 // Then instances of specific detector digitizers are created:
34 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager)
35 // and the I/O configured (you have to specify input files
36 // and an output file). The manager connects appropriate trees from
37 // the input files according a combination returned by AliMergeCombi
38 // class. It creates TreeD in the output and runs once per
39 // event Digitize method of all existing AliDetDigitizers
40 // (without any option). AliDetDigitizers ask manager
41 // for a TTree with input (manager->GetInputTreeS(Int_t i),
42 // merge all inputs, digitize it, and save it in the TreeD
43 // obtained by manager->GetTreeD(). Output events are stored with
44 // numbers from 0, this default can be changed by
45 // manager->SetFirstOutputEventNr(Int_t) method. The particle numbers
46 // in the output are shifted by MASK, which is taken from manager.
48 // The default output is to the signal file (stream 0). This can be
49 // changed with the SetOutputFile(TString fn) method.
51 // Single input file is permitted. Maximum kMaxStreamsToMerge can be merged.
52 // Input from the memory (on-the-fly merging) is not yet
53 // supported, as well as access to the input data by invoking methods
54 // on the output data.
56 // Access to the some data is via gAlice for now (supposing the
57 // same geometry in all input files), gAlice is taken from the first
58 // input file on the first stream.
60 // Example with MUON digitizer, no merging, just digitization
62 // AliRunDigitizer * manager = new AliRunDigitizer(1,1);
63 // manager->SetInputStream(0,"galice.root");
64 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
67 // Example with MUON digitizer, merge all events from
68 // galice.root (signal) file with events from bgr.root
69 // (background) file. Number of merged events is
70 // min(number of events in galice.root, number of events in bgr.root)
72 // AliRunDigitizer * manager = new AliRunDigitizer(2,1);
73 // manager->SetInputStream(0,"galice.root");
74 // manager->SetInputStream(1,"bgr.root");
75 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
78 // Example with MUON digitizer, save digits in a new file digits.root,
79 // process only 1 event
81 // AliRunDigitizer * manager = new AliRunDigitizer(2,1);
82 // manager->SetInputStream(0,"galice.root");
83 // manager->SetInputStream(1,"bgr.root");
84 // manager->SetOutputFile("digits.root");
85 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
86 // manager->SetNrOfEventsToWrite(1);
89 ////////////////////////////////////////////////////////////////////////
93 #include <Riostream.h>
103 #include "AliDigitizer.h"
104 #include "AliMergeCombi.h"
106 #include "AliRunDigitizer.h"
107 #include "AliStream.h"
109 ClassImp(AliRunDigitizer)
111 //_______________________________________________________________________
112 AliRunDigitizer::AliRunDigitizer():
118 fNrOfEventsToWrite(0),
119 fNrOfEventsWritten(0),
120 fCopyTreesFromInput(0),
128 fTreeDTPCBaseName(0),
129 fTreeTPCSBaseName(0),
132 fCombinationFileName(0),
136 // default ctor, where no new objects can be created
137 // do not use this ctor, it is supplied only for root needs
139 for (Int_t i=0;i<kMaxStreamsToMerge;i++) {
140 fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
145 //_______________________________________________________________________
146 AliRunDigitizer::AliRunDigitizer(const AliRunDigitizer& dig):
153 fNrOfEventsToWrite(0),
154 fNrOfEventsWritten(0),
155 fCopyTreesFromInput(0),
163 fTreeDTPCBaseName(0),
164 fTreeTPCSBaseName(0),
167 fCombinationFileName(0),
176 //_______________________________________________________________________
177 void AliRunDigitizer::Copy(AliRunDigitizer&) const
179 Fatal("Copy","Not installed\n");
183 //_______________________________________________________________________
184 AliRunDigitizer::AliRunDigitizer(Int_t nInputStreams, Int_t sperb):
185 TTask("AliRunDigitizer","The manager for Merging"),
186 fkMASKSTEP(10000000),
191 fNrOfEventsToWrite(-1),
192 fNrOfEventsWritten(0),
193 fCopyTreesFromInput(-1),
198 fNinputs(nInputStreams),
200 fInputStreams(new TClonesArray("AliStream",nInputStreams)),
201 fTreeDTPCBaseName("TreeD_75x40_100x60_150x60_"),
202 fTreeTPCSBaseName("TreeS_75x40_100x60_150x60_"),
203 fCombi(new AliMergeCombi(nInputStreams,sperb)),
204 fCombination(kMaxStreamsToMerge),
205 fCombinationFileName(0),
209 // ctor which should be used to create a manager for merging/digitization
211 if (nInputStreams == 0) {
212 Error("AliRunDigitizer","Specify nr of input streams");
215 for (Int_t i=0;i<kMaxStreamsToMerge;i++) {
216 fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
220 for (Int_t i=1;i<kMaxStreamsToMerge;i++) {
221 fkMASK[i] = fkMASK[i-1] + fkMASKSTEP;
224 TClonesArray &lInputStreams = *fInputStreams;
225 // the first Input is open RW to be output as well
226 new(lInputStreams[0]) AliStream("UPDATE");
227 for (Int_t i=1;i<nInputStreams;i++) {
228 new(lInputStreams[i]) AliStream("READ");
231 for (Int_t i=0; i<kMaxStreamsToMerge; i++) fInputFiles[i]=0;
234 //_______________________________________________________________________
235 AliRunDigitizer::~AliRunDigitizer() {
238 // do not delete subtasks, let the creator delete them
239 if (GetListOfTasks())
240 GetListOfTasks()->Clear("nodelete");
242 delete fInputStreams;
246 //_______________________________________________________________________
247 void AliRunDigitizer::AddDigitizer(AliDigitizer *digitizer)
249 // add digitizer to the list of active digitizers
250 this->Add(digitizer);
253 //_______________________________________________________________________
254 void AliRunDigitizer::SetInputStream(Int_t i, const char *inputFile)
257 // Sets the name of the input file
259 if (i > fInputStreams->GetLast()) {
260 Error("SetInputStream","Input stream number too high");
263 static_cast<AliStream*>(fInputStreams->At(i))->AddFile(inputFile);
266 //_______________________________________________________________________
267 void AliRunDigitizer::Digitize(Option_t* option)
269 // get a new combination of inputs, connect input trees and loop
270 // over all digitizers
272 // take gAlice from the first input file. It is needed to access
274 // If gAlice is already in memory, use it
276 if (!static_cast<AliStream*>(fInputStreams->At(0))->ImportgAlice()) {
277 cerr<<"gAlice object not found in the first file of "
278 <<"the 1st stream"<<endl;
283 cerr<<"False from InitGlobal"<<endl;
286 Int_t eventsCreated = 0;
287 // loop until there is anything on the input in case fNrOfEventsToWrite < 0
288 while ((eventsCreated++ < fNrOfEventsToWrite) || (fNrOfEventsToWrite < 0)) {
289 if (!ConnectInputTrees()) break;
291 // loop over all registered digitizers and let them do the work
292 ExecuteTasks(option);
299 //_______________________________________________________________________
300 Bool_t AliRunDigitizer::ConnectInputTrees()
302 // fill arrays fArrayTreeS, fArrayTreeH and fArrayTreeTPCS with
303 // pointers to the correct events according fCombination values
304 // null pointers can be in the output, AliDigitizer has to check it
309 Int_t eventNr[kMaxStreamsToMerge], delta[kMaxStreamsToMerge];
310 fCombi->Combination(eventNr, delta);
311 for (Int_t i=0;i<fNinputs;i++) {
313 AliStream *iStream = static_cast<AliStream*>(fInputStreams->At(i));
314 if (!iStream->NextEventInStream(serialNr)) return kFALSE;
315 fInputFiles[i]=iStream->CurrentFile();
316 sprintf(treeName,"TreeS%d",serialNr);
317 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
318 if (fArrayTreeS[i]) {
319 delete fArrayTreeS[i];
322 fArrayTreeS[i] = tree;
323 sprintf(treeName,"TreeH%d",serialNr);
324 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
325 if (fArrayTreeH[i]) {
326 delete fArrayTreeH[i];
329 fArrayTreeH[i] = tree;
330 sprintf(treeName,"%s%d",fTreeTPCSBaseName,serialNr);
331 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
332 if (fArrayTreeTPCS[i]) {
333 delete fArrayTreeTPCS[i];
334 fArrayTreeTPCS[i] = 0;
336 fArrayTreeTPCS[i] = tree;
337 sprintf(treeName,"TreeS%d_TRD",serialNr);
338 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
339 if (fArrayTreeTRDS[i]) {
340 delete fArrayTreeTRDS[i];
341 fArrayTreeTRDS[i] = 0;
343 fArrayTreeTRDS[i] = tree;
344 } else if (delta[i] != 0) {
345 Error("ConnectInputTrees","Only delta 0 or 1 is implemented");
352 //_______________________________________________________________________
353 Bool_t AliRunDigitizer::InitGlobal()
355 // called once before Digitize() is called, initialize digitizers and output
357 TList* subTasks = GetListOfTasks();
359 // subTasks->ForEach(AliDigitizer,Init)();
360 TIter next(subTasks);
361 while (AliDigitizer * dig = (AliDigitizer *) next())
367 //_______________________________________________________________________
368 void AliRunDigitizer::SetOutputFile(TString fn)
370 // the output will be to separate file, not to the signal file
371 fOutputFileName = fn;
372 (static_cast<AliStream*>(fInputStreams->At(0)))->ChangeMode("READ");
376 //_______________________________________________________________________
377 Bool_t AliRunDigitizer::InitOutputGlobal()
379 // Creates the output file, called by InitEvent()
382 fn = fOutputDirName + '/' + fOutputFileName;
383 fOutput = new TFile(fn,"update");
385 cerr<<"AliRunDigitizer::InitOutputGlobal(): file "<<fn.Data()<<" was opened"<<endl;
387 if (fOutput) return kTRUE;
388 Error("InitOutputGlobal","Could not create output file.");
393 //_______________________________________________________________________
394 void AliRunDigitizer::InitEvent()
396 // Creates TreeDxx in the output file, called from Digitize() once for
397 // each event. xx = fEvent
400 cerr<<"AliRunDigitizer::InitEvent: fEvent = "<<fEvent<<endl;
402 // if fOutputFileName was not given, write output to signal file
403 if (fOutputFileName == "") {
404 fOutput = (static_cast<AliStream*>(fInputStreams->At(0)))->CurrentFile();
408 sprintf(treeName,"TreeD%d",fEvent);
409 fTreeD = static_cast<TTree*>(fOutput->Get(treeName));
411 fTreeD = new TTree(treeName,"Digits");
412 fTreeD->Write(0,TObject::kOverwrite);
415 // tree for ITS fast points
416 sprintf(treeName,"TreeR%d",fEvent);
417 fTreeR = static_cast<TTree*>(fOutput->Get(treeName));
419 fTreeR = new TTree(treeName,"Reconstruction");
420 fTreeR->Write(0,TObject::kOverwrite);
423 // special tree for TPC
424 sprintf(treeName,"%s%d",fTreeDTPCBaseName,fEvent);
425 fTreeDTPC = static_cast<TTree*>(fOutput->Get(treeName));
427 fTreeDTPC = new TTree(treeName,"TPC_Digits");
428 fTreeDTPC->Write(0,TObject::kOverwrite);
431 // special tree for TRD
432 sprintf(treeName,"TreeD%d_TRD",fEvent);
433 fTreeDTRD = static_cast<TTree*>(fOutput->Get(treeName));
435 fTreeDTRD = new TTree(treeName,"TRD_Digits");
436 fTreeDTRD->Write(0,TObject::kOverwrite);
441 //_______________________________________________________________________
442 void AliRunDigitizer::FinishEvent()
444 // called at the end of loop over digitizers
448 if (fCopyTreesFromInput > -1) {
450 i = fCopyTreesFromInput;
451 sprintf(treeName,"TreeK%d",fCombination[i]);
452 fInputFiles[i]->Get(treeName)->Clone()->Write();
453 sprintf(treeName,"TreeH%d",fCombination[i]);
454 fInputFiles[i]->Get(treeName)->Clone()->Write();
457 fNrOfEventsWritten++;
476 //_______________________________________________________________________
477 void AliRunDigitizer::FinishGlobal()
479 // called at the end of Exec
480 // save unique objects to the output file
483 this->Write(0,TObject::kOverwrite);
484 if (fCopyTreesFromInput > -1) {
485 fInputFiles[fCopyTreesFromInput]->Get("TE")->Clone()->Write();
492 //_______________________________________________________________________
493 Int_t AliRunDigitizer::GetNParticles(Int_t event) const
495 // return number of particles in all input files for a given
496 // event (as numbered in the output file)
497 // return -1 if some file cannot be accessed
501 for (Int_t i = 0; i < fNinputs; i++) {
502 sumI = GetNParticles(GetInputEventNumber(event,i), i);
503 if (sumI < 0) return -1;
509 //_______________________________________________________________________
510 Int_t AliRunDigitizer::GetNParticles(Int_t /* event */, Int_t /* input */) const
512 // return number of particles in input file input for a given
513 // event (as numbered in this input file)
514 // return -1 if some error
516 // Must be revised in the version with AliStream
521 TFile *file = ConnectInputFile(input);
523 Error("GetNParticles","Cannot open input file");
527 // find the header and get Nprimaries and Nsecondaries
528 TTree* tE = (TTree *)file->Get("TE") ;
530 Error("GetNParticles","input file does not contain TE");
535 tE->SetBranchAddress("Header", &header);
536 if (!tE->GetEntry(event)) {
537 Error("GetNParticles","event %d not found",event);
541 cerr<<"Nprimary: "<< header->GetNprimary()<<endl;
542 cerr<<"Nsecondary: "<<header->GetNsecondary()<<endl;
544 return header->GetNprimary() + header->GetNsecondary();
548 //_______________________________________________________________________
549 Int_t* AliRunDigitizer::GetInputEventNumbers(Int_t event) const
551 // return pointer to an int array with input event numbers which were
552 // merged in the output event event
554 // simplified for now, implement later
555 Int_t * a = new Int_t[kMaxStreamsToMerge];
556 for (Int_t i = 0; i < fNinputs; i++) {
562 //_______________________________________________________________________
563 Int_t AliRunDigitizer::GetInputEventNumber(Int_t event, Int_t /* input */) const
565 // return an event number of an eventInput from input file input
566 // which was merged to create output event event
568 // simplified for now, implement later
572 //_______________________________________________________________________
573 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t event) const
575 // return pointer to particle with index i (index with mask)
578 Int_t input = i/fkMASKSTEP;
579 return GetParticle(i,input,GetInputEventNumber(event,input));
582 //_______________________________________________________________________
583 TParticle* AliRunDigitizer::GetParticle(Int_t /* i */, Int_t /* input */,
584 Int_t /* event */) const
586 // return pointer to particle with index i in the input file input
587 // (index without mask)
588 // event is the event number in the file input
589 // return 0 i fit does not exist
591 // Must be revised in the version with AliStream
595 TFile *file = ConnectInputFile(input);
597 Error("GetParticle","Cannot open input file");
601 // find the header and get Nprimaries and Nsecondaries
602 TTree* tE = (TTree *)file->Get("TE") ;
604 Error("GetParticle","input file does not contain TE");
609 tE->SetBranchAddress("Header", &header);
610 if (!tE->GetEntry(event)) {
611 Error("GetParticle","event %d not found",event);
617 sprintf(treeName,"TreeK%d",event);
618 TTree* tK = static_cast<TTree*>(file->Get(treeName));
620 Error("GetParticle","input file does not contain TreeK%d",event);
623 TParticle *particleBuffer;
625 tK->SetBranchAddress("Particles", &particleBuffer);
628 // algorithmic way of getting entry index
629 // (primary particles are filled after secondaries)
631 if (i<header->GetNprimary())
632 entry = i+header->GetNsecondary();
634 entry = i-header->GetNprimary();
635 Int_t bytesRead = tK->GetEntry(entry);
636 // new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
638 return particleBuffer;
643 //_______________________________________________________________________
644 void AliRunDigitizer::ExecuteTask(Option_t* option)
646 // overwrite ExecuteTask to do Digitize only
648 if (!IsActive()) return;
650 fHasExecuted = kTRUE;
654 //_______________________________________________________________________
655 TString AliRunDigitizer::GetInputFileName(const Int_t input, const Int_t order) const
657 // returns file name of the order-th file in the input stream input
658 // returns empty string if such file does not exist
659 // first input stream is 0
660 // first file in the input stream is 0
661 TString fileName("");
662 if (input >= fNinputs) return fileName;
663 AliStream * stream = static_cast<AliStream*>(fInputStreams->At(input));
664 if (order > stream->GetNInputFiles()) return fileName;
665 fileName = stream->GetFileName(order);