]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
Correction in the memory managemen. Adding new getter
[u/mrichter/AliRoot.git] / STEER / AliSimulation.cxx
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 passed as argument to the       //
34 // AliSimulation constructor or can be specified by                          //
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 //                                                                           //
64 // Backgound events can be merged by calling                                 //
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.      //
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.                          //
73 //                                                                           //
74 // The methods RunSimulation, RunSDigitization, RunDigitization and          //
75 // RunHitsDigitization can be used to run only parts of the full simulation  //
76 // chain.                                                                    //
77 //                                                                           //
78 ///////////////////////////////////////////////////////////////////////////////
79
80
81 #include "AliSimulation.h"
82 #include "AliRunLoader.h"
83 #include "AliRun.h"
84 #include "AliModule.h"
85 #include "AliGenerator.h"
86 #include "AliRunDigitizer.h"
87 #include "AliDigitizer.h"
88 #include <TObjString.h>
89
90
91 ClassImp(AliSimulation)
92
93
94 //_____________________________________________________________________________
95 AliSimulation::AliSimulation(const char* configFileName,
96                              const char* name, const char* title) :
97   TNamed(name, title),
98
99   fRunGeneration(kTRUE),
100   fRunSimulation(kTRUE),
101   fMakeSDigits("ALL"),
102   fMakeDigits("ALL"),
103   fMakeDigitsFromHits(""),
104   fStopOnError(kFALSE),
105
106   fNEvents(1),
107   fConfigFileName(configFileName),
108   fGAliceFileName("galice.root"),
109   fBkgrdFileNames(NULL),
110   fRegionOfInterest(kTRUE)
111 {
112 // create simulation object with default parameters
113
114 }
115
116 //_____________________________________________________________________________
117 AliSimulation::AliSimulation(const AliSimulation& sim) :
118   TNamed(sim),
119
120   fRunGeneration(sim.fRunGeneration),
121   fRunSimulation(sim.fRunSimulation),
122   fMakeSDigits(sim.fMakeSDigits),
123   fMakeDigits(sim.fMakeDigits),
124   fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
125   fStopOnError(sim.fStopOnError),
126
127   fNEvents(sim.fNEvents),
128   fConfigFileName(sim.fConfigFileName),
129   fGAliceFileName(sim.fGAliceFileName),
130   fBkgrdFileNames(NULL),
131   fRegionOfInterest(sim.fRegionOfInterest)
132 {
133 // copy constructor
134
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());
139   }
140 }
141
142 //_____________________________________________________________________________
143 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
144 {
145 // assignment operator
146
147   this->~AliSimulation();
148   new(this) AliSimulation(sim);
149   return *this;
150 }
151
152 //_____________________________________________________________________________
153 AliSimulation::~AliSimulation()
154 {
155 // clean up
156
157   if (fBkgrdFileNames) {
158     fBkgrdFileNames->Delete();
159     delete fBkgrdFileNames;
160   }
161 }
162
163
164 //_____________________________________________________________________________
165 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
166 {
167 // set the number of events for one run
168
169   fNEvents = nEvents;
170 }
171
172 //_____________________________________________________________________________
173 void AliSimulation::SetConfigFile(const char* fileName)
174 {
175 // set the name of the config file
176
177   fConfigFileName = fileName;
178 }
179
180 //_____________________________________________________________________________
181 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
182 {
183 // add a file with background events for merging
184
185   TObjString* fileNameStr = new TObjString(fileName);
186   fileNameStr->SetUniqueID(nSignalPerBkgrd);
187   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
188   fBkgrdFileNames->Add(fileNameStr);
189 }
190
191
192 //_____________________________________________________________________________
193 Bool_t AliSimulation::Run(Int_t nEvents)
194 {
195 // run the generation, simulation and digitization
196
197   if (nEvents > 0) fNEvents = nEvents;
198
199   // generation and simulation -> hits
200   if (fRunGeneration) {
201     if (!RunSimulation()) if (fStopOnError) return kFALSE;
202   }
203
204   // hits -> summable digits
205   if (!fMakeSDigits.IsNull()) {
206     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
207   }
208
209   // summable digits -> digits
210   if (!fMakeDigits.IsNull()) {
211     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
212       if (fStopOnError) return kFALSE;
213     }
214   }
215
216   // hits -> digits
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());
223     }
224     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
225       if (fStopOnError) return kFALSE;
226     }
227   }
228
229   return kTRUE;
230 }
231
232 //_____________________________________________________________________________
233 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
234 {
235 // run the generation and simulation
236
237   TStopwatch stopwatch;
238   stopwatch.Start();
239
240   if (!gAlice) {
241     Error("RunSimulation", "no gAlice object. Restart aliroot and try again.");
242     return kFALSE;
243   }
244   if (gAlice->Modules()->GetEntries() > 0) {
245     Error("RunSimulation", 
246           "gAlice was already run. Restart aliroot and try again.");
247     return kFALSE;
248   }
249
250   Info("RunSimulation", "initializing gAlice with config file %s",
251        fConfigFileName.Data());
252   gAlice->Init(fConfigFileName.Data());
253   AliRunLoader* runLoader = gAlice->GetRunLoader();
254   if (!runLoader) {
255     Error("RunSimulation", "gAlice has no run loader object. "
256           "Check your config file: %s", fConfigFileName.Data());
257     return kFALSE;
258   }
259   fGAliceFileName = runLoader->GetFileName();
260
261   if (!fRunSimulation) {
262     if (!gAlice->Generator()) {
263       Error("RunSimulation", "gAlice has no generator object. "
264             "Check your config file: %s", fConfigFileName.Data());
265       return kFALSE;
266     }
267     gAlice->Generator()->SetTrackingFlag(0);
268   }
269
270   Info("RunSimulation", "running gAlice");
271   if (nEvents <= 0) nEvents = fNEvents;
272   gAlice->Run(nEvents);
273
274   delete runLoader;
275
276   Info("RunSimulation", "execution time:");
277   stopwatch.Print();
278
279   return kTRUE;
280 }
281
282 //_____________________________________________________________________________
283 Bool_t AliSimulation::RunSDigitization(const char* detectors)
284 {
285 // run the digitization and produce summable digits
286
287   TStopwatch stopwatch;
288   stopwatch.Start();
289
290   AliRunLoader* runLoader = LoadRun();
291   if (!runLoader) return kFALSE;
292
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", 
300            det->GetName());
301       det->Hits2SDigits();
302     }
303   }
304
305   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
306     Error("RunSDigitization", "the following detectors were not found: %s", 
307           detStr.Data());
308     if (fStopOnError) return kFALSE;
309   }
310
311   delete runLoader;
312
313   Info("RunSDigitization", "execution time:");
314   stopwatch.Print();
315
316   return kTRUE;
317 }
318
319
320 //_____________________________________________________________________________
321 Bool_t AliSimulation::RunDigitization(const char* detectors, 
322                                       const char* excludeDetectors)
323 {
324 // run the digitization and produce digits from sdigits
325
326   TStopwatch stopwatch;
327   stopwatch.Start();
328
329   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
330   if (gAlice) delete gAlice;
331   gAlice = NULL;
332
333   Int_t nStreams = 1;
334   Int_t signalPerBkgrd = 1;
335   if (fBkgrdFileNames) {
336     nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
337     if (nStreams > 1) signalPerBkgrd = fBkgrdFileNames->At(0)->GetUniqueID();
338   }
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);
345   }
346
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);
359       if (!digitizer) {
360         Error("RunDigitization", "no digitizer for %s", det->GetName());
361         if (fStopOnError) return kFALSE;
362       } else {
363         digitizer->SetRegionOfInterest(fRegionOfInterest);
364       }
365     }
366   }
367
368   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
369     Error("RunDigitization", "the following detectors were not found: %s", 
370           detStr.Data());
371     if (fStopOnError) return kFALSE;
372   }
373
374   if (!manager->GetListOfTasks()->IsEmpty()) {
375     Info("RunDigitization", "executing digitization");
376     manager->Exec("");
377   }
378
379   delete manager;
380
381   Info("RunDigitization", "execution time:");
382   stopwatch.Print();
383
384   return kTRUE;
385 }
386
387 //_____________________________________________________________________________
388 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
389 {
390 // run the digitization and produce digits from hits
391
392   TStopwatch stopwatch;
393   stopwatch.Start();
394
395   AliRunLoader* runLoader = LoadRun();
396   if (!runLoader) return kFALSE;
397
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", 
405            det->GetName());
406       det->Hits2Digits();
407     }
408   }
409
410   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
411     Error("RunHitsDigitization", "the following detectors were not found: %s", 
412           detStr.Data());
413     if (fStopOnError) return kFALSE;
414   }
415
416   delete runLoader;
417
418   Info("RunHitsDigitization", "execution time:");
419   stopwatch.Print();
420
421   return kTRUE;
422 }
423
424
425 //_____________________________________________________________________________
426 AliRunLoader* AliSimulation::LoadRun() const
427 {
428 // delete existing run loaders, open a new one and load gAlice
429
430   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
431   AliRunLoader* runLoader = 
432     AliRunLoader::Open(fGAliceFileName.Data(), 
433                        AliConfig::fgkDefaultEventFolderName, "UPDATE");
434   if (!runLoader) {
435     Error("LoadRun", "no run loader found in file %s", 
436           fGAliceFileName.Data());
437     return NULL;
438   }
439   runLoader->LoadgAlice();
440   gAlice = runLoader->GetAliRun();
441   if (!gAlice) {
442     Error("LoadRun", "no gAlice object found in file %s", 
443           fGAliceFileName.Data());
444     return NULL;
445   }
446   return runLoader;
447 }
448
449 //_____________________________________________________________________________
450 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
451 {
452 // check whether detName is contained in detectors
453 // if yes, it is removed from detectors
454
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 ")) {
460     detectors = "ALL";
461     return kTRUE;
462   }
463
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, "");
471     result = kTRUE;
472   }
473
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);
478
479   return result;
480 }