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.20 2002/10/14 14:57:32 hristov
19 Merging the VirtualMC branch to the main development branch (HEAD)
21 Revision 1.13.6.2 2002/07/24 10:08:13 alibrary
24 Revision 1.19 2002/07/19 12:46:05 hristov
25 Write file instead of closing it
27 Revision 1.18 2002/07/17 08:59:39 jchudoba
28 Do not delete subtasks when AliRunDigitizer is deleted. Owner should delete them itself.
30 Revision 1.17 2002/07/16 13:47:53 jchudoba
31 Add methods to get access to names of files used in merging.
33 Revision 1.16 2002/06/07 09:18:47 jchudoba
34 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.
36 Revision 1.15 2002/04/09 13:38:47 jchudoba
37 Add const to the filename argument
39 Revision 1.14 2002/04/04 09:28:04 jchudoba
40 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.
42 Revision 1.13 2002/02/13 09:03:32 jchudoba
43 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.
45 Revision 1.12 2001/12/10 16:40:52 jchudoba
46 Import gAlice from the signal file before InitGlobal() to allow detectors to use it during initialization
48 Revision 1.11 2001/12/03 07:10:13 jchudoba
49 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
51 Revision 1.10 2001/11/15 11:07:25 jchudoba
52 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.
54 Revision 1.9 2001/11/15 09:00:11 jchudoba
55 Add special treatment for TPC and TRD, they use different trees than other detectors
57 Revision 1.8 2001/10/21 18:38:43 hristov
58 Several pointers were set to zero in the default constructors to avoid memory management problems
60 Revision 1.7 2001/10/04 15:56:07 jchudoba
63 Revision 1.4 2001/09/19 06:23:50 jchudoba
64 Move some tasks to AliStream and AliMergeCombi classes
66 Revision 1.3 2001/07/30 14:04:18 jchudoba
67 correct bug in the initialization
69 Revision 1.2 2001/07/28 10:44:32 hristov
70 Loop variable declared once; typos corrected
72 Revision 1.1 2001/07/27 12:59:00 jchudoba
73 Manager class for merging/digitization
77 ////////////////////////////////////////////////////////////////////////
79 // AliRunDigitizer.cxx
81 // Manager object for merging/digitization
83 // Instance of this class manages the digitization and/or merging of
84 // Sdigits into Digits.
86 // Only one instance of this class is created in the macro:
87 // AliRunDigitizer * manager =
88 // new AliRunDigitizer(nInputStreams,SPERB);
89 // where nInputStreams is number of input streams and SPERB is
90 // signals per background variable, which determines how combinations
91 // of signal and background events are generated.
92 // Then instances of specific detector digitizers are created:
93 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager)
94 // and the I/O configured (you have to specify input files
95 // and an output file). The manager connects appropriate trees from
96 // the input files according a combination returned by AliMergeCombi
97 // class. It creates TreeD in the output and runs once per
98 // event Digitize method of all existing AliDetDigitizers
99 // (without any option). AliDetDigitizers ask manager
100 // for a TTree with input (manager->GetInputTreeS(Int_t i),
101 // merge all inputs, digitize it, and save it in the TreeD
102 // obtained by manager->GetTreeD(). Output events are stored with
103 // numbers from 0, this default can be changed by
104 // manager->SetFirstOutputEventNr(Int_t) method. The particle numbers
105 // in the output are shifted by MASK, which is taken from manager.
107 // The default output is to the signal file (stream 0). This can be
108 // changed with the SetOutputFile(TString fn) method.
110 // Single input file is permitted. Maximum kMaxStreamsToMerge can be merged.
111 // Input from the memory (on-the-fly merging) is not yet
112 // supported, as well as access to the input data by invoking methods
113 // on the output data.
115 // Access to the some data is via gAlice for now (supposing the
116 // same geometry in all input files), gAlice is taken from the first
117 // input file on the first stream.
119 // Example with MUON digitizer, no merging, just digitization
121 // AliRunDigitizer * manager = new AliRunDigitizer(1,1);
122 // manager->SetInputStream(0,"galice.root");
123 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
124 // manager->Exec("");
126 // Example with MUON digitizer, merge all events from
127 // galice.root (signal) file with events from bgr.root
128 // (background) file. Number of merged events is
129 // min(number of events in galice.root, number of events in bgr.root)
131 // AliRunDigitizer * manager = new AliRunDigitizer(2,1);
132 // manager->SetInputStream(0,"galice.root");
133 // manager->SetInputStream(1,"bgr.root");
134 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
135 // manager->Exec("");
137 // Example with MUON digitizer, save digits in a new file digits.root,
138 // process only 1 event
140 // AliRunDigitizer * manager = new AliRunDigitizer(2,1);
141 // manager->SetInputStream(0,"galice.root");
142 // manager->SetInputStream(1,"bgr.root");
143 // manager->SetOutputFile("digits.root");
144 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
145 // manager->SetNrOfEventsToWrite(1);
146 // manager->Exec("");
148 ////////////////////////////////////////////////////////////////////////
152 #include <Riostream.h>
162 #include "AliRunDigitizer.h"
163 #include "AliDigitizer.h"
165 #include "AliHeader.h"
166 #include "TParticle.h"
167 #include "AliStream.h"
168 #include "AliMergeCombi.h"
170 ClassImp(AliRunDigitizer)
172 ////////////////////////////////////////////////////////////////////////
173 AliRunDigitizer::AliRunDigitizer()
175 // root requires default ctor, where no new objects can be created
176 // do not use this ctor, it is supplied only for root needs
178 // just set all pointers - data members to 0
185 for (Int_t i=0;i<kMaxStreamsToMerge;i++) {
186 fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
193 ////////////////////////////////////////////////////////////////////////
194 AliRunDigitizer::AliRunDigitizer(Int_t nInputStreams, Int_t sperb) : TTask("AliRunDigitizer","The manager for Merging")
196 // ctor which should be used to create a manager for merging/digitization
197 if (nInputStreams == 0) {
198 Error("AliRunDigitizer","Specify nr of input streams");
202 fNinputs = nInputStreams;
203 fOutputFileName = "";
204 fOutputDirName = ".";
205 fCombination.Set(kMaxStreamsToMerge);
206 for (i=0;i<kMaxStreamsToMerge;i++) {
207 fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
210 fkMASKSTEP = 10000000;
212 for (i=1;i<kMaxStreamsToMerge;i++) {
213 fkMASK[i] = fkMASK[i-1] + fkMASKSTEP;
215 fInputStreams = new TClonesArray("AliStream",nInputStreams);
216 TClonesArray &lInputStreams = *fInputStreams;
217 // the first Input is open RW to be output as well
218 new(lInputStreams[0]) AliStream("UPDATE");
219 for (i=1;i<nInputStreams;i++) {
220 new(lInputStreams[i]) AliStream("READ");
224 fNrOfEventsToWrite = -1;
225 fNrOfEventsWritten = 0;
226 fCopyTreesFromInput = -1;
227 fCombi = new AliMergeCombi(nInputStreams,sperb);
233 fTreeDTPCBaseName = "TreeD_75x40_100x60_150x60_";
234 fTreeTPCSBaseName = "TreeS_75x40_100x60_150x60_";
236 for (i=0; i<kMaxStreamsToMerge; i++) fInputFiles[i]=0;
239 ////////////////////////////////////////////////////////////////////////
241 AliRunDigitizer::~AliRunDigitizer() {
244 // do not delete subtasks, let the creator delete them
245 if (GetListOfTasks())
246 GetListOfTasks()->Clear("nodelete");
249 delete fInputStreams;
258 ////////////////////////////////////////////////////////////////////////
259 void AliRunDigitizer::AddDigitizer(AliDigitizer *digitizer)
261 // add digitizer to the list of active digitizers
262 this->Add(digitizer);
264 ////////////////////////////////////////////////////////////////////////
265 void AliRunDigitizer::SetInputStream(Int_t i, const char *inputFile)
267 if (i > fInputStreams->GetLast()) {
268 Error("SetInputStream","Input stream number too high");
271 static_cast<AliStream*>(fInputStreams->At(i))->AddFile(inputFile);
274 ////////////////////////////////////////////////////////////////////////
275 void AliRunDigitizer::Digitize(Option_t* option)
277 // get a new combination of inputs, connect input trees and loop
278 // over all digitizers
280 // take gAlice from the first input file. It is needed to access
282 // If gAlice is already in memory, use it
284 if (!static_cast<AliStream*>(fInputStreams->At(0))->ImportgAlice()) {
285 cerr<<"gAlice object not found in the first file of "
286 <<"the 1st stream"<<endl;
291 cerr<<"False from InitGlobal"<<endl;
294 Int_t eventsCreated = 0;
295 // loop until there is anything on the input in case fNrOfEventsToWrite < 0
296 while ((eventsCreated++ < fNrOfEventsToWrite) || (fNrOfEventsToWrite < 0)) {
297 if (!ConnectInputTrees()) break;
299 // loop over all registered digitizers and let them do the work
300 ExecuteTasks(option);
307 ////////////////////////////////////////////////////////////////////////
308 Bool_t AliRunDigitizer::ConnectInputTrees()
310 // fill arrays fArrayTreeS, fArrayTreeH and fArrayTreeTPCS with
311 // pointers to the correct events according fCombination values
312 // null pointers can be in the output, AliDigitizer has to check it
317 Int_t eventNr[kMaxStreamsToMerge], delta[kMaxStreamsToMerge];
318 fCombi->Combination(eventNr, delta);
319 for (Int_t i=0;i<fNinputs;i++) {
321 AliStream *iStream = static_cast<AliStream*>(fInputStreams->At(i));
322 if (!iStream->NextEventInStream(serialNr)) return kFALSE;
323 fInputFiles[i]=iStream->CurrentFile();
324 sprintf(treeName,"TreeS%d",serialNr);
325 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
326 if (fArrayTreeS[i]) {
327 delete fArrayTreeS[i];
330 fArrayTreeS[i] = tree;
331 sprintf(treeName,"TreeH%d",serialNr);
332 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
333 if (fArrayTreeH[i]) {
334 delete fArrayTreeH[i];
337 fArrayTreeH[i] = tree;
338 sprintf(treeName,"%s%d",fTreeTPCSBaseName,serialNr);
339 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
340 if (fArrayTreeTPCS[i]) {
341 delete fArrayTreeTPCS[i];
342 fArrayTreeTPCS[i] = 0;
344 fArrayTreeTPCS[i] = tree;
345 sprintf(treeName,"TreeS%d_TRD",serialNr);
346 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
347 if (fArrayTreeTRDS[i]) {
348 delete fArrayTreeTRDS[i];
349 fArrayTreeTRDS[i] = 0;
351 fArrayTreeTRDS[i] = tree;
352 } else if (delta[i] != 0) {
353 Error("ConnectInputTrees","Only delta 0 or 1 is implemented");
360 ////////////////////////////////////////////////////////////////////////
361 Bool_t AliRunDigitizer::InitGlobal()
363 // called once before Digitize() is called, initialize digitizers and output
365 TList* subTasks = this->GetListOfTasks();
367 subTasks->ForEach(AliDigitizer,Init)();
372 ////////////////////////////////////////////////////////////////////////
374 void AliRunDigitizer::SetOutputFile(TString fn)
375 // the output will be to separate file, not to the signal file
377 fOutputFileName = fn;
378 (static_cast<AliStream*>(fInputStreams->At(0)))->ChangeMode("READ");
382 ////////////////////////////////////////////////////////////////////////
383 Bool_t AliRunDigitizer::InitOutputGlobal()
385 // Creates the output file, called by InitEvent()
388 fn = fOutputDirName + '/' + fOutputFileName;
389 fOutput = new TFile(fn,"update");
391 cerr<<"AliRunDigitizer::InitOutputGlobal(): file "<<fn.Data()<<" was opened"<<endl;
393 if (fOutput) return kTRUE;
394 Error("InitOutputGlobal","Could not create output file.");
399 ////////////////////////////////////////////////////////////////////////
400 void AliRunDigitizer::InitEvent()
402 // Creates TreeDxx in the output file, called from Digitize() once for
403 // each event. xx = fEvent
406 cerr<<"AliRunDigitizer::InitEvent: fEvent = "<<fEvent<<endl;
408 // if fOutputFileName was not given, write output to signal file
409 if (fOutputFileName == "") {
410 fOutput = (static_cast<AliStream*>(fInputStreams->At(0)))->CurrentFile();
414 sprintf(treeName,"TreeD%d",fEvent);
415 fTreeD = static_cast<TTree*>(fOutput->Get(treeName));
417 fTreeD = new TTree(treeName,"Digits");
418 fTreeD->Write(0,TObject::kOverwrite);
421 // tree for ITS fast points
422 sprintf(treeName,"TreeR%d",fEvent);
423 fTreeR = static_cast<TTree*>(fOutput->Get(treeName));
425 fTreeR = new TTree(treeName,"Reconstruction");
426 fTreeR->Write(0,TObject::kOverwrite);
429 // special tree for TPC
430 sprintf(treeName,"%s%d",fTreeDTPCBaseName,fEvent);
431 fTreeDTPC = static_cast<TTree*>(fOutput->Get(treeName));
433 fTreeDTPC = new TTree(treeName,"TPC_Digits");
434 fTreeDTPC->Write(0,TObject::kOverwrite);
437 // special tree for TRD
438 sprintf(treeName,"TreeD%d_TRD",fEvent);
439 fTreeDTRD = static_cast<TTree*>(fOutput->Get(treeName));
441 fTreeDTRD = new TTree(treeName,"TRD_Digits");
442 fTreeDTRD->Write(0,TObject::kOverwrite);
447 ////////////////////////////////////////////////////////////////////////
448 void AliRunDigitizer::FinishEvent()
450 // called at the end of loop over digitizers
454 if (fCopyTreesFromInput > -1) {
456 i = fCopyTreesFromInput;
457 sprintf(treeName,"TreeK%d",fCombination[i]);
458 fInputFiles[i]->Get(treeName)->Clone()->Write();
459 sprintf(treeName,"TreeH%d",fCombination[i]);
460 fInputFiles[i]->Get(treeName)->Clone()->Write();
463 fNrOfEventsWritten++;
481 ////////////////////////////////////////////////////////////////////////
482 void AliRunDigitizer::FinishGlobal()
484 // called at the end of Exec
485 // save unique objects to the output file
488 this->Write(0,TObject::kOverwrite);
489 if (fCopyTreesFromInput > -1) {
490 fInputFiles[fCopyTreesFromInput]->Get("TE")->Clone()->Write();
497 ////////////////////////////////////////////////////////////////////////
498 Int_t AliRunDigitizer::GetNParticles(Int_t event) const
500 // return number of particles in all input files for a given
501 // event (as numbered in the output file)
502 // return -1 if some file cannot be accessed
506 for (Int_t i = 0; i < fNinputs; i++) {
507 sumI = GetNParticles(GetInputEventNumber(event,i), i);
508 if (sumI < 0) return -1;
514 ////////////////////////////////////////////////////////////////////////
515 Int_t AliRunDigitizer::GetNParticles(Int_t event, Int_t input) const
517 // return number of particles in input file input for a given
518 // event (as numbered in this input file)
519 // return -1 if some error
521 // Must be revised in the version with AliStream
526 TFile *file = ConnectInputFile(input);
528 Error("GetNParticles","Cannot open input file");
532 // find the header and get Nprimaries and Nsecondaries
533 TTree* tE = (TTree *)file->Get("TE") ;
535 Error("GetNParticles","input file does not contain TE");
540 tE->SetBranchAddress("Header", &header);
541 if (!tE->GetEntry(event)) {
542 Error("GetNParticles","event %d not found",event);
546 cerr<<"Nprimary: "<< header->GetNprimary()<<endl;
547 cerr<<"Nsecondary: "<<header->GetNsecondary()<<endl;
549 return header->GetNprimary() + header->GetNsecondary();
553 ////////////////////////////////////////////////////////////////////////
554 Int_t* AliRunDigitizer::GetInputEventNumbers(Int_t event) const
556 // return pointer to an int array with input event numbers which were
557 // merged in the output event event
559 // simplified for now, implement later
560 Int_t * a = new Int_t[kMaxStreamsToMerge];
561 for (Int_t i = 0; i < fNinputs; i++) {
566 ////////////////////////////////////////////////////////////////////////
567 Int_t AliRunDigitizer::GetInputEventNumber(Int_t event, Int_t input) const
569 // return an event number of an eventInput from input file input
570 // which was merged to create output event event
572 // simplified for now, implement later
575 ////////////////////////////////////////////////////////////////////////
576 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t event) const
578 // return pointer to particle with index i (index with mask)
581 Int_t input = i/fkMASKSTEP;
582 return GetParticle(i,input,GetInputEventNumber(event,input));
585 ////////////////////////////////////////////////////////////////////////
586 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t input, Int_t event) const
588 // return pointer to particle with index i in the input file input
589 // (index without mask)
590 // event is the event number in the file input
591 // return 0 i fit does not exist
593 // Must be revised in the version with AliStream
597 TFile *file = ConnectInputFile(input);
599 Error("GetParticle","Cannot open input file");
603 // find the header and get Nprimaries and Nsecondaries
604 TTree* tE = (TTree *)file->Get("TE") ;
606 Error("GetParticle","input file does not contain TE");
611 tE->SetBranchAddress("Header", &header);
612 if (!tE->GetEntry(event)) {
613 Error("GetParticle","event %d not found",event);
619 sprintf(treeName,"TreeK%d",event);
620 TTree* tK = static_cast<TTree*>(file->Get(treeName));
622 Error("GetParticle","input file does not contain TreeK%d",event);
625 TParticle *particleBuffer;
627 tK->SetBranchAddress("Particles", &particleBuffer);
630 // algorithmic way of getting entry index
631 // (primary particles are filled after secondaries)
633 if (i<header->GetNprimary())
634 entry = i+header->GetNsecondary();
636 entry = i-header->GetNprimary();
637 Int_t bytesRead = tK->GetEntry(entry);
638 // new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
640 return particleBuffer;
645 ////////////////////////////////////////////////////////////////////////
646 void AliRunDigitizer::ExecuteTask(Option_t* option)
648 // overwrite ExecuteTask to do Digitize only
650 if (!IsActive()) return;
652 fHasExecuted = kTRUE;
655 ////////////////////////////////////////////////////////////////////////
656 TString AliRunDigitizer::GetInputFileName(const Int_t input, const Int_t order) const
658 // returns file name of the order-th file in the input stream input
659 // returns empty string if such file does not exist
660 // first input stream is 0
661 // first file in the input stream is 0
662 TString fileName("");
663 if (input >= fNinputs) return fileName;
664 AliStream * stream = static_cast<AliStream*>(fInputStreams->At(input));
665 if (order > stream->GetNInputFiles()) return fileName;
666 fileName = stream->GetFileName(order);
669 ////////////////////////////////////////////////////////////////////////