]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliSimulation.cxx
Code for simulation, sdigitization and digitization moved from macros to compiled...
[u/mrichter/AliRoot.git] / STEER / AliSimulation.cxx
CommitLineData
85a5290f 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/* $Id$ */
17
18///////////////////////////////////////////////////////////////////////////////
19// //
20// class for running generation, simulation and digitization //
21// //
22// Hits, sdigits and digits are created for all detectors by typing: //
23// //
24// AliSimulation sim; //
25// sim.Run(); //
26// //
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 //
29// can be set by //
30// //
31// sim.SetNumberOfEvents(n); //
32// //
33// The name of the configuration file can be specified by //
34// //
35// sim.SetConfigFile("..."); //
36// //
37// The generation of particles and the simulation of detector hits can be //
38// switched on or off by //
39// //
40// sim.SetRunGeneration(kTRUE); // generation of primary particles //
41// sim.SetRunSimulation(kFALSE); // but no tracking //
42// //
43// For which detectors sdigits and digits will be created, can be steered //
44// by //
45// //
46// sim.SetMakeSDigits("ALL"); // make sdigits for all detectors //
47// sim.SetMakeDigits("ITS TPC"); // make digits only for ITS and TPC //
48// //
49// The argument is a (case sensitive) string with the names of the //
50// detectors separated by a space. An empty string ("") can be used to //
51// disable the creation of sdigits or digits. The special string "ALL" //
52// selects all available detectors. This is the default. //
53// //
54// The creation of digits from hits instead of from sdigits can be selected //
55// by //
56// //
57// sim.SetMakeDigitsFromHits("TRD"); //
58// //
59// The argument is again a string with the selected detectors. Be aware that //
60// this feature is not available for all detectors and that merging is not //
61// possible, when digits are created directly from hits. //
62// //
63// Backgound events can be merged by calling //
64// //
65// sim.MergeWith("background/galice.root", 2); //
66// //
67// The first argument is the file name of the background galice file. The //
68// second argument is the number of signal events per background event. //
69// The default value for this is 1. MergeWith can be called several times //
70// to merge more than two event streams. It is assumed that the sdigits //
71// were already produced for the background events. //
72// //
73///////////////////////////////////////////////////////////////////////////////
74
75
76#include "AliSimulation.h"
77#include "AliRunLoader.h"
78#include "AliRun.h"
79#include "AliModule.h"
80#include "AliGenerator.h"
81#include "AliRunDigitizer.h"
82#include <TObjString.h>
83
84
85ClassImp(AliSimulation)
86
87
88//_____________________________________________________________________________
89AliSimulation::AliSimulation(const char* name, const char* title) :
90 TNamed(name, title)
91{
92// create simulation object with default parameters
93
94 Init();
95}
96
97//_____________________________________________________________________________
98AliSimulation::AliSimulation(const AliSimulation& sim) :
99 TNamed(sim)
100{
101// copy constructor
102
103 fRunGeneration = sim.fRunGeneration;
104 fRunSimulation = sim.fRunSimulation;
105 fMakeSDigits = sim.fMakeSDigits;
106 fMakeDigits = sim.fMakeDigits;
107 fMakeDigitsFromHits = sim.fMakeDigitsFromHits;
108 fStopOnError = sim.fStopOnError;
109
110 fNEvents = sim.fNEvents;
111 fConfigFileName = sim.fConfigFileName;
112 fGAliceFileName = sim.fGAliceFileName;
113 fBkgrdFileNames = new TObjArray;
114 for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
115 if (!sim.fBkgrdFileNames->At(i)) continue;
116 fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
117 }
118
119 fRunLoader = NULL;
120}
121
122//_____________________________________________________________________________
123AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
124{
125// assignment operator
126
127 this->~AliSimulation();
128 new(this) AliSimulation(sim);
129 return *this;
130}
131
132//_____________________________________________________________________________
133AliSimulation::~AliSimulation()
134{
135// clean up
136
137 fBkgrdFileNames->Delete();
138 delete fBkgrdFileNames;
139}
140
141//_____________________________________________________________________________
142void AliSimulation::Init()
143{
144// set default parameters
145
146 fRunGeneration = kTRUE;
147 fRunSimulation = kTRUE;
148 fMakeSDigits = "ALL";
149 fMakeDigits = "ALL";
150 fMakeDigitsFromHits = "";
151 fStopOnError = kFALSE;
152
153 fNEvents = 1;
154 fConfigFileName = "Config.C";
155 fGAliceFileName = "galice.root";
156 fBkgrdFileNames = new TObjArray;
157
158 fRunLoader = NULL;
159}
160
161
162//_____________________________________________________________________________
163void AliSimulation::SetNumberOfEvents(Int_t nEvents)
164{
165// set the number of events for one run
166
167 fNEvents = nEvents;
168}
169
170//_____________________________________________________________________________
171void AliSimulation::SetConfigFile(const char* fileName)
172{
173// set the name of the config file
174
175 fConfigFileName = fileName;
176}
177
178//_____________________________________________________________________________
179void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
180{
181// add a file with background events for merging
182
183 TObjString* fileNameStr = new TObjString(fileName);
184 fileNameStr->SetUniqueID(nSignalPerBkgrd);
185 fBkgrdFileNames->Add(fileNameStr);
186}
187
188
189//_____________________________________________________________________________
190Bool_t AliSimulation::Run(Int_t nEvents)
191{
192// run the generation, simulation and digitization
193
194 if (nEvents > 0) fNEvents = nEvents;
195
196 // generation and simulation -> hits
197 if (fRunGeneration) {
198 if (!gAlice) {
199 Error("Run", "no gAlice object. Restart aliroot and try again.");
200 return kFALSE;
201 }
202 if (gAlice->Modules()->GetEntries() > 0) {
203 Error("Run", "gAlice was already run. Restart aliroot and try again.");
204 return kFALSE;
205 }
206 if (!RunSimulation()) if (fStopOnError) return kFALSE;
207 }
208
209 // reopen the run loader
210 if (fRunLoader) delete fRunLoader;
211 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
212 if (!fRunLoader) {
213 Error("Run", "no run loader found in file %s",
214 fGAliceFileName.Data());
215 return kFALSE;
216 }
217 fRunLoader->LoadgAlice();
218 gAlice = fRunLoader->GetAliRun();
219 if (!gAlice) {
220 Error("GetLoadersAndDetectors", "no gAlice object found in file %s",
221 fGAliceFileName.Data());
222 return kFALSE;
223 }
224
225 // hits -> summable digits
226 if (!fMakeSDigits.IsNull()) {
227 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
228 }
229
230 // summable digits -> digits
231 if (!fMakeDigits.IsNull()) {
232 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
233 if (fStopOnError) return kFALSE;
234 }
235 }
236
237 // hits -> digits
238 if (!fMakeDigitsFromHits.IsNull()) {
239 if (fBkgrdFileNames->GetEntriesFast() > 0) {
240 Warning("Run", "Merging and direct creation of digits from hits "
241 "was selected for some detectors. "
242 "No merging will be done for the following detectors: %s",
243 fMakeDigitsFromHits.Data());
244 }
245 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
246 if (fStopOnError) return kFALSE;
247 }
248 }
249
250 return kTRUE;
251}
252
253//_____________________________________________________________________________
254Bool_t AliSimulation::RunSimulation()
255{
256// run the generation and simulation
257
258 Info("RunSimulation", "initializing gAlice with config file %s",
259 fConfigFileName.Data());
260 gAlice->Init(fConfigFileName.Data());
261 fRunLoader = gAlice->GetRunLoader();
262 if (!fRunLoader) {
263 Error("RunSimulation", "gAlice has no run loader object. "
264 "Check your config file: %s", fConfigFileName.Data());
265 return kFALSE;
266 }
267 fGAliceFileName = fRunLoader->GetFileName();
268
269 if (!fRunSimulation) {
270 if (!gAlice->Generator()) {
271 Error("RunSimulation", "gAlice has no generator object. "
272 "Check your config file: %s", fConfigFileName.Data());
273 return kFALSE;
274 }
275 gAlice->Generator()->SetTrackingFlag(0);
276 }
277
278 Info("Run", "running gAlice");
279 gAlice->Run(fNEvents);
280
281 return kTRUE;
282}
283
284//_____________________________________________________________________________
285Bool_t AliSimulation::RunSDigitization(const TString& detectors)
286{
287// run the digitization and produce summable digits
288
289 TString detStr = detectors;
290 TObjArray* detArray = gAlice->Detectors();
291 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
292 AliModule* det = (AliModule*) detArray->At(iDet);
293 if (!det || !det->IsActive()) continue;
294 if (IsSelected(det->GetName(), detStr)) {
295 Info("RunSDigitization", "creating summable digits for %s",
296 det->GetName());
297 det->Hits2SDigits();
298 }
299 }
300
301 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
302 Error("RunSDigitization", "the following detectors were not found: %s",
303 detStr.Data());
304 if (fStopOnError) return kFALSE;
305 }
306
307 return kTRUE;
308}
309
310
311//_____________________________________________________________________________
312Bool_t AliSimulation::RunDigitization(const TString& detectors,
313 const TString& excludeDetectors)
314{
315// run the digitization and produce digits from sdigits
316
317 Int_t nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
318 Int_t signalPerBkgrd = 1;
319 if (nStreams > 1) signalPerBkgrd = fBkgrdFileNames->At(0)->GetUniqueID();
320 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
321 manager->SetInputStream(0, fGAliceFileName.Data());
322 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
323 const char* fileName = ((TObjString*)
324 (fBkgrdFileNames->At(iStream-1)))->GetName();
325 manager->SetInputStream(iStream, fileName);
326 }
327
328 TString detStr = detectors;
329 TString detExcl = excludeDetectors;
330 TObjArray* detArray = gAlice->Detectors();
331 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
332 AliModule* det = (AliModule*) detArray->At(iDet);
333 if (!det || !det->IsActive()) continue;
334 if (IsSelected(det->GetName(), detStr) &&
335 !IsSelected(det->GetName(), detExcl)) {
336 if (!det->CreateDigitizer(manager)) {
337 Error("RunDigitization", "no digitizer for %s", det->GetName());
338 if (fStopOnError) return kFALSE;
339 }
340 }
341 }
342
343 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
344 Error("RunDigitization", "the following detectors were not found: %s",
345 detStr.Data());
346 if (fStopOnError) return kFALSE;
347 }
348
349 if (!manager->GetListOfTasks()->IsEmpty()) {
350 Info("RunDigitization", "executing digitization");
351 manager->Exec("");
352 }
353 delete manager;
354
355 return kTRUE;
356}
357
358//_____________________________________________________________________________
359Bool_t AliSimulation::RunHitsDigitization(const TString& detectors)
360{
361// run the digitization and produce digits from hits
362
363 TString detStr = detectors;
364 TObjArray* detArray = gAlice->Detectors();
365 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
366 AliModule* det = (AliModule*) detArray->At(iDet);
367 if (!det || !det->IsActive()) continue;
368 if (IsSelected(det->GetName(), detStr)) {
369 Info("RunHitsDigitization", "creating digits from hits for %s",
370 det->GetName());
371 det->Hits2Digits();
372 }
373 }
374
375 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
376 Error("RunHitsDigitization", "the following detectors were not found: %s",
377 detStr.Data());
378 if (fStopOnError) return kFALSE;
379 }
380
381 return kTRUE;
382}
383
384
385//_____________________________________________________________________________
386Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
387{
388// check whether detName is contained in detectors
389// if yes, it is removed from detectors
390
391 // check if all detectors are selected
392 if ((detectors.CompareTo("ALL") == 0) ||
393 detectors.BeginsWith("ALL ") ||
394 detectors.EndsWith(" ALL") ||
395 detectors.Contains(" ALL ")) {
396 detectors = "ALL";
397 return kTRUE;
398 }
399
400 // search for the given detector
401 Bool_t result = kFALSE;
402 if ((detectors.CompareTo(detName) == 0) ||
403 detectors.BeginsWith(detName+" ") ||
404 detectors.EndsWith(" "+detName) ||
405 detectors.Contains(" "+detName+" ")) {
406 detectors.ReplaceAll(detName, "");
407 result = kTRUE;
408 }
409
410 // clean up the detectors string
411 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
412 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
413 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
414
415 return result;
416}