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.10 2001/11/15 11:07:25 jchudoba
19 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.
21 Revision 1.9 2001/11/15 09:00:11 jchudoba
22 Add special treatment for TPC and TRD, they use different trees than other detectors
24 Revision 1.8 2001/10/21 18:38:43 hristov
25 Several pointers were set to zero in the default constructors to avoid memory management problems
27 Revision 1.7 2001/10/04 15:56:07 jchudoba
30 Revision 1.4 2001/09/19 06:23:50 jchudoba
31 Move some tasks to AliStream and AliMergeCombi classes
33 Revision 1.3 2001/07/30 14:04:18 jchudoba
34 correct bug in the initialization
36 Revision 1.2 2001/07/28 10:44:32 hristov
37 Loop variable declared once; typos corrected
39 Revision 1.1 2001/07/27 12:59:00 jchudoba
40 Manager class for merging/digitization
44 ////////////////////////////////////////////////////////////////////////
46 // AliRunDigitizer.cxx
48 // Manager object for merging/digitization
50 // Instance of this class manages the digitization and/or merging of
51 // Sdigits into Digits.
53 // Only one instance of this class is created in the macro:
54 // AliRunDigitizer * manager =
55 // new AliRunDigitizer(nInputStreams,SPERB);
56 // where nInputStreams is number of input streams and SPERB is
57 // signals per background variable, which determines how combinations
58 // of signal and background events are generated.
59 // Then instances of specific detector digitizers are created:
60 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager)
61 // and the I/O configured (you have to specify input files
62 // and an output file). The manager connects appropriate trees from
63 // the input files according a combination returned by AliMergeCombi
64 // class. It creates TreeD in the output and runs once per
65 // event Digitize method of all existing AliDetDigitizers
66 // (without any option). AliDetDigitizers ask manager
67 // for a TTree with input (manager->GetInputTreeS(Int_t i),
68 // merge all inputs, digitize it, and save it in the TreeD
69 // obtained by manager->GetTreeD(). Output events are stored with
70 // numbers from 0, this default can be changed by
71 // manager->SetFirstOutputEventNr(Int_t) method. The particle numbers
72 // in the output are shifted by MASK, which is taken from manager.
74 // The default output is to the signal file (stream 0). This can be
75 // changed with the SetOutputFile(TString fn) method.
77 // Single input file is permitted. Maximum kMaxStreamsToMerge can be merged.
78 // Input from the memory (on-the-fly merging) is not yet
79 // supported, as well as access to the input data by invoking methods
80 // on the output data.
82 // Access to the some data is via gAlice for now (supposing the
83 // same geometry in all input files), gAlice is taken from the first
84 // input file on the first stream.
86 // Example with MUON digitizer, no merging, just digitization
88 // AliRunDigitizer * manager = new AliRunDigitizer(1,1);
89 // manager->SetInputStream(0,"galice.root");
90 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
93 // Example with MUON digitizer, merge all events from
94 // galice.root (signal) file with events from bgr.root
95 // (background) file. Number of merged events is
96 // min(number of events in galice.root, number of events in bgr.root)
98 // AliRunDigitizer * manager = new AliRunDigitizer(2,1);
99 // manager->SetInputStream(0,"galice.root");
100 // manager->SetInputStream(1,"bgr.root");
101 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
102 // manager->Exec("");
104 // Example with MUON digitizer, save digits in a new file digits.root,
105 // process only 1 event
107 // AliRunDigitizer * manager = new AliRunDigitizer(2,1);
108 // manager->SetInputStream(0,"galice.root");
109 // manager->SetInputStream(1,"bgr.root");
110 // manager->SetOutputFile("digits.root");
111 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
112 // manager->SetNrOfEventsToWrite(1);
113 // manager->Exec("");
115 ////////////////////////////////////////////////////////////////////////
119 #include <iostream.h>
129 #include "AliRunDigitizer.h"
130 #include "AliDigitizer.h"
132 #include "AliHeader.h"
133 #include "TParticle.h"
134 #include "AliStream.h"
135 #include "AliMergeCombi.h"
137 ClassImp(AliRunDigitizer)
139 ////////////////////////////////////////////////////////////////////////
140 AliRunDigitizer::AliRunDigitizer()
142 // root requires default ctor, where no new objects can be created
143 // do not use this ctor, it is supplied only for root needs
145 // just set all pointers - data members to 0
151 for (Int_t i=0;i<kMaxStreamsToMerge;i++) {
152 fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
159 ////////////////////////////////////////////////////////////////////////
160 AliRunDigitizer::AliRunDigitizer(Int_t nInputStreams, Int_t sperb) : TTask("AliRunDigitizer","The manager for Merging")
162 // ctor which should be used to create a manager for merging/digitization
163 if (nInputStreams == 0) {
164 Error("AliRunDigitizer","Specify nr of input streams");
168 fNinputs = nInputStreams;
169 fOutputFileName = "";
170 fOutputDirName = ".";
171 fCombination.Set(kMaxStreamsToMerge);
172 for (i=0;i<kMaxStreamsToMerge;i++) {
173 fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
176 fkMASKSTEP = 10000000;
178 for (i=1;i<kMaxStreamsToMerge;i++) {
179 fkMASK[i] = fkMASK[i-1] + fkMASKSTEP;
181 fInputStreams = new TClonesArray("AliStream",nInputStreams);
182 TClonesArray &lInputStreams = *fInputStreams;
183 // the first Input is open RW to be output as well
184 new(lInputStreams[0]) AliStream("UPDATE");
185 for (i=1;i<nInputStreams;i++) {
186 new(lInputStreams[i]) AliStream("READ");
190 fNrOfEventsToWrite = -1;
191 fNrOfEventsWritten = 0;
192 fCopyTreesFromInput = -1;
193 fCombi = new AliMergeCombi(nInputStreams,sperb);
199 for (i=0; i<kMaxStreamsToMerge; i++) fInputFiles[i]=0;
202 ////////////////////////////////////////////////////////////////////////
204 AliRunDigitizer::~AliRunDigitizer() {
208 delete fInputStreams;
217 ////////////////////////////////////////////////////////////////////////
218 void AliRunDigitizer::AddDigitizer(AliDigitizer *digitizer)
220 // add digitizer to the list of active digitizers
221 this->Add(digitizer);
223 ////////////////////////////////////////////////////////////////////////
224 void AliRunDigitizer::SetInputStream(Int_t i, char *inputFile)
226 if (i > fInputStreams->GetLast()) {
227 Error("SetInputStream","Input stream number too high");
230 static_cast<AliStream*>(fInputStreams->At(i))->AddFile(inputFile);
233 ////////////////////////////////////////////////////////////////////////
234 void AliRunDigitizer::Digitize(Option_t* option)
236 // get a new combination of inputs, connect input trees and loop
237 // over all digitizers
240 cerr<<"False from InitGlobal"<<endl;
243 // take gAlice from the first input file. It is needed to access
245 if (!static_cast<AliStream*>(fInputStreams->At(0))->ImportgAlice()) {
246 cerr<<"gAlice object not found in the first file of "
247 <<"the 1st stream"<<endl;
250 Int_t eventsCreated = 0;
251 // loop until there is anything on the input in case fNrOfEventsToWrite < 0
252 while ((eventsCreated++ < fNrOfEventsToWrite) || (fNrOfEventsToWrite < 0)) {
253 if (!ConnectInputTrees()) break;
255 // loop over all registered digitizers and let them do the work
263 ////////////////////////////////////////////////////////////////////////
264 Bool_t AliRunDigitizer::ConnectInputTrees()
266 // fill arrays fArrayTreeS, fArrayTreeH and fArrayTreeTPCS with
267 // pointers to the correct events according fCombination values
268 // null pointers can be in the output, AliDigitizer has to check it
273 Int_t eventNr[kMaxStreamsToMerge], delta[kMaxStreamsToMerge];
274 fCombi->Combination(eventNr, delta);
275 for (Int_t i=0;i<fNinputs;i++) {
277 AliStream *iStream = static_cast<AliStream*>(fInputStreams->At(i));
278 if (!iStream->NextEventInStream(serialNr)) return kFALSE;
279 fInputFiles[i]=iStream->CurrentFile();
280 sprintf(treeName,"TreeS%d",serialNr);
281 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
282 fArrayTreeS[i] = tree;
283 sprintf(treeName,"TreeH%d",serialNr);
284 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
285 fArrayTreeH[i] = tree;
286 sprintf(treeName,"TreeS_75x40_100x60_%d",serialNr);
287 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
288 fArrayTreeTPCS[i] = tree;
289 sprintf(treeName,"TreeS%d_TRD",serialNr);
290 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
291 fArrayTreeTRDS[i] = tree;
292 } else if (delta[i] != 0) {
293 Error("ConnectInputTrees","Only delta 0 or 1 is implemented");
300 ////////////////////////////////////////////////////////////////////////
301 Bool_t AliRunDigitizer::InitGlobal()
303 // called once before Digitize() is called, initialize digitizers and output
305 TList* subTasks = this->GetListOfTasks();
307 subTasks->ForEach(AliDigitizer,Init)();
312 ////////////////////////////////////////////////////////////////////////
314 void AliRunDigitizer::SetOutputFile(TString fn)
315 // the output will be to separate file, not to the signal file
317 fOutputFileName = fn;
318 (static_cast<AliStream*>(fInputStreams->At(0)))->ChangeMode("READ");
322 ////////////////////////////////////////////////////////////////////////
323 Bool_t AliRunDigitizer::InitOutputGlobal()
325 // Creates the output file, called by InitEvent()
328 fn = fOutputDirName + '/' + fOutputFileName;
329 fOutput = new TFile(fn,"recreate");
331 cerr<<"AliRunDigitizer::InitOutputGlobal(): file "<<fn.Data()<<" was opened"<<endl;
333 if (fOutput) return kTRUE;
334 Error("InitOutputGlobal","Could not create output file.");
339 ////////////////////////////////////////////////////////////////////////
340 void AliRunDigitizer::InitEvent()
342 // Creates TreeDxx in the output file, called from Digitize() once for
343 // each event. xx = fEvent
346 cerr<<"AliRunDigitizer::InitEvent: fEvent = "<<fEvent<<endl;
348 // if fOutputFileName was not given, write output to signal file
349 if (fOutputFileName == "") {
350 fOutput = (static_cast<AliStream*>(fInputStreams->At(0)))->CurrentFile();
354 sprintf(treeName,"TreeD%d",fEvent);
355 fTreeD = static_cast<TTree*>(fOutput->Get(treeName));
357 fTreeD = new TTree(treeName,"Digits");
358 fTreeD->Write(0,TObject::kOverwrite);
361 // special tree for TPC
362 sprintf(treeName,"TreeD_75x40_100x60_%d",fEvent);
363 fTreeDTPC = static_cast<TTree*>(fOutput->Get(treeName));
365 fTreeDTPC = new TTree(treeName,"TPC_Digits");
366 fTreeDTPC->Write(0,TObject::kOverwrite);
369 // special tree for TRD
370 sprintf(treeName,"TreeD%d_TRD",fEvent);
371 fTreeDTRD = static_cast<TTree*>(fOutput->Get(treeName));
373 fTreeDTRD = new TTree(treeName,"TRD_Digits");
374 fTreeDTRD->Write(0,TObject::kOverwrite);
379 ////////////////////////////////////////////////////////////////////////
380 void AliRunDigitizer::FinishEvent()
382 // called at the end of loop over digitizers
385 if (fCopyTreesFromInput > -1) {
387 Int_t i = fCopyTreesFromInput;
388 sprintf(treeName,"TreeK%d",fCombination[i]);
389 fInputFiles[i]->Get(treeName)->Clone()->Write();
390 sprintf(treeName,"TreeH%d",fCombination[i]);
391 fInputFiles[i]->Get(treeName)->Clone()->Write();
394 fNrOfEventsWritten++;
408 ////////////////////////////////////////////////////////////////////////
409 void AliRunDigitizer::FinishGlobal()
411 // called at the end of Exec
412 // save unique objects to the output file
416 if (fCopyTreesFromInput > -1) {
417 fInputFiles[fCopyTreesFromInput]->Get("TE")->Clone()->Write();
424 ////////////////////////////////////////////////////////////////////////
425 Int_t AliRunDigitizer::GetNParticles(Int_t event) const
427 // return number of particles in all input files for a given
428 // event (as numbered in the output file)
429 // return -1 if some file cannot be accessed
433 for (Int_t i = 0; i < fNinputs; i++) {
434 sumI = GetNParticles(GetInputEventNumber(event,i), i);
435 if (sumI < 0) return -1;
441 ////////////////////////////////////////////////////////////////////////
442 Int_t AliRunDigitizer::GetNParticles(Int_t event, Int_t input) const
444 // return number of particles in input file input for a given
445 // event (as numbered in this input file)
446 // return -1 if some error
448 // Must be revised in the version with AliStream
453 TFile *file = ConnectInputFile(input);
455 Error("GetNParticles","Cannot open input file");
459 // find the header and get Nprimaries and Nsecondaries
460 TTree* tE = (TTree *)file->Get("TE") ;
462 Error("GetNParticles","input file does not contain TE");
467 tE->SetBranchAddress("Header", &header);
468 if (!tE->GetEntry(event)) {
469 Error("GetNParticles","event %d not found",event);
473 cerr<<"Nprimary: "<< header->GetNprimary()<<endl;
474 cerr<<"Nsecondary: "<<header->GetNsecondary()<<endl;
476 return header->GetNprimary() + header->GetNsecondary();
480 ////////////////////////////////////////////////////////////////////////
481 Int_t* AliRunDigitizer::GetInputEventNumbers(Int_t event) const
483 // return pointer to an int array with input event numbers which were
484 // merged in the output event event
486 // simplified for now, implement later
487 Int_t * a = new Int_t[kMaxStreamsToMerge];
488 for (Int_t i = 0; i < fNinputs; i++) {
493 ////////////////////////////////////////////////////////////////////////
494 Int_t AliRunDigitizer::GetInputEventNumber(Int_t event, Int_t input) const
496 // return an event number of an eventInput from input file input
497 // which was merged to create output event event
499 // simplified for now, implement later
502 ////////////////////////////////////////////////////////////////////////
503 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t event) const
505 // return pointer to particle with index i (index with mask)
508 Int_t input = i/fkMASKSTEP;
509 return GetParticle(i,input,GetInputEventNumber(event,input));
512 ////////////////////////////////////////////////////////////////////////
513 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t input, Int_t event) const
515 // return pointer to particle with index i in the input file input
516 // (index without mask)
517 // event is the event number in the file input
518 // return 0 i fit does not exist
520 // Must be revised in the version with AliStream
524 TFile *file = ConnectInputFile(input);
526 Error("GetParticle","Cannot open input file");
530 // find the header and get Nprimaries and Nsecondaries
531 TTree* tE = (TTree *)file->Get("TE") ;
533 Error("GetParticle","input file does not contain TE");
538 tE->SetBranchAddress("Header", &header);
539 if (!tE->GetEntry(event)) {
540 Error("GetParticle","event %d not found",event);
546 sprintf(treeName,"TreeK%d",event);
547 TTree* tK = static_cast<TTree*>(file->Get(treeName));
549 Error("GetParticle","input file does not contain TreeK%d",event);
552 TParticle *particleBuffer;
554 tK->SetBranchAddress("Particles", &particleBuffer);
557 // algorithmic way of getting entry index
558 // (primary particles are filled after secondaries)
560 if (i<header->GetNprimary())
561 entry = i+header->GetNsecondary();
563 entry = i-header->GetNprimary();
564 Int_t bytesRead = tK->GetEntry(entry);
565 // new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
567 return particleBuffer;
572 ////////////////////////////////////////////////////////////////////////
573 void AliRunDigitizer::ExecuteTask(Option_t* option)
575 // overwrite ExecuteTask to do Digitize only
577 if (!IsActive()) return;
579 fHasExecuted = kTRUE;