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 // class for running generation, simulation and digitization //
22 // Hits, sdigits and digits are created for all detectors by typing: //
24 // AliSimulation sim; //
27 // The Run method returns kTRUE in case of successful execution. //
28 // The number of events can be given as argument to the Run method or it //
31 // sim.SetNumberOfEvents(n); //
33 // The name of the configuration file can be passed as argument to the //
34 // AliSimulation constructor or can be specified by //
36 // sim.SetConfigFile("..."); //
38 // The generation of particles and the simulation of detector hits can be //
39 // switched on or off by //
41 // sim.SetRunGeneration(kTRUE); // generation of primary particles //
42 // sim.SetRunSimulation(kFALSE); // but no tracking //
44 // For which detectors sdigits and digits will be created, can be steered //
47 // sim.SetMakeSDigits("ALL"); // make sdigits for all detectors //
48 // sim.SetMakeDigits("ITS TPC"); // make digits only for ITS and TPC //
50 // The argument is a (case sensitive) string with the names of the //
51 // detectors separated by a space. An empty string ("") can be used to //
52 // disable the creation of sdigits or digits. The special string "ALL" //
53 // selects all available detectors. This is the default. //
55 // The creation of digits from hits instead of from sdigits can be selected //
58 // sim.SetMakeDigitsFromHits("TRD"); //
60 // The argument is again a string with the selected detectors. Be aware that //
61 // this feature is not available for all detectors and that merging is not //
62 // possible, when digits are created directly from hits. //
64 // Backgound events can be merged by calling //
66 // sim.MergeWith("background/galice.root", 2); //
68 // The first argument is the file name of the background galice file. The //
69 // second argument is the number of signal events per background event. //
70 // The default value for this is 1. MergeWith can be called several times //
71 // to merge more than two event streams. It is assumed that the sdigits //
72 // were already produced for the background events. //
74 // The methods RunSimulation, RunSDigitization, RunDigitization and //
75 // RunHitsDigitization can be used to run only parts of the full simulation //
78 ///////////////////////////////////////////////////////////////////////////////
81 #include "AliSimulation.h"
82 #include "AliRunLoader.h"
84 #include "AliModule.h"
85 #include "AliGenerator.h"
86 #include "AliRunDigitizer.h"
87 #include "AliDigitizer.h"
88 #include <TObjString.h>
91 ClassImp(AliSimulation)
94 //_____________________________________________________________________________
95 AliSimulation::AliSimulation(const char* configFileName,
96 const char* name, const char* title) :
99 fRunGeneration(kTRUE),
100 fRunSimulation(kTRUE),
103 fMakeDigitsFromHits(""),
104 fStopOnError(kFALSE),
107 fConfigFileName(configFileName),
108 fGAliceFileName("galice.root"),
109 fBkgrdFileNames(NULL),
110 fRegionOfInterest(kTRUE)
112 // create simulation object with default parameters
116 //_____________________________________________________________________________
117 AliSimulation::AliSimulation(const AliSimulation& sim) :
120 fRunGeneration(sim.fRunGeneration),
121 fRunSimulation(sim.fRunSimulation),
122 fMakeSDigits(sim.fMakeSDigits),
123 fMakeDigits(sim.fMakeDigits),
124 fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
125 fStopOnError(sim.fStopOnError),
127 fNEvents(sim.fNEvents),
128 fConfigFileName(sim.fConfigFileName),
129 fGAliceFileName(sim.fGAliceFileName),
130 fBkgrdFileNames(NULL),
131 fRegionOfInterest(sim.fRegionOfInterest)
135 fBkgrdFileNames = new TObjArray;
136 for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
137 if (!sim.fBkgrdFileNames->At(i)) continue;
138 fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
142 //_____________________________________________________________________________
143 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
145 // assignment operator
147 this->~AliSimulation();
148 new(this) AliSimulation(sim);
152 //_____________________________________________________________________________
153 AliSimulation::~AliSimulation()
157 if (fBkgrdFileNames) {
158 fBkgrdFileNames->Delete();
159 delete fBkgrdFileNames;
164 //_____________________________________________________________________________
165 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
167 // set the number of events for one run
172 //_____________________________________________________________________________
173 void AliSimulation::SetConfigFile(const char* fileName)
175 // set the name of the config file
177 fConfigFileName = fileName;
180 //_____________________________________________________________________________
181 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
183 // add a file with background events for merging
185 TObjString* fileNameStr = new TObjString(fileName);
186 fileNameStr->SetUniqueID(nSignalPerBkgrd);
187 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
188 fBkgrdFileNames->Add(fileNameStr);
192 //_____________________________________________________________________________
193 Bool_t AliSimulation::Run(Int_t nEvents)
195 // run the generation, simulation and digitization
197 if (nEvents > 0) fNEvents = nEvents;
199 // generation and simulation -> hits
200 if (fRunGeneration) {
201 if (!RunSimulation()) if (fStopOnError) return kFALSE;
204 // hits -> summable digits
205 if (!fMakeSDigits.IsNull()) {
206 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
209 // summable digits -> digits
210 if (!fMakeDigits.IsNull()) {
211 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
212 if (fStopOnError) return kFALSE;
217 if (!fMakeDigitsFromHits.IsNull()) {
218 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
219 Warning("Run", "Merging and direct creation of digits from hits "
220 "was selected for some detectors. "
221 "No merging will be done for the following detectors: %s",
222 fMakeDigitsFromHits.Data());
224 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
225 if (fStopOnError) return kFALSE;
232 //_____________________________________________________________________________
233 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
235 // run the generation and simulation
237 TStopwatch stopwatch;
241 Error("RunSimulation", "no gAlice object. Restart aliroot and try again.");
244 if (gAlice->Modules()->GetEntries() > 0) {
245 Error("RunSimulation",
246 "gAlice was already run. Restart aliroot and try again.");
250 Info("RunSimulation", "initializing gAlice with config file %s",
251 fConfigFileName.Data());
252 gAlice->Init(fConfigFileName.Data());
253 AliRunLoader* runLoader = gAlice->GetRunLoader();
255 Error("RunSimulation", "gAlice has no run loader object. "
256 "Check your config file: %s", fConfigFileName.Data());
259 fGAliceFileName = runLoader->GetFileName();
261 if (!fRunSimulation) {
262 if (!gAlice->Generator()) {
263 Error("RunSimulation", "gAlice has no generator object. "
264 "Check your config file: %s", fConfigFileName.Data());
267 gAlice->Generator()->SetTrackingFlag(0);
270 Info("RunSimulation", "running gAlice");
271 if (nEvents <= 0) nEvents = fNEvents;
272 gAlice->Run(nEvents);
276 Info("RunSimulation", "execution time:");
282 //_____________________________________________________________________________
283 Bool_t AliSimulation::RunSDigitization(const char* detectors)
285 // run the digitization and produce summable digits
287 TStopwatch stopwatch;
290 AliRunLoader* runLoader = LoadRun();
291 if (!runLoader) return kFALSE;
293 TString detStr = detectors;
294 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
295 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
296 AliModule* det = (AliModule*) detArray->At(iDet);
297 if (!det || !det->IsActive()) continue;
298 if (IsSelected(det->GetName(), detStr)) {
299 Info("RunSDigitization", "creating summable digits for %s",
305 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
306 Error("RunSDigitization", "the following detectors were not found: %s",
308 if (fStopOnError) return kFALSE;
313 Info("RunSDigitization", "execution time:");
320 //_____________________________________________________________________________
321 Bool_t AliSimulation::RunDigitization(const char* detectors,
322 const char* excludeDetectors)
324 // run the digitization and produce digits from sdigits
326 TStopwatch stopwatch;
329 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
330 if (gAlice) delete gAlice;
334 Int_t signalPerBkgrd = 1;
335 if (fBkgrdFileNames) {
336 nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
337 if (nStreams > 1) signalPerBkgrd = fBkgrdFileNames->At(0)->GetUniqueID();
339 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
340 manager->SetInputStream(0, fGAliceFileName.Data());
341 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
342 const char* fileName = ((TObjString*)
343 (fBkgrdFileNames->At(iStream-1)))->GetName();
344 manager->SetInputStream(iStream, fileName);
347 TString detStr = detectors;
348 TString detExcl = excludeDetectors;
349 manager->GetInputStream(0)->ImportgAlice();
350 AliRunLoader* runLoader =
351 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
352 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
353 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
354 AliModule* det = (AliModule*) detArray->At(iDet);
355 if (!det || !det->IsActive()) continue;
356 if (IsSelected(det->GetName(), detStr) &&
357 !IsSelected(det->GetName(), detExcl)) {
358 AliDigitizer* digitizer = det->CreateDigitizer(manager);
360 Error("RunDigitization", "no digitizer for %s", det->GetName());
361 if (fStopOnError) return kFALSE;
363 digitizer->SetRegionOfInterest(fRegionOfInterest);
368 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
369 Error("RunDigitization", "the following detectors were not found: %s",
371 if (fStopOnError) return kFALSE;
374 if (!manager->GetListOfTasks()->IsEmpty()) {
375 Info("RunDigitization", "executing digitization");
381 Info("RunDigitization", "execution time:");
387 //_____________________________________________________________________________
388 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
390 // run the digitization and produce digits from hits
392 TStopwatch stopwatch;
395 AliRunLoader* runLoader = LoadRun();
396 if (!runLoader) return kFALSE;
398 TString detStr = detectors;
399 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
400 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
401 AliModule* det = (AliModule*) detArray->At(iDet);
402 if (!det || !det->IsActive()) continue;
403 if (IsSelected(det->GetName(), detStr)) {
404 Info("RunHitsDigitization", "creating digits from hits for %s",
410 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
411 Error("RunHitsDigitization", "the following detectors were not found: %s",
413 if (fStopOnError) return kFALSE;
418 Info("RunHitsDigitization", "execution time:");
425 //_____________________________________________________________________________
426 AliRunLoader* AliSimulation::LoadRun() const
428 // delete existing run loaders, open a new one and load gAlice
430 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
431 AliRunLoader* runLoader =
432 AliRunLoader::Open(fGAliceFileName.Data(),
433 AliConfig::fgkDefaultEventFolderName, "UPDATE");
435 Error("LoadRun", "no run loader found in file %s",
436 fGAliceFileName.Data());
439 runLoader->LoadgAlice();
440 gAlice = runLoader->GetAliRun();
442 Error("LoadRun", "no gAlice object found in file %s",
443 fGAliceFileName.Data());
449 //_____________________________________________________________________________
450 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
452 // check whether detName is contained in detectors
453 // if yes, it is removed from detectors
455 // check if all detectors are selected
456 if ((detectors.CompareTo("ALL") == 0) ||
457 detectors.BeginsWith("ALL ") ||
458 detectors.EndsWith(" ALL") ||
459 detectors.Contains(" ALL ")) {
464 // search for the given detector
465 Bool_t result = kFALSE;
466 if ((detectors.CompareTo(detName) == 0) ||
467 detectors.BeginsWith(detName+" ") ||
468 detectors.EndsWith(" "+detName) ||
469 detectors.Contains(" "+detName+" ")) {
470 detectors.ReplaceAll(detName, "");
474 // clean up the detectors string
475 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
476 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
477 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);