]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliSimulation.cxx
Splitting of TRD library (T.Kuhr)
[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// //
95601147 33// The name of the configuration file can be passed as argument to the //
34// AliSimulation constructor or can be specified by //
85a5290f 35// //
36// sim.SetConfigFile("..."); //
37// //
38// The generation of particles and the simulation of detector hits can be //
39// switched on or off by //
40// //
41// sim.SetRunGeneration(kTRUE); // generation of primary particles //
42// sim.SetRunSimulation(kFALSE); // but no tracking //
43// //
44// For which detectors sdigits and digits will be created, can be steered //
45// by //
46// //
47// sim.SetMakeSDigits("ALL"); // make sdigits for all detectors //
48// sim.SetMakeDigits("ITS TPC"); // make digits only for ITS and TPC //
49// //
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. //
54// //
55// The creation of digits from hits instead of from sdigits can be selected //
56// by //
57// //
58// sim.SetMakeDigitsFromHits("TRD"); //
59// //
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. //
63// //
05526d44 64// Background events can be merged by calling //
85a5290f 65// //
66// sim.MergeWith("background/galice.root", 2); //
67// //
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. //
05526d44 70// By default this number is calculated from the number of available //
71// background events. MergeWith can be called several times to merge more //
72// than two event streams. It is assumed that the sdigits were already //
73// produced for the background events. //
85a5290f 74// //
0421c3d1 75// The output of raw data can be switched on by calling //
76// //
77// sim.SetWriteRawData("MUON"); // write raw data for MUON //
78// //
79// The methods RunSimulation, RunSDigitization, RunDigitization, //
80// RunHitsDigitization and WriteRawData can be used to run only parts of //
81// the full simulation chain. //
95601147 82// //
85a5290f 83///////////////////////////////////////////////////////////////////////////////
84
85a5290f 85#include <TObjString.h>
fd46e2d2 86#include <TStopwatch.h>
af7ba10c 87#include <TSystem.h>
85a5290f 88
af7ba10c 89#include "AliDigitizer.h"
90#include "AliGenerator.h"
91#include "AliModule.h"
92#include "AliRun.h"
93#include "AliRunDigitizer.h"
94#include "AliRunLoader.h"
95#include "AliSimulation.h"
96#include "AliVertexGenFile.h"
85a5290f 97
98ClassImp(AliSimulation)
99
100
101//_____________________________________________________________________________
95601147 102AliSimulation::AliSimulation(const char* configFileName,
103 const char* name, const char* title) :
104 TNamed(name, title),
105
106 fRunGeneration(kTRUE),
107 fRunSimulation(kTRUE),
108 fMakeSDigits("ALL"),
109 fMakeDigits("ALL"),
110 fMakeDigitsFromHits(""),
0421c3d1 111 fWriteRawData(""),
95601147 112 fStopOnError(kFALSE),
113
114 fNEvents(1),
115 fConfigFileName(configFileName),
116 fGAliceFileName("galice.root"),
117 fBkgrdFileNames(NULL),
04bae0a0 118 fUseBkgrdVertex(kTRUE),
95601147 119 fRegionOfInterest(kTRUE)
85a5290f 120{
121// create simulation object with default parameters
122
0421c3d1 123 SetGAliceFile("galice.root");
85a5290f 124}
125
126//_____________________________________________________________________________
127AliSimulation::AliSimulation(const AliSimulation& sim) :
95601147 128 TNamed(sim),
129
130 fRunGeneration(sim.fRunGeneration),
131 fRunSimulation(sim.fRunSimulation),
132 fMakeSDigits(sim.fMakeSDigits),
133 fMakeDigits(sim.fMakeDigits),
134 fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
0421c3d1 135 fWriteRawData(sim.fWriteRawData),
95601147 136 fStopOnError(sim.fStopOnError),
137
138 fNEvents(sim.fNEvents),
139 fConfigFileName(sim.fConfigFileName),
140 fGAliceFileName(sim.fGAliceFileName),
141 fBkgrdFileNames(NULL),
04bae0a0 142 fUseBkgrdVertex(sim.fUseBkgrdVertex),
95601147 143 fRegionOfInterest(sim.fRegionOfInterest)
85a5290f 144{
145// copy constructor
146
85a5290f 147 fBkgrdFileNames = new TObjArray;
148 for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
149 if (!sim.fBkgrdFileNames->At(i)) continue;
150 fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
151 }
85a5290f 152}
153
154//_____________________________________________________________________________
155AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
156{
157// assignment operator
158
159 this->~AliSimulation();
160 new(this) AliSimulation(sim);
161 return *this;
162}
163
164//_____________________________________________________________________________
165AliSimulation::~AliSimulation()
166{
167// clean up
168
95601147 169 if (fBkgrdFileNames) {
170 fBkgrdFileNames->Delete();
171 delete fBkgrdFileNames;
172 }
85a5290f 173}
174
175
176//_____________________________________________________________________________
177void AliSimulation::SetNumberOfEvents(Int_t nEvents)
178{
179// set the number of events for one run
180
181 fNEvents = nEvents;
182}
183
184//_____________________________________________________________________________
185void AliSimulation::SetConfigFile(const char* fileName)
186{
187// set the name of the config file
188
189 fConfigFileName = fileName;
190}
191
0421c3d1 192//_____________________________________________________________________________
193void AliSimulation::SetGAliceFile(const char* fileName)
194{
195// set the name of the galice file
196// the path is converted to an absolute one if it is relative
197
198 fGAliceFileName = fileName;
199 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
200 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
201 fGAliceFileName);
202 fGAliceFileName = absFileName;
203 delete[] absFileName;
204 }
205}
206
85a5290f 207//_____________________________________________________________________________
208void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
209{
210// add a file with background events for merging
211
212 TObjString* fileNameStr = new TObjString(fileName);
213 fileNameStr->SetUniqueID(nSignalPerBkgrd);
95601147 214 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
85a5290f 215 fBkgrdFileNames->Add(fileNameStr);
216}
217
218
219//_____________________________________________________________________________
220Bool_t AliSimulation::Run(Int_t nEvents)
221{
222// run the generation, simulation and digitization
223
224 if (nEvents > 0) fNEvents = nEvents;
225
226 // generation and simulation -> hits
227 if (fRunGeneration) {
85a5290f 228 if (!RunSimulation()) if (fStopOnError) return kFALSE;
229 }
230
85a5290f 231 // hits -> summable digits
232 if (!fMakeSDigits.IsNull()) {
233 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
234 }
235
236 // summable digits -> digits
237 if (!fMakeDigits.IsNull()) {
238 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
239 if (fStopOnError) return kFALSE;
240 }
241 }
242
243 // hits -> digits
244 if (!fMakeDigitsFromHits.IsNull()) {
95601147 245 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
85a5290f 246 Warning("Run", "Merging and direct creation of digits from hits "
247 "was selected for some detectors. "
248 "No merging will be done for the following detectors: %s",
249 fMakeDigitsFromHits.Data());
250 }
251 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
252 if (fStopOnError) return kFALSE;
253 }
254 }
255
0421c3d1 256 // digits -> raw data
257 if (!fWriteRawData.IsNull()) {
258 if (!WriteRawData(fWriteRawData)) {
259 if (fStopOnError) return kFALSE;
260 }
261 }
262
85a5290f 263 return kTRUE;
264}
265
266//_____________________________________________________________________________
95601147 267Bool_t AliSimulation::RunSimulation(Int_t nEvents)
85a5290f 268{
269// run the generation and simulation
270
4df28b43 271 TStopwatch stopwatch;
272 stopwatch.Start();
273
95601147 274 if (!gAlice) {
275 Error("RunSimulation", "no gAlice object. Restart aliroot and try again.");
276 return kFALSE;
277 }
278 if (gAlice->Modules()->GetEntries() > 0) {
279 Error("RunSimulation",
280 "gAlice was already run. Restart aliroot and try again.");
281 return kFALSE;
282 }
283
85a5290f 284 Info("RunSimulation", "initializing gAlice with config file %s",
285 fConfigFileName.Data());
286 gAlice->Init(fConfigFileName.Data());
95601147 287 AliRunLoader* runLoader = gAlice->GetRunLoader();
288 if (!runLoader) {
85a5290f 289 Error("RunSimulation", "gAlice has no run loader object. "
290 "Check your config file: %s", fConfigFileName.Data());
291 return kFALSE;
292 }
0421c3d1 293 SetGAliceFile(runLoader->GetFileName());
85a5290f 294
04bae0a0 295 if (!gAlice->Generator()) {
296 Error("RunSimulation", "gAlice has no generator object. "
297 "Check your config file: %s", fConfigFileName.Data());
298 return kFALSE;
299 }
05526d44 300 if (nEvents <= 0) nEvents = fNEvents;
04bae0a0 301
302 // get vertex from background file in case of merging
303 if (fUseBkgrdVertex &&
304 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
05526d44 305 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
04bae0a0 306 const char* fileName = ((TObjString*)
307 (fBkgrdFileNames->At(0)))->GetName();
308 Info("RunSimulation", "The vertex will be taken from the background "
309 "file %s with nSignalPerBackground = %d",
310 fileName, signalPerBkgrd);
311 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
312 gAlice->Generator()->SetVertexGenerator(vtxGen);
313 }
314
85a5290f 315 if (!fRunSimulation) {
85a5290f 316 gAlice->Generator()->SetTrackingFlag(0);
317 }
318
4df28b43 319 Info("RunSimulation", "running gAlice");
95601147 320 gAlice->Run(nEvents);
321
322 delete runLoader;
85a5290f 323
4df28b43 324 Info("RunSimulation", "execution time:");
325 stopwatch.Print();
326
85a5290f 327 return kTRUE;
328}
329
330//_____________________________________________________________________________
95601147 331Bool_t AliSimulation::RunSDigitization(const char* detectors)
85a5290f 332{
333// run the digitization and produce summable digits
334
4df28b43 335 TStopwatch stopwatch;
336 stopwatch.Start();
337
95601147 338 AliRunLoader* runLoader = LoadRun();
339 if (!runLoader) return kFALSE;
340
85a5290f 341 TString detStr = detectors;
95601147 342 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
85a5290f 343 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
344 AliModule* det = (AliModule*) detArray->At(iDet);
345 if (!det || !det->IsActive()) continue;
346 if (IsSelected(det->GetName(), detStr)) {
347 Info("RunSDigitization", "creating summable digits for %s",
348 det->GetName());
de76655b 349 TStopwatch stopwatchDet;
350 stopwatchDet.Start();
85a5290f 351 det->Hits2SDigits();
de76655b 352 Info("RunSDigitization", "execution time for %s:", det->GetName());
353 stopwatchDet.Print();
85a5290f 354 }
355 }
356
357 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
358 Error("RunSDigitization", "the following detectors were not found: %s",
359 detStr.Data());
360 if (fStopOnError) return kFALSE;
361 }
362
95601147 363 delete runLoader;
364
4df28b43 365 Info("RunSDigitization", "execution time:");
366 stopwatch.Print();
367
85a5290f 368 return kTRUE;
369}
370
371
372//_____________________________________________________________________________
95601147 373Bool_t AliSimulation::RunDigitization(const char* detectors,
374 const char* excludeDetectors)
85a5290f 375{
376// run the digitization and produce digits from sdigits
377
4df28b43 378 TStopwatch stopwatch;
379 stopwatch.Start();
380
95601147 381 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
382 if (gAlice) delete gAlice;
383 gAlice = NULL;
384
385 Int_t nStreams = 1;
05526d44 386 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
387 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
85a5290f 388 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
389 manager->SetInputStream(0, fGAliceFileName.Data());
390 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
391 const char* fileName = ((TObjString*)
392 (fBkgrdFileNames->At(iStream-1)))->GetName();
393 manager->SetInputStream(iStream, fileName);
394 }
395
396 TString detStr = detectors;
397 TString detExcl = excludeDetectors;
95601147 398 manager->GetInputStream(0)->ImportgAlice();
399 AliRunLoader* runLoader =
400 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
401 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
85a5290f 402 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
403 AliModule* det = (AliModule*) detArray->At(iDet);
404 if (!det || !det->IsActive()) continue;
405 if (IsSelected(det->GetName(), detStr) &&
406 !IsSelected(det->GetName(), detExcl)) {
4df28b43 407 AliDigitizer* digitizer = det->CreateDigitizer(manager);
408 if (!digitizer) {
85a5290f 409 Error("RunDigitization", "no digitizer for %s", det->GetName());
410 if (fStopOnError) return kFALSE;
4df28b43 411 } else {
412 digitizer->SetRegionOfInterest(fRegionOfInterest);
85a5290f 413 }
414 }
415 }
416
417 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
418 Error("RunDigitization", "the following detectors were not found: %s",
419 detStr.Data());
420 if (fStopOnError) return kFALSE;
421 }
422
423 if (!manager->GetListOfTasks()->IsEmpty()) {
424 Info("RunDigitization", "executing digitization");
425 manager->Exec("");
426 }
95601147 427
85a5290f 428 delete manager;
429
4df28b43 430 Info("RunDigitization", "execution time:");
431 stopwatch.Print();
432
85a5290f 433 return kTRUE;
434}
435
436//_____________________________________________________________________________
95601147 437Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
85a5290f 438{
439// run the digitization and produce digits from hits
440
4df28b43 441 TStopwatch stopwatch;
442 stopwatch.Start();
443
95601147 444 AliRunLoader* runLoader = LoadRun();
445 if (!runLoader) return kFALSE;
446
85a5290f 447 TString detStr = detectors;
95601147 448 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
85a5290f 449 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
450 AliModule* det = (AliModule*) detArray->At(iDet);
451 if (!det || !det->IsActive()) continue;
452 if (IsSelected(det->GetName(), detStr)) {
453 Info("RunHitsDigitization", "creating digits from hits for %s",
454 det->GetName());
455 det->Hits2Digits();
456 }
457 }
458
459 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
460 Error("RunHitsDigitization", "the following detectors were not found: %s",
461 detStr.Data());
462 if (fStopOnError) return kFALSE;
463 }
464
95601147 465 delete runLoader;
466
4df28b43 467 Info("RunHitsDigitization", "execution time:");
468 stopwatch.Print();
469
85a5290f 470 return kTRUE;
471}
472
0421c3d1 473//_____________________________________________________________________________
474Bool_t AliSimulation::WriteRawData(const char* detectors)
475{
476// convert the digits to raw data
477
478 TStopwatch stopwatch;
479 stopwatch.Start();
480
481 AliRunLoader* runLoader = LoadRun();
482 if (!runLoader) return kFALSE;
483
484 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
485 Info("WriteRawData", "processing event %d", iEvent);
486 runLoader->GetEvent(iEvent);
487 TString baseDir = gSystem->WorkingDirectory();
488 char dirName[256];
489 sprintf(dirName, "raw%d", iEvent);
490 gSystem->MakeDirectory(dirName);
491 if (!gSystem->ChangeDirectory(dirName)) {
492 Error("WriteRawData", "couldn't change to directory %s", dirName);
493 if (fStopOnError) return kFALSE; else continue;
494 }
495
496 TString detStr = detectors;
497 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
498 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
499 AliModule* det = (AliModule*) detArray->At(iDet);
500 if (!det || !det->IsActive()) continue;
501 if (IsSelected(det->GetName(), detStr)) {
502 Info("WriteRawData", "creating raw data from digits for %s",
503 det->GetName());
504 det->Digits2Raw();
505 }
506 }
507
508 gSystem->ChangeDirectory(baseDir);
509 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
510 Error("WriteRawData", "the following detectors were not found: %s",
511 detStr.Data());
512 if (fStopOnError) return kFALSE;
513 }
514 }
515
516 delete runLoader;
517
518 Info("WriteRawData", "execution time:");
519 stopwatch.Print();
520
521 return kTRUE;
522}
523
85a5290f 524
95601147 525//_____________________________________________________________________________
526AliRunLoader* AliSimulation::LoadRun() const
527{
528// delete existing run loaders, open a new one and load gAlice
529
530 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
531 AliRunLoader* runLoader =
532 AliRunLoader::Open(fGAliceFileName.Data(),
e191bb57 533 AliConfig::GetDefaultEventFolderName(), "UPDATE");
95601147 534 if (!runLoader) {
535 Error("LoadRun", "no run loader found in file %s",
536 fGAliceFileName.Data());
537 return NULL;
538 }
539 runLoader->LoadgAlice();
540 gAlice = runLoader->GetAliRun();
541 if (!gAlice) {
542 Error("LoadRun", "no gAlice object found in file %s",
543 fGAliceFileName.Data());
544 return NULL;
545 }
546 return runLoader;
547}
548
85a5290f 549//_____________________________________________________________________________
05526d44 550Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
551{
552// get or calculate the number of signal events per background event
553
554 if (!fBkgrdFileNames) return 1;
555 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
556 if (nBkgrdFiles == 0) return 1;
557
558 // get the number of signal events
559 if (nEvents <= 0) {
560 AliRunLoader* runLoader =
561 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
562 if (!runLoader) return 1;
563 nEvents = runLoader->GetNumberOfEvents();
564 delete runLoader;
565 }
566
567 Int_t result = 0;
568 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
569 // get the number of background events
570 const char* fileName = ((TObjString*)
571 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
572 AliRunLoader* runLoader =
573 AliRunLoader::Open(fileName, "BKGRD");
574 if (!runLoader) continue;
575 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
576 delete runLoader;
577
578 // get or calculate the number of signal per background events
579 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
580 if (nSignalPerBkgrd <= 0) {
581 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
582 } else if (result && (result != nSignalPerBkgrd)) {
583 Info("GetNSignalPerBkgrd", "the number of signal events per "
584 "background event will be changed from %d to %d for stream %d",
585 nSignalPerBkgrd, result, iBkgrdFile+1);
586 nSignalPerBkgrd = result;
587 }
588
589 if (!result) result = nSignalPerBkgrd;
590 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
591 Warning("GetNSignalPerBkgrd", "not enough background events (%d) for "
592 "%d signal events using %d signal per background events for "
593 "stream %d",
594 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1);
595 }
596 }
597
598 return result;
599}
600
601//_____________________________________________________________________________
85a5290f 602Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
603{
604// check whether detName is contained in detectors
605// if yes, it is removed from detectors
606
607 // check if all detectors are selected
608 if ((detectors.CompareTo("ALL") == 0) ||
609 detectors.BeginsWith("ALL ") ||
610 detectors.EndsWith(" ALL") ||
611 detectors.Contains(" ALL ")) {
612 detectors = "ALL";
613 return kTRUE;
614 }
615
616 // search for the given detector
617 Bool_t result = kFALSE;
618 if ((detectors.CompareTo(detName) == 0) ||
619 detectors.BeginsWith(detName+" ") ||
620 detectors.EndsWith(" "+detName) ||
621 detectors.Contains(" "+detName+" ")) {
622 detectors.ReplaceAll(detName, "");
623 result = kTRUE;
624 }
625
626 // clean up the detectors string
627 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
628 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
629 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
630
631 return result;
632}