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 Revision 1.4 2001/09/19 06:23:50 jchudoba
19 Move some tasks to AliStream and AliMergeCombi classes
21 Revision 1.3 2001/07/30 14:04:18 jchudoba
22 correct bug in the initialization
24 Revision 1.2 2001/07/28 10:44:32 hristov
25 Loop variable declared once; typos corrected
27 Revision 1.1 2001/07/27 12:59:00 jchudoba
28 Manager class for merging/digitization
32 ////////////////////////////////////////////////////////////////////////
34 // AliRunDigitizer.cxx
36 // Manager object for merging/digitization
38 // Instance of this class manages the digitization and/or merging of
39 // Sdigits into Digits.
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.
62 // The default output is to the signal file (stream 0). This can be
63 // changed with the SetOutputFile(TString fn) method.
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.
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.
74 // Example with MUON digitizer:
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);
84 ////////////////////////////////////////////////////////////////////////
98 #include "AliRunDigitizer.h"
99 #include "AliDigitizer.h"
101 #include "AliHeader.h"
102 #include "TParticle.h"
103 #include "AliStream.h"
104 #include "AliMergeCombi.h"
106 ClassImp(AliRunDigitizer)
108 ////////////////////////////////////////////////////////////////////////
109 AliRunDigitizer::AliRunDigitizer(Int_t nInputStreams, Int_t sperb) : TTask("AliRunDigitizer","The manager for Merging")
112 if (nInputStreams == 0) {
113 Error("AliRunDigitizer","Specify nr of input streams");
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;
125 fkMASKSTEP = 10000000;
127 for (i=1;i<MAXSTREAMSTOMERGE;i++) {
128 fkMASK[i] = fkMASK[i-1] + fkMASKSTEP;
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();
139 fNrOfEventsToWrite = -1;
140 fNrOfEventsWritten = 0;
141 fCopyTreesFromInput = -1;
142 fCombi = new AliMergeCombi(nInputStreams,sperb);
145 cerr<<"AliRunDigitizer::AliRunDigitizer() called"<<endl;
148 ////////////////////////////////////////////////////////////////////////
150 AliRunDigitizer::~AliRunDigitizer() {
154 delete fInputStreams;
163 ////////////////////////////////////////////////////////////////////////
164 void AliRunDigitizer::AddDigitizer(AliDigitizer *digitizer)
166 // add digitizer to the list of active digitizers
167 this->Add(digitizer);
169 ////////////////////////////////////////////////////////////////////////
170 void AliRunDigitizer::SetInputStream(Int_t i, char *inputFile)
172 if (i > fInputStreams->GetLast()) {
173 Error("SetInputStream","Input stream number too high");
176 static_cast<AliStream*>(fInputStreams->At(i))->AddFile(inputFile);
179 ////////////////////////////////////////////////////////////////////////
180 void AliRunDigitizer::Digitize(Option_t* option)
182 // get a new combination of inputs, connect input trees and loop
183 // over all digitizers
186 cerr<<"False from InitGlobal"<<endl;
189 // take gAlice from the first input file. It is needed to access
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;
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;
201 // loop over all registered digitizers and let them do the work
209 ////////////////////////////////////////////////////////////////////////
210 Bool_t AliRunDigitizer::ConnectInputTrees()
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
219 Int_t eventNr[MAXSTREAMSTOMERGE], delta[MAXSTREAMSTOMERGE];
220 fCombi->Combination(eventNr, delta);
221 for (Int_t i=0;i<fNinputs;i++) {
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");
243 ////////////////////////////////////////////////////////////////////////
244 Bool_t AliRunDigitizer::InitGlobal()
246 // called once before Digitize() is called, initialize digitizers and output
248 TList* subTasks = this->GetListOfTasks();
250 subTasks->ForEach(AliDigitizer,Init)();
255 ////////////////////////////////////////////////////////////////////////
257 void AliRunDigitizer::SetOutputFile(TString fn)
258 // the output will be to separate file, not to the signal file
260 fOutputFileName = fn;
261 (static_cast<AliStream*>(fInputStreams->At(0)))->ChangeMode("READ");
265 ////////////////////////////////////////////////////////////////////////
266 Bool_t AliRunDigitizer::InitOutputGlobal()
268 // Creates the output file, called by InitEvent()
271 fn = fOutputDirName + '/' + fOutputFileName;
272 fOutput = new TFile(fn,"recreate");
274 cerr<<"AliRunDigitizer::InitOutputGlobal(): file "<<fn.Data()<<" was opened"<<endl;
276 if (fOutput) return kTRUE;
277 Error("InitOutputGlobal","Could not create output file.");
282 ////////////////////////////////////////////////////////////////////////
283 void AliRunDigitizer::InitEvent()
285 // Creates TreeDxx in the output file, called from Digitize() once for
286 // each event. xx = fEvent
289 cerr<<"AliRunDigitizer::InitEvent: fEvent = "<<fEvent<<endl;
291 // if fOutputFileName was not given, write output to signal file
292 if (fOutputFileName == "") {
293 fOutput = (static_cast<AliStream*>(fInputStreams->At(0)))->CurrentFile();
297 sprintf(treeName,"TreeD%d",fEvent);
298 fTreeD = static_cast<TTree*>(fOutput->Get(treeName));
300 fTreeD = new TTree(treeName,"Digits");
301 fTreeD->Write(0,TObject::kOverwrite);
305 ////////////////////////////////////////////////////////////////////////
306 void AliRunDigitizer::FinishEvent()
308 // called at the end of loop over digitizers
311 if (fCopyTreesFromInput > -1) {
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();
320 fNrOfEventsWritten++;
326 ////////////////////////////////////////////////////////////////////////
327 void AliRunDigitizer::FinishGlobal()
329 // called at the end of Exec
330 // save unique objects to the output file
334 if (fCopyTreesFromInput > -1) {
335 fInputFiles[fCopyTreesFromInput]->Get("TE")->Clone()->Write();
342 ////////////////////////////////////////////////////////////////////////
343 Int_t AliRunDigitizer::GetNParticles(Int_t event)
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
351 for (Int_t i = 0; i < fNinputs; i++) {
352 sumI = GetNParticles(GetInputEventNumber(event,i), i);
353 if (sumI < 0) return -1;
359 ////////////////////////////////////////////////////////////////////////
360 Int_t AliRunDigitizer::GetNParticles(Int_t event, Int_t input)
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
366 // Must be revised in the version with AliStream
371 TFile *file = ConnectInputFile(input);
373 Error("GetNParticles","Cannot open input file");
377 // find the header and get Nprimaries and Nsecondaries
378 TTree* tE = (TTree *)file->Get("TE") ;
380 Error("GetNParticles","input file does not contain TE");
385 tE->SetBranchAddress("Header", &header);
386 if (!tE->GetEntry(event)) {
387 Error("GetNParticles","event %d not found",event);
391 cerr<<"Nprimary: "<< header->GetNprimary()<<endl;
392 cerr<<"Nsecondary: "<<header->GetNsecondary()<<endl;
394 return header->GetNprimary() + header->GetNsecondary();
398 ////////////////////////////////////////////////////////////////////////
399 Int_t* AliRunDigitizer::GetInputEventNumbers(Int_t event)
401 // return pointer to an int array with input event numbers which were
402 // merged in the output event event
404 // simplified for now, implement later
405 Int_t * a = new Int_t[MAXSTREAMSTOMERGE];
406 for (Int_t i = 0; i < fNinputs; i++) {
411 ////////////////////////////////////////////////////////////////////////
412 Int_t AliRunDigitizer::GetInputEventNumber(Int_t event, Int_t input)
414 // return an event number of an eventInput from input file input
415 // which was merged to create output event event
417 // simplified for now, implement later
420 ////////////////////////////////////////////////////////////////////////
421 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t event)
423 // return pointer to particle with index i (index with mask)
426 Int_t input = i/fkMASKSTEP;
427 return GetParticle(i,input,GetInputEventNumber(event,input));
430 ////////////////////////////////////////////////////////////////////////
431 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t input, Int_t event)
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
438 // Must be revised in the version with AliStream
442 TFile *file = ConnectInputFile(input);
444 Error("GetParticle","Cannot open input file");
448 // find the header and get Nprimaries and Nsecondaries
449 TTree* tE = (TTree *)file->Get("TE") ;
451 Error("GetParticle","input file does not contain TE");
456 tE->SetBranchAddress("Header", &header);
457 if (!tE->GetEntry(event)) {
458 Error("GetParticle","event %d not found",event);
464 sprintf(treeName,"TreeK%d",event);
465 TTree* tK = static_cast<TTree*>(file->Get(treeName));
467 Error("GetParticle","input file does not contain TreeK%d",event);
470 TParticle *particleBuffer;
472 tK->SetBranchAddress("Particles", &particleBuffer);
475 // algorithmic way of getting entry index
476 // (primary particles are filled after secondaries)
478 if (i<header->GetNprimary())
479 entry = i+header->GetNsecondary();
481 entry = i-header->GetNprimary();
482 Int_t bytesRead = tK->GetEntry(entry);
483 // new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
485 return particleBuffer;
490 ////////////////////////////////////////////////////////////////////////
491 void AliRunDigitizer::ExecuteTask(Option_t* option)
493 // overwrite ExecuteTask to do Digitize only
495 if (!IsActive()) return;
497 fHasExecuted = kTRUE;