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.21 2002/10/22 15:02:15 alibrary
19 Introducing Riostream.h
21 Revision 1.20 2002/10/14 14:57:32 hristov
22 Merging the VirtualMC branch to the main development branch (HEAD)
24 Revision 1.13.6.2 2002/07/24 10:08:13 alibrary
27 Revision 1.19 2002/07/19 12:46:05 hristov
28 Write file instead of closing it
30 Revision 1.18 2002/07/17 08:59:39 jchudoba
31 Do not delete subtasks when AliRunDigitizer is deleted. Owner should delete them itself.
33 Revision 1.17 2002/07/16 13:47:53 jchudoba
34 Add methods to get access to names of files used in merging.
36 Revision 1.16 2002/06/07 09:18:47 jchudoba
37 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.
39 Revision 1.15 2002/04/09 13:38:47 jchudoba
40 Add const to the filename argument
42 Revision 1.14 2002/04/04 09:28:04 jchudoba
43 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.
45 Revision 1.13 2002/02/13 09:03:32 jchudoba
46 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.
48 Revision 1.12 2001/12/10 16:40:52 jchudoba
49 Import gAlice from the signal file before InitGlobal() to allow detectors to use it during initialization
51 Revision 1.11 2001/12/03 07:10:13 jchudoba
52 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
54 Revision 1.10 2001/11/15 11:07:25 jchudoba
55 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.
57 Revision 1.9 2001/11/15 09:00:11 jchudoba
58 Add special treatment for TPC and TRD, they use different trees than other detectors
60 Revision 1.8 2001/10/21 18:38:43 hristov
61 Several pointers were set to zero in the default constructors to avoid memory management problems
63 Revision 1.7 2001/10/04 15:56:07 jchudoba
66 Revision 1.4 2001/09/19 06:23:50 jchudoba
67 Move some tasks to AliStream and AliMergeCombi classes
69 Revision 1.3 2001/07/30 14:04:18 jchudoba
70 correct bug in the initialization
72 Revision 1.2 2001/07/28 10:44:32 hristov
73 Loop variable declared once; typos corrected
75 Revision 1.1 2001/07/27 12:59:00 jchudoba
76 Manager class for merging/digitization
80 ////////////////////////////////////////////////////////////////////////
82 // AliRunDigitizer.cxx
84 // Manager object for merging/digitization
86 // Instance of this class manages the digitization and/or merging of
87 // Sdigits into Digits.
89 // Only one instance of this class is created in the macro:
90 // AliRunDigitizer * manager =
91 // new AliRunDigitizer(nInputStreams,SPERB);
92 // where nInputStreams is number of input streams and SPERB is
93 // signals per background variable, which determines how combinations
94 // of signal and background events are generated.
95 // Then instances of specific detector digitizers are created:
96 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager)
97 // and the I/O configured (you have to specify input files
98 // and an output file). The manager connects appropriate trees from
99 // the input files according a combination returned by AliMergeCombi
100 // class. It creates TreeD in the output and runs once per
101 // event Digitize method of all existing AliDetDigitizers
102 // (without any option). AliDetDigitizers ask manager
103 // for a TTree with input (manager->GetInputTreeS(Int_t i),
104 // merge all inputs, digitize it, and save it in the TreeD
105 // obtained by manager->GetTreeD(). Output events are stored with
106 // numbers from 0, this default can be changed by
107 // manager->SetFirstOutputEventNr(Int_t) method. The particle numbers
108 // in the output are shifted by MASK, which is taken from manager.
110 // The default output is to the signal file (stream 0). This can be
111 // changed with the SetOutputFile(TString fn) method.
113 // Single input file is permitted. Maximum kMaxStreamsToMerge can be merged.
114 // Input from the memory (on-the-fly merging) is not yet
115 // supported, as well as access to the input data by invoking methods
116 // on the output data.
118 // Access to the some data is via gAlice for now (supposing the
119 // same geometry in all input files), gAlice is taken from the first
120 // input file on the first stream.
122 // Example with MUON digitizer, no merging, just digitization
124 // AliRunDigitizer * manager = new AliRunDigitizer(1,1);
125 // manager->SetInputStream(0,"galice.root");
126 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
127 // manager->Exec("");
129 // Example with MUON digitizer, merge all events from
130 // galice.root (signal) file with events from bgr.root
131 // (background) file. Number of merged events is
132 // min(number of events in galice.root, number of events in bgr.root)
134 // AliRunDigitizer * manager = new AliRunDigitizer(2,1);
135 // manager->SetInputStream(0,"galice.root");
136 // manager->SetInputStream(1,"bgr.root");
137 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
138 // manager->Exec("");
140 // Example with MUON digitizer, save digits in a new file digits.root,
141 // process only 1 event
143 // AliRunDigitizer * manager = new AliRunDigitizer(2,1);
144 // manager->SetInputStream(0,"galice.root");
145 // manager->SetInputStream(1,"bgr.root");
146 // manager->SetOutputFile("digits.root");
147 // AliMUONDigitizer *dMUON = new AliMUONDigitizer(manager);
148 // manager->SetNrOfEventsToWrite(1);
149 // manager->Exec("");
151 ////////////////////////////////////////////////////////////////////////
155 #include <Riostream.h>
161 #include "TParticle.h"
166 #include "AliDigitizer.h"
167 #include "AliHeader.h"
168 #include "AliMergeCombi.h"
170 #include "AliRunDigitizer.h"
171 #include "AliStream.h"
173 ClassImp(AliRunDigitizer)
175 //_______________________________________________________________________
176 AliRunDigitizer::AliRunDigitizer():
182 fNrOfEventsToWrite(0),
183 fNrOfEventsWritten(0),
184 fCopyTreesFromInput(0),
192 fTreeDTPCBaseName(0),
193 fTreeTPCSBaseName(0),
196 fCombinationFileName(0),
200 // default ctor, where no new objects can be created
201 // do not use this ctor, it is supplied only for root needs
203 for (Int_t i=0;i<kMaxStreamsToMerge;i++) {
204 fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
209 //_______________________________________________________________________
210 AliRunDigitizer::AliRunDigitizer(const AliRunDigitizer& dig):
217 fNrOfEventsToWrite(0),
218 fNrOfEventsWritten(0),
219 fCopyTreesFromInput(0),
227 fTreeDTPCBaseName(0),
228 fTreeTPCSBaseName(0),
231 fCombinationFileName(0),
240 //_______________________________________________________________________
241 void AliRunDigitizer::Copy(AliRunDigitizer&) const
243 Fatal("Copy","Not installed\n");
247 //_______________________________________________________________________
248 AliRunDigitizer::AliRunDigitizer(Int_t nInputStreams, Int_t sperb):
249 TTask("AliRunDigitizer","The manager for Merging"),
250 fkMASKSTEP(10000000),
255 fNrOfEventsToWrite(-1),
256 fNrOfEventsWritten(0),
257 fCopyTreesFromInput(-1),
262 fNinputs(nInputStreams),
264 fInputStreams(new TClonesArray("AliStream",nInputStreams)),
265 fTreeDTPCBaseName("TreeD_75x40_100x60_150x60_"),
266 fTreeTPCSBaseName("TreeS_75x40_100x60_150x60_"),
267 fCombi(new AliMergeCombi(nInputStreams,sperb)),
268 fCombination(kMaxStreamsToMerge),
269 fCombinationFileName(0),
273 // ctor which should be used to create a manager for merging/digitization
275 if (nInputStreams == 0) {
276 Error("AliRunDigitizer","Specify nr of input streams");
279 for (Int_t i=0;i<kMaxStreamsToMerge;i++) {
280 fArrayTreeS[i]=fArrayTreeH[i]=fArrayTreeTPCS[i]=fArrayTreeTRDS[i]=NULL;
284 for (Int_t i=1;i<kMaxStreamsToMerge;i++) {
285 fkMASK[i] = fkMASK[i-1] + fkMASKSTEP;
288 TClonesArray &lInputStreams = *fInputStreams;
289 // the first Input is open RW to be output as well
290 new(lInputStreams[0]) AliStream("UPDATE");
291 for (Int_t i=1;i<nInputStreams;i++) {
292 new(lInputStreams[i]) AliStream("READ");
295 for (Int_t i=0; i<kMaxStreamsToMerge; i++) fInputFiles[i]=0;
298 //_______________________________________________________________________
299 AliRunDigitizer::~AliRunDigitizer() {
302 // do not delete subtasks, let the creator delete them
303 if (GetListOfTasks())
304 GetListOfTasks()->Clear("nodelete");
306 delete fInputStreams;
310 //_______________________________________________________________________
311 void AliRunDigitizer::AddDigitizer(AliDigitizer *digitizer)
313 // add digitizer to the list of active digitizers
314 this->Add(digitizer);
317 //_______________________________________________________________________
318 void AliRunDigitizer::SetInputStream(Int_t i, const char *inputFile)
320 if (i > fInputStreams->GetLast()) {
321 Error("SetInputStream","Input stream number too high");
324 static_cast<AliStream*>(fInputStreams->At(i))->AddFile(inputFile);
327 //_______________________________________________________________________
328 void AliRunDigitizer::Digitize(Option_t* option)
330 // get a new combination of inputs, connect input trees and loop
331 // over all digitizers
333 // take gAlice from the first input file. It is needed to access
335 // If gAlice is already in memory, use it
337 if (!static_cast<AliStream*>(fInputStreams->At(0))->ImportgAlice()) {
338 cerr<<"gAlice object not found in the first file of "
339 <<"the 1st stream"<<endl;
344 cerr<<"False from InitGlobal"<<endl;
347 Int_t eventsCreated = 0;
348 // loop until there is anything on the input in case fNrOfEventsToWrite < 0
349 while ((eventsCreated++ < fNrOfEventsToWrite) || (fNrOfEventsToWrite < 0)) {
350 if (!ConnectInputTrees()) break;
352 // loop over all registered digitizers and let them do the work
353 ExecuteTasks(option);
360 //_______________________________________________________________________
361 Bool_t AliRunDigitizer::ConnectInputTrees()
363 // fill arrays fArrayTreeS, fArrayTreeH and fArrayTreeTPCS with
364 // pointers to the correct events according fCombination values
365 // null pointers can be in the output, AliDigitizer has to check it
370 Int_t eventNr[kMaxStreamsToMerge], delta[kMaxStreamsToMerge];
371 fCombi->Combination(eventNr, delta);
372 for (Int_t i=0;i<fNinputs;i++) {
374 AliStream *iStream = static_cast<AliStream*>(fInputStreams->At(i));
375 if (!iStream->NextEventInStream(serialNr)) return kFALSE;
376 fInputFiles[i]=iStream->CurrentFile();
377 sprintf(treeName,"TreeS%d",serialNr);
378 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
379 if (fArrayTreeS[i]) {
380 delete fArrayTreeS[i];
383 fArrayTreeS[i] = tree;
384 sprintf(treeName,"TreeH%d",serialNr);
385 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
386 if (fArrayTreeH[i]) {
387 delete fArrayTreeH[i];
390 fArrayTreeH[i] = tree;
391 sprintf(treeName,"%s%d",fTreeTPCSBaseName,serialNr);
392 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
393 if (fArrayTreeTPCS[i]) {
394 delete fArrayTreeTPCS[i];
395 fArrayTreeTPCS[i] = 0;
397 fArrayTreeTPCS[i] = tree;
398 sprintf(treeName,"TreeS%d_TRD",serialNr);
399 tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
400 if (fArrayTreeTRDS[i]) {
401 delete fArrayTreeTRDS[i];
402 fArrayTreeTRDS[i] = 0;
404 fArrayTreeTRDS[i] = tree;
405 } else if (delta[i] != 0) {
406 Error("ConnectInputTrees","Only delta 0 or 1 is implemented");
413 //_______________________________________________________________________
414 Bool_t AliRunDigitizer::InitGlobal()
416 // called once before Digitize() is called, initialize digitizers and output
418 TList* subTasks = GetListOfTasks();
420 subTasks->ForEach(AliDigitizer,Init)();
425 //_______________________________________________________________________
426 void AliRunDigitizer::SetOutputFile(TString fn)
427 // the output will be to separate file, not to the signal file
429 fOutputFileName = fn;
430 (static_cast<AliStream*>(fInputStreams->At(0)))->ChangeMode("READ");
434 //_______________________________________________________________________
435 Bool_t AliRunDigitizer::InitOutputGlobal()
437 // Creates the output file, called by InitEvent()
440 fn = fOutputDirName + '/' + fOutputFileName;
441 fOutput = new TFile(fn,"update");
443 cerr<<"AliRunDigitizer::InitOutputGlobal(): file "<<fn.Data()<<" was opened"<<endl;
445 if (fOutput) return kTRUE;
446 Error("InitOutputGlobal","Could not create output file.");
451 //_______________________________________________________________________
452 void AliRunDigitizer::InitEvent()
454 // Creates TreeDxx in the output file, called from Digitize() once for
455 // each event. xx = fEvent
458 cerr<<"AliRunDigitizer::InitEvent: fEvent = "<<fEvent<<endl;
460 // if fOutputFileName was not given, write output to signal file
461 if (fOutputFileName == "") {
462 fOutput = (static_cast<AliStream*>(fInputStreams->At(0)))->CurrentFile();
466 sprintf(treeName,"TreeD%d",fEvent);
467 fTreeD = static_cast<TTree*>(fOutput->Get(treeName));
469 fTreeD = new TTree(treeName,"Digits");
470 fTreeD->Write(0,TObject::kOverwrite);
473 // tree for ITS fast points
474 sprintf(treeName,"TreeR%d",fEvent);
475 fTreeR = static_cast<TTree*>(fOutput->Get(treeName));
477 fTreeR = new TTree(treeName,"Reconstruction");
478 fTreeR->Write(0,TObject::kOverwrite);
481 // special tree for TPC
482 sprintf(treeName,"%s%d",fTreeDTPCBaseName,fEvent);
483 fTreeDTPC = static_cast<TTree*>(fOutput->Get(treeName));
485 fTreeDTPC = new TTree(treeName,"TPC_Digits");
486 fTreeDTPC->Write(0,TObject::kOverwrite);
489 // special tree for TRD
490 sprintf(treeName,"TreeD%d_TRD",fEvent);
491 fTreeDTRD = static_cast<TTree*>(fOutput->Get(treeName));
493 fTreeDTRD = new TTree(treeName,"TRD_Digits");
494 fTreeDTRD->Write(0,TObject::kOverwrite);
499 //_______________________________________________________________________
500 void AliRunDigitizer::FinishEvent()
502 // called at the end of loop over digitizers
506 if (fCopyTreesFromInput > -1) {
508 i = fCopyTreesFromInput;
509 sprintf(treeName,"TreeK%d",fCombination[i]);
510 fInputFiles[i]->Get(treeName)->Clone()->Write();
511 sprintf(treeName,"TreeH%d",fCombination[i]);
512 fInputFiles[i]->Get(treeName)->Clone()->Write();
515 fNrOfEventsWritten++;
534 //_______________________________________________________________________
535 void AliRunDigitizer::FinishGlobal()
537 // called at the end of Exec
538 // save unique objects to the output file
541 this->Write(0,TObject::kOverwrite);
542 if (fCopyTreesFromInput > -1) {
543 fInputFiles[fCopyTreesFromInput]->Get("TE")->Clone()->Write();
550 //_______________________________________________________________________
551 Int_t AliRunDigitizer::GetNParticles(Int_t event) const
553 // return number of particles in all input files for a given
554 // event (as numbered in the output file)
555 // return -1 if some file cannot be accessed
559 for (Int_t i = 0; i < fNinputs; i++) {
560 sumI = GetNParticles(GetInputEventNumber(event,i), i);
561 if (sumI < 0) return -1;
567 //_______________________________________________________________________
568 Int_t AliRunDigitizer::GetNParticles(Int_t /* event */, Int_t /* input */) const
570 // return number of particles in input file input for a given
571 // event (as numbered in this input file)
572 // return -1 if some error
574 // Must be revised in the version with AliStream
579 TFile *file = ConnectInputFile(input);
581 Error("GetNParticles","Cannot open input file");
585 // find the header and get Nprimaries and Nsecondaries
586 TTree* tE = (TTree *)file->Get("TE") ;
588 Error("GetNParticles","input file does not contain TE");
593 tE->SetBranchAddress("Header", &header);
594 if (!tE->GetEntry(event)) {
595 Error("GetNParticles","event %d not found",event);
599 cerr<<"Nprimary: "<< header->GetNprimary()<<endl;
600 cerr<<"Nsecondary: "<<header->GetNsecondary()<<endl;
602 return header->GetNprimary() + header->GetNsecondary();
606 //_______________________________________________________________________
607 Int_t* AliRunDigitizer::GetInputEventNumbers(Int_t event) const
609 // return pointer to an int array with input event numbers which were
610 // merged in the output event event
612 // simplified for now, implement later
613 Int_t * a = new Int_t[kMaxStreamsToMerge];
614 for (Int_t i = 0; i < fNinputs; i++) {
620 //_______________________________________________________________________
621 Int_t AliRunDigitizer::GetInputEventNumber(Int_t event, Int_t /* input */) const
623 // return an event number of an eventInput from input file input
624 // which was merged to create output event event
626 // simplified for now, implement later
630 //_______________________________________________________________________
631 TParticle* AliRunDigitizer::GetParticle(Int_t i, Int_t event) const
633 // return pointer to particle with index i (index with mask)
636 Int_t input = i/fkMASKSTEP;
637 return GetParticle(i,input,GetInputEventNumber(event,input));
640 //_______________________________________________________________________
641 TParticle* AliRunDigitizer::GetParticle(Int_t /* i */, Int_t /* input */,
642 Int_t /* event */) const
644 // return pointer to particle with index i in the input file input
645 // (index without mask)
646 // event is the event number in the file input
647 // return 0 i fit does not exist
649 // Must be revised in the version with AliStream
653 TFile *file = ConnectInputFile(input);
655 Error("GetParticle","Cannot open input file");
659 // find the header and get Nprimaries and Nsecondaries
660 TTree* tE = (TTree *)file->Get("TE") ;
662 Error("GetParticle","input file does not contain TE");
667 tE->SetBranchAddress("Header", &header);
668 if (!tE->GetEntry(event)) {
669 Error("GetParticle","event %d not found",event);
675 sprintf(treeName,"TreeK%d",event);
676 TTree* tK = static_cast<TTree*>(file->Get(treeName));
678 Error("GetParticle","input file does not contain TreeK%d",event);
681 TParticle *particleBuffer;
683 tK->SetBranchAddress("Particles", &particleBuffer);
686 // algorithmic way of getting entry index
687 // (primary particles are filled after secondaries)
689 if (i<header->GetNprimary())
690 entry = i+header->GetNsecondary();
692 entry = i-header->GetNprimary();
693 Int_t bytesRead = tK->GetEntry(entry);
694 // new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
696 return particleBuffer;
701 //_______________________________________________________________________
702 void AliRunDigitizer::ExecuteTask(Option_t* option)
704 // overwrite ExecuteTask to do Digitize only
706 if (!IsActive()) return;
708 fHasExecuted = kTRUE;
712 //_______________________________________________________________________
713 TString AliRunDigitizer::GetInputFileName(const Int_t input, const Int_t order) const
715 // returns file name of the order-th file in the input stream input
716 // returns empty string if such file does not exist
717 // first input stream is 0
718 // first file in the input stream is 0
719 TString fileName("");
720 if (input >= fNinputs) return fileName;
721 AliStream * stream = static_cast<AliStream*>(fInputStreams->At(input));
722 if (order > stream->GetNInputFiles()) return fileName;
723 fileName = stream->GetFileName(order);