]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
Remove old geometries
[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 // Background 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 // 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.                                       //
74 //                                                                           //
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 default output format of the raw data are DDL files. They are         //
80 // converted to a DATE file, if a file name is given as second argument.     //
81 // For this conversion the program "dateStream" is required. If the file     //
82 // name has the extension ".root", the DATE file is converted to a root      //
83 // file. The program "alimdc" is used for this purpose. For the conversion   //
84 // to DATE and root format the two conversion programs have to be installed. //
85 // Only the raw data in the final format is kept if the third argument is    //
86 // kTRUE.                                                                    //
87 //                                                                           //
88 // The methods RunSimulation, RunSDigitization, RunDigitization,             //
89 // RunHitsDigitization and WriteRawData can be used to run only parts of     //
90 // the full simulation chain. The creation of raw data DDL files and their   //
91 // conversion to the DATE or root format can be run directly by calling      //
92 // the methods WriteRawFiles, ConvertRawFilesToDate and ConvertDateToRoot.   //
93 //                                                                           //
94 // The default number of events per file, which is usually set in the        //
95 // config file, can be changed for individual detectors and data types       //
96 // by calling                                                                //
97 //                                                                           //
98 //   sim.SetEventsPerFile("PHOS", "Reconstructed Points", 3);                //
99 //                                                                           //
100 // The first argument is the detector, the second one the data type and the  //
101 // last one the number of events per file. Valid data types are "Hits",      //
102 // "Summable Digits", "Digits", "Reconstructed Points" and "Tracks".         //
103 // The number of events per file has to be set before the simulation of      //
104 // hits. Otherwise it has no effect.                                         //
105 //                                                                           //
106 ///////////////////////////////////////////////////////////////////////////////
107
108 #include <TGeoManager.h>
109 #include <TObjString.h>
110 #include <TStopwatch.h>
111 #include <TSystem.h>
112 #include <TFile.h>
113
114 #include "AliCDBStorage.h"
115 #include "AliCDBEntry.h"
116 #include "AliCDBManager.h"
117 #include "AliAlignObj.h"
118 #include "AliCentralTrigger.h"
119 #include "AliDAQConfig.h"
120 #include "AliDigitizer.h"
121 #include "AliGenerator.h"
122 #include "AliLog.h"
123 #include "AliModule.h"
124 #include "AliRun.h"
125 #include "AliRunDigitizer.h"
126 #include "AliRunLoader.h"
127 #include "AliSimulation.h"
128 #include "AliVertexGenFile.h"
129 #include "AliCentralTrigger.h"
130 #include "AliCTPRawData.h"
131
132 ClassImp(AliSimulation)
133
134
135 //_____________________________________________________________________________
136 AliSimulation::AliSimulation(const char* configFileName, const char* cdbUri,
137                              const char* name, const char* title) :
138   TNamed(name, title),
139
140   fRunGeneration(kTRUE),
141   fRunSimulation(kTRUE),
142   fLoadAlignFromCDB(kTRUE),
143   fLoadAlignData("ALL"),
144   fMakeSDigits("ALL"),
145   fMakeDigits("ALL"),
146   fMakeTrigger(""),
147   fMakeDigitsFromHits(""),
148   fWriteRawData(""),
149   fRawDataFileName(""),
150   fDeleteIntermediateFiles(kFALSE),
151   fStopOnError(kFALSE),
152
153   fNEvents(1),
154   fConfigFileName(configFileName),
155   fGAliceFileName("galice.root"),
156   fEventsPerFile(),
157   fBkgrdFileNames(NULL),
158   fAlignObjArray(NULL),
159   fUseBkgrdVertex(kTRUE),
160   fRegionOfInterest(kFALSE),
161   fCDBUri(cdbUri)
162 {
163 // create simulation object with default parameters
164
165   SetGAliceFile("galice.root");
166 }
167
168 //_____________________________________________________________________________
169 AliSimulation::AliSimulation(const AliSimulation& sim) :
170   TNamed(sim),
171
172   fRunGeneration(sim.fRunGeneration),
173   fRunSimulation(sim.fRunSimulation),
174   fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
175   fLoadAlignData(sim.fLoadAlignData),
176   fMakeSDigits(sim.fMakeSDigits),
177   fMakeDigits(sim.fMakeDigits),
178   fMakeTrigger(sim.fMakeTrigger),
179   fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
180   fWriteRawData(sim.fWriteRawData),
181   fRawDataFileName(""),
182   fDeleteIntermediateFiles(kFALSE),
183   fStopOnError(sim.fStopOnError),
184
185   fNEvents(sim.fNEvents),
186   fConfigFileName(sim.fConfigFileName),
187   fGAliceFileName(sim.fGAliceFileName),
188   fEventsPerFile(),
189   fBkgrdFileNames(NULL),
190   fAlignObjArray(NULL),
191   fUseBkgrdVertex(sim.fUseBkgrdVertex),
192   fRegionOfInterest(sim.fRegionOfInterest),
193   fCDBUri(sim.fCDBUri)
194 {
195 // copy constructor
196
197   for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
198     if (!sim.fEventsPerFile[i]) continue;
199     fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
200   }
201
202   fBkgrdFileNames = new TObjArray;
203   for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
204     if (!sim.fBkgrdFileNames->At(i)) continue;
205     fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
206   }
207 }
208
209 //_____________________________________________________________________________
210 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
211 {
212 // assignment operator
213
214   this->~AliSimulation();
215   new(this) AliSimulation(sim);
216   return *this;
217 }
218
219 //_____________________________________________________________________________
220 AliSimulation::~AliSimulation()
221 {
222 // clean up
223
224   fEventsPerFile.Delete();
225 //  if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
226 //  delete fAlignObjArray; fAlignObjArray=0;
227
228   if (fBkgrdFileNames) {
229     fBkgrdFileNames->Delete();
230     delete fBkgrdFileNames;
231   }
232 }
233
234
235 //_____________________________________________________________________________
236 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
237 {
238 // set the number of events for one run
239
240   fNEvents = nEvents;
241 }
242
243 //_____________________________________________________________________________
244 void AliSimulation::InitCDBStorage(const char* uri)
245 {
246 // activate a default CDB storage
247 // First check if we have any CDB storage set, because it is used 
248 // to retrieve the calibration and alignment constants
249
250   AliCDBManager* man = AliCDBManager::Instance();
251   if (!man->IsDefaultStorageSet())
252   {
253     AliWarningClass("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
254     AliWarningClass("Default CDB storage not yet set");
255     AliWarningClass(Form("Using default storage declared in AliSimulation: %s",uri));
256     AliWarningClass("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
257     SetDefaultStorage(uri);
258   }  
259   
260 }
261
262 //_____________________________________________________________________________
263 void AliSimulation::SetDefaultStorage(const char* uri) {
264 // activate a default CDB storage 
265
266    AliCDBManager::Instance()->SetDefaultStorage(uri);
267
268 }
269
270 //_____________________________________________________________________________
271 void AliSimulation::SetSpecificStorage(const char* detName, const char* uri) {
272 // activate a detector-specific CDB storage 
273
274    AliCDBManager::Instance()->SetSpecificStorage(detName, uri);
275
276 }
277
278 //_____________________________________________________________________________
279 void AliSimulation::SetConfigFile(const char* fileName)
280 {
281 // set the name of the config file
282
283   fConfigFileName = fileName;
284 }
285
286 //_____________________________________________________________________________
287 void AliSimulation::SetGAliceFile(const char* fileName)
288 {
289 // set the name of the galice file
290 // the path is converted to an absolute one if it is relative
291
292   fGAliceFileName = fileName;
293   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
294     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
295                                                 fGAliceFileName);
296     fGAliceFileName = absFileName;
297     delete[] absFileName;
298   }
299
300   AliDebug(2, Form("galice file name set to %s", fileName));
301 }
302
303 //_____________________________________________________________________________
304 void AliSimulation::SetEventsPerFile(const char* detector, const char* type, 
305                                      Int_t nEvents)
306 {
307 // set the number of events per file for the given detector and data type
308 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
309
310   TNamed* obj = new TNamed(detector, type);
311   obj->SetUniqueID(nEvents);
312   fEventsPerFile.Add(obj);
313 }
314
315 //_____________________________________________________________________________
316 Bool_t AliSimulation::ApplyAlignObjsToGeom(TObjArray* alObjArray)
317 {
318   // Read collection of alignment objects (AliAlignObj derived) saved
319   // in the TClonesArray ClArrayName and apply them to the geometry
320   // manager singleton.
321   //
322   alObjArray->Sort();
323   Int_t nvols = alObjArray->GetEntriesFast();
324
325   Bool_t flag = kTRUE;
326
327   for(Int_t j=0; j<nvols; j++)
328     {
329       AliAlignObj* alobj = (AliAlignObj*) alObjArray->UncheckedAt(j);
330       if (alobj->ApplyToGeometry() == kFALSE) flag = kFALSE;
331     }
332
333   if (AliDebugLevelClass() >= 1) {
334     gGeoManager->CheckOverlaps(20);
335     TObjArray* ovexlist = gGeoManager->GetListOfOverlaps();
336     if(ovexlist->GetEntriesFast()){  
337       AliError("The application of alignment objects to the geometry caused huge overlaps/extrusions!");
338    }
339   }
340
341   return flag;
342
343 }
344
345 //_____________________________________________________________________________
346 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* fileName, const char* clArrayName)
347 {
348   // read collection of alignment objects (AliAlignObj derived) saved
349   // in the TClonesArray ClArrayName in the file fileName and apply
350   // them to the TGeo geometry passed as argument
351   //
352
353   TFile* inFile = TFile::Open(fileName,"READ");
354   if (!inFile || !inFile->IsOpen()) {
355     AliErrorClass(Form("Could not open file %s !",fileName));
356     return kFALSE;
357   }
358
359   TClonesArray* alObjArray = ((TClonesArray*) inFile->Get(clArrayName));
360   inFile->Close();
361   if (!alObjArray) {
362     AliErrorClass(Form("Could not get array (%s) from file (%s) !",clArrayName,fileName));
363     return kFALSE;
364   }
365
366   return ApplyAlignObjsToGeom(alObjArray);
367
368 }
369
370 //_____________________________________________________________________________
371 Bool_t AliSimulation::ApplyAlignObjsToGeom(AliCDBParam* param, AliCDBId& Id)
372 {
373   // read collection of alignment objects (AliAlignObj derived) saved
374   // in the TClonesArray ClArrayName in the AliCDBEntry identified by
375   // param (to get the AliCDBStorage) and Id; apply the alignment objects
376   // to the TGeo geometry passed as argument
377   //
378
379   AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage(param);
380   AliCDBEntry* entry = storage->Get(Id);
381   TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
382
383   return ApplyAlignObjsToGeom(AlObjArray);
384
385 }
386
387 //_____________________________________________________________________________
388 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* uri, const char* path, Int_t runnum, Int_t version, Int_t sversion)
389 {
390   // read collection of alignment objects (AliAlignObj derived) saved
391   // in the TClonesArray ClArrayName in the AliCDBEntry identified by
392   // param (to get the AliCDBStorage) and Id; apply the alignment objects
393   // to the TGeo geometry passed as argument
394   //
395
396   AliCDBParam* param = AliCDBManager::Instance()->CreateParameter(uri);
397   AliCDBId id(path, runnum, runnum, version, sversion);
398
399   return ApplyAlignObjsToGeom(param, id);
400
401 }
402
403 //_____________________________________________________________________________
404 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* detName, Int_t runnum, Int_t version, Int_t sversion)
405 {
406   // read collection of alignment objects (AliAlignObj derived) saved
407   // in the TClonesArray ClArrayName in the AliCDBEntry identified by
408   // param (to get the AliCDBStorage) and Id; apply the alignment objects
409   // to the TGeo geometry passed as argument
410   //
411
412   InitCDBStorage("local://$ALICE_ROOT");
413   AliCDBPath path(detName,"Align","Data");
414   AliCDBEntry* entry = AliCDBManager::Instance()->Get(path.GetPath(),runnum,version,sversion);
415
416   if(!entry) return kFALSE;
417   TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
418
419   return ApplyAlignObjsToGeom(AlObjArray);
420 }
421
422 //_____________________________________________________________________________
423 Bool_t AliSimulation::SetAlignObjArraySingleDet(const char* detName)
424 {
425   // Fills array of single detector's alignable objects from CDB
426   
427   AliDebug(2, Form("Loading alignment data for detector: %s",detName));
428   
429   AliCDBEntry *entry;
430         
431   AliCDBPath path(detName,"Align","Data");
432         
433   entry=AliCDBManager::Instance()->Get(path.GetPath());
434   if(!entry){ 
435         AliDebug(2,Form("Couldn't load alignment data for detector %s",detName));
436         return kFALSE;
437   }
438   entry->SetOwner(1);
439   TClonesArray *alignArray = (TClonesArray*) entry->GetObject();        
440   alignArray->SetOwner(0);
441   AliDebug(2,Form("Found %d alignment objects for %s",
442                         alignArray->GetEntries(),detName));
443
444   AliAlignObj *alignObj=0;
445   TIter iter(alignArray);
446         
447   // loop over align objects in detector
448   while( ( alignObj=(AliAlignObj *) iter.Next() ) ){
449         fAlignObjArray->Add(alignObj);
450   }
451   // delete entry --- Don't delete, it is cached!
452         
453   AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
454   return kTRUE;
455
456 }
457
458 //_____________________________________________________________________________
459 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
460 {
461   // Read the alignment objects from CDB.
462   // Each detector is supposed to have the
463   // alignment objects in DET/Align/Data CDB path.
464   // All the detector objects are then collected,
465   // sorted by geometry level (starting from ALIC) and
466   // then applied to the TGeo geometry.
467   // Finally an overlaps check is performed.
468
469   Bool_t delRunLoader = kFALSE;
470   if (!runLoader) {
471     runLoader = LoadRun("READ");
472     if (!runLoader) return kFALSE;
473     delRunLoader = kTRUE;
474   }
475
476   // Load alignment data from CDB and fill fAlignObjArray 
477   if(fLoadAlignFromCDB){
478         if(!fAlignObjArray) fAlignObjArray = new TObjArray();
479         
480         //fAlignObjArray->RemoveAll(); 
481         fAlignObjArray->Clear();        
482         fAlignObjArray->SetOwner(0);
483  
484         TString detStr = fLoadAlignData;
485         TString dataNotLoaded="";
486         TString dataLoaded="";
487   
488         TObjArray* detArray = runLoader->GetAliRun()->Detectors();
489         for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
490                 AliModule* det = (AliModule*) detArray->At(iDet);
491                 if (!det || !det->IsActive()) continue;
492                 if (IsSelected(det->GetName(), detStr)) {
493                         if(!SetAlignObjArraySingleDet(det->GetName())){
494                                 dataNotLoaded += det->GetName();
495                                 dataNotLoaded += " ";
496                         } else {
497                                 dataLoaded += det->GetName();
498                                 dataLoaded += " ";
499                         }
500                 }
501         } // end loop over detectors
502   
503         if ((detStr.CompareTo("ALL") == 0)) detStr = "";
504         dataNotLoaded += detStr;
505         AliInfo(Form("Alignment data loaded for: %s",
506                           dataLoaded.Data()));
507         AliInfo(Form("Didn't/couldn't load alignment data for: %s",
508                           dataNotLoaded.Data()));
509   } // fLoadAlignFromCDB flag
510  
511   // Check if the array with alignment objects was
512   // provided by the user. If yes, apply the objects
513   // to the present TGeo geometry
514   if (fAlignObjArray) {
515     if (gGeoManager && gGeoManager->IsClosed()) {
516       if (ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
517         AliError("The misalignment of one or more volumes failed!"
518                  "Compare the list of simulated detectors and the list of detector alignment data!");
519         if (delRunLoader) delete runLoader;
520         return kFALSE;
521       }
522     }
523     else {
524       AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
525       if (delRunLoader) delete runLoader;
526       return kFALSE;
527     }
528   }
529
530   if (delRunLoader) delete runLoader;
531
532   return kTRUE;
533 }
534
535
536 //_____________________________________________________________________________
537 Bool_t AliSimulation::SetRunNumber()
538 {
539   // Set the CDB manager run number
540   // The run number is retrieved from gAlice
541
542   if(AliCDBManager::Instance()->GetRun() < 0) {
543     AliRunLoader* runLoader = LoadRun("READ");
544     if (!runLoader) return kFALSE;
545     else {
546       AliCDBManager::Instance()->SetRun(runLoader->GetAliRun()->GetRunNumber());
547       AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
548       delete runLoader;
549     }
550   }
551   return kTRUE;
552 }
553
554 //_____________________________________________________________________________
555 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
556 {
557 // add a file with background events for merging
558
559   TObjString* fileNameStr = new TObjString(fileName);
560   fileNameStr->SetUniqueID(nSignalPerBkgrd);
561   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
562   fBkgrdFileNames->Add(fileNameStr);
563 }
564
565
566 //_____________________________________________________________________________
567 Bool_t AliSimulation::Run(Int_t nEvents)
568 {
569 // run the generation, simulation and digitization
570
571   InitCDBStorage(fCDBUri);
572
573   if (nEvents > 0) fNEvents = nEvents;
574
575   // generation and simulation -> hits
576   if (fRunGeneration) {
577     if (!RunSimulation()) if (fStopOnError) return kFALSE;
578   }
579
580   // Set run number in CDBManager (if it is not already set in RunSimulation)
581   if (!SetRunNumber()) if (fStopOnError) return kFALSE;
582
583   // Load and misalign the geometry
584   if (!gGeoManager) {
585     TGeoManager::Import("geometry.root");
586     if (!gGeoManager) if (fStopOnError) return kFALSE;
587     if (!MisalignGeometry()) if (fStopOnError) return kFALSE;
588   }
589
590   // hits -> summable digits
591   if (!fMakeSDigits.IsNull()) {
592     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
593   }
594
595   // summable digits -> digits
596   if (!fMakeDigits.IsNull()) {
597     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
598       if (fStopOnError) return kFALSE;
599     }
600   }
601
602   // hits -> digits
603   if (!fMakeDigitsFromHits.IsNull()) {
604     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
605       AliWarning(Form("Merging and direct creation of digits from hits " 
606                  "was selected for some detectors. "
607                  "No merging will be done for the following detectors: %s",
608                  fMakeDigitsFromHits.Data()));
609     }
610     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
611       if (fStopOnError) return kFALSE;
612     }
613   }
614
615   // digits -> trigger
616   if (!RunTrigger(fMakeTrigger)) {
617     if (fStopOnError) return kFALSE;
618   }
619
620   // digits -> raw data
621   if (!fWriteRawData.IsNull()) {
622     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
623                       fDeleteIntermediateFiles)) {
624       if (fStopOnError) return kFALSE;
625     }
626   }
627
628   return kTRUE;
629 }
630
631 //_____________________________________________________________________________
632 Bool_t AliSimulation::RunTrigger(const char* descriptors)
633 {
634   // run the trigger
635
636    TStopwatch stopwatch;
637    stopwatch.Start();
638
639    AliRunLoader* runLoader = LoadRun("READ");
640    if (!runLoader) return kFALSE;
641    TString des = descriptors;
642
643    if (des.IsNull()) {
644      if (gAlice->GetTriggerDescriptor() != "") {
645        des = gAlice->GetTriggerDescriptor();
646      }
647      else {
648        AliWarning("No trigger descriptor is specified. Skipping the trigger simulation...");
649        return kTRUE;
650      }
651    }
652
653    runLoader->MakeTree( "CT" );
654    AliCentralTrigger* aCTP = runLoader->GetTrigger();
655   // Load Descriptors
656    aCTP->LoadDescriptor( des );
657
658   // digits -> trigger
659    if( !aCTP->RunTrigger( runLoader ) ) {
660       if (fStopOnError) {
661     //  delete aCTP;
662          return kFALSE;
663       }
664    }
665
666    AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
667            stopwatch.RealTime(),stopwatch.CpuTime()));
668
669    delete runLoader;
670
671    return kTRUE;
672 }
673
674 //_____________________________________________________________________________
675 Bool_t AliSimulation::WriteTriggerRawData()
676 {
677   // Writes the CTP (trigger) DDL raw data
678   // Details of the format are given in the
679   // trigger TDR - pages 134 and 135.
680   AliCTPRawData writer;
681   writer.RawData();
682
683   return kTRUE;
684 }
685
686 //_____________________________________________________________________________
687 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
688 {
689 // run the generation and simulation
690
691   TStopwatch stopwatch;
692   stopwatch.Start();
693
694   if (!gAlice) {
695     AliError("no gAlice object. Restart aliroot and try again.");
696     return kFALSE;
697   }
698   if (gAlice->Modules()->GetEntries() > 0) {
699     AliError("gAlice was already run. Restart aliroot and try again.");
700     return kFALSE;
701   }
702
703   AliInfo(Form("initializing gAlice with config file %s",
704           fConfigFileName.Data()));
705   StdoutToAliInfo(StderrToAliError(
706     gAlice->Init(fConfigFileName.Data());
707   ););
708
709   // Get the trigger descriptor string
710   // Either from AliSimulation or from
711   // gAlice
712   if (fMakeTrigger.IsNull()) {
713     if (gAlice->GetTriggerDescriptor() != "")
714       fMakeTrigger = gAlice->GetTriggerDescriptor();
715   }
716   else
717     gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
718
719   // Set run number in CDBManager
720   AliCDBManager::Instance()->SetRun(gAlice->GetRunNumber());
721   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
722
723   AliRunLoader* runLoader = gAlice->GetRunLoader();
724   if (!runLoader) {
725              AliError(Form("gAlice has no run loader object. "
726                              "Check your config file: %s", fConfigFileName.Data()));
727              return kFALSE;
728   }
729   SetGAliceFile(runLoader->GetFileName());
730  
731   // Export ideal geometry 
732   if (gGeoManager) gGeoManager->Export("geometry.root");
733
734   // Misalign geometry
735   if (!MisalignGeometry(runLoader)) return kFALSE;
736
737   // Temporary fix by A.Gheata
738   // Could be removed with the next Root version (>5.11)
739   if (gGeoManager) {
740     TIter next(gGeoManager->GetListOfVolumes());
741     TGeoVolume *vol;
742     while ((vol = (TGeoVolume *)next())) {
743       if (vol->GetVoxels()) {
744         if (vol->GetVoxels()->NeedRebuild()) {
745           vol->GetVoxels()->Voxelize();
746           vol->FindOverlaps();
747         }
748       }
749     }
750   }
751
752   // Export (mis)aligned geometry 
753   if (gGeoManager) gGeoManager->Export("misaligned_geometry.root");
754
755 //   AliRunLoader* runLoader = gAlice->GetRunLoader();
756 //   if (!runLoader) {
757 //     AliError(Form("gAlice has no run loader object. "
758 //                   "Check your config file: %s", fConfigFileName.Data()));
759 //     return kFALSE;
760 //   }
761 //   SetGAliceFile(runLoader->GetFileName());
762
763   if (!gAlice->Generator()) {
764     AliError(Form("gAlice has no generator object. "
765                   "Check your config file: %s", fConfigFileName.Data()));
766     return kFALSE;
767   }
768   if (nEvents <= 0) nEvents = fNEvents;
769
770   // get vertex from background file in case of merging
771   if (fUseBkgrdVertex &&
772       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
773     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
774     const char* fileName = ((TObjString*)
775                             (fBkgrdFileNames->At(0)))->GetName();
776     AliInfo(Form("The vertex will be taken from the background "
777                  "file %s with nSignalPerBackground = %d", 
778                  fileName, signalPerBkgrd));
779     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
780     gAlice->Generator()->SetVertexGenerator(vtxGen);
781   }
782
783   if (!fRunSimulation) {
784     gAlice->Generator()->SetTrackingFlag(0);
785   }
786
787   // set the number of events per file for given detectors and data types
788   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
789     if (!fEventsPerFile[i]) continue;
790     const char* detName = fEventsPerFile[i]->GetName();
791     const char* typeName = fEventsPerFile[i]->GetTitle();
792     TString loaderName(detName);
793     loaderName += "Loader";
794     AliLoader* loader = runLoader->GetLoader(loaderName);
795     if (!loader) {
796       AliError(Form("RunSimulation", "no loader for %s found\n"
797                     "Number of events per file not set for %s %s", 
798                     detName, typeName, detName));
799       continue;
800     }
801     AliDataLoader* dataLoader = 
802       loader->GetDataLoader(typeName);
803     if (!dataLoader) {
804       AliError(Form("no data loader for %s found\n"
805                     "Number of events per file not set for %s %s", 
806                     typeName, detName, typeName));
807       continue;
808     }
809     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
810     AliDebug(1, Form("number of events per file set to %d for %s %s",
811                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
812   }
813
814   AliInfo("running gAlice");
815   StdoutToAliInfo(StderrToAliError(
816     gAlice->Run(nEvents);
817   ););
818
819   delete runLoader;
820
821   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
822                stopwatch.RealTime(),stopwatch.CpuTime()));
823
824   return kTRUE;
825 }
826
827 //_____________________________________________________________________________
828 Bool_t AliSimulation::RunSDigitization(const char* detectors)
829 {
830 // run the digitization and produce summable digits
831
832   TStopwatch stopwatch;
833   stopwatch.Start();
834
835   AliRunLoader* runLoader = LoadRun();
836   if (!runLoader) return kFALSE;
837
838   TString detStr = detectors;
839   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
840   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
841     AliModule* det = (AliModule*) detArray->At(iDet);
842     if (!det || !det->IsActive()) continue;
843     if (IsSelected(det->GetName(), detStr)) {
844       AliInfo(Form("creating summable digits for %s", det->GetName()));
845       TStopwatch stopwatchDet;
846       stopwatchDet.Start();
847       det->Hits2SDigits();
848       AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
849            det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
850     }
851   }
852
853   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
854     AliError(Form("the following detectors were not found: %s",
855                   detStr.Data()));
856     if (fStopOnError) return kFALSE;
857   }
858
859   delete runLoader;
860
861   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
862            stopwatch.RealTime(),stopwatch.CpuTime()));
863
864   return kTRUE;
865 }
866
867
868 //_____________________________________________________________________________
869 Bool_t AliSimulation::RunDigitization(const char* detectors, 
870                                       const char* excludeDetectors)
871 {
872 // run the digitization and produce digits from sdigits
873
874   TStopwatch stopwatch;
875   stopwatch.Start();
876
877   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
878   if (gAlice) delete gAlice;
879   gAlice = NULL;
880
881   Int_t nStreams = 1;
882   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
883   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
884   AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
885   manager->SetInputStream(0, fGAliceFileName.Data());
886   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
887     const char* fileName = ((TObjString*)
888                             (fBkgrdFileNames->At(iStream-1)))->GetName();
889     manager->SetInputStream(iStream, fileName);
890   }
891
892   TString detStr = detectors;
893   TString detExcl = excludeDetectors;
894   manager->GetInputStream(0)->ImportgAlice();
895   AliRunLoader* runLoader = 
896     AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
897   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
898   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
899     AliModule* det = (AliModule*) detArray->At(iDet);
900     if (!det || !det->IsActive()) continue;
901     if (IsSelected(det->GetName(), detStr) && 
902         !IsSelected(det->GetName(), detExcl)) {
903       AliDigitizer* digitizer = det->CreateDigitizer(manager);
904       if (!digitizer) {
905         AliError(Form("no digitizer for %s", det->GetName()));
906         if (fStopOnError) return kFALSE;
907       } else {
908         digitizer->SetRegionOfInterest(fRegionOfInterest);
909       }
910     }
911   }
912
913   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
914     AliError(Form("the following detectors were not found: %s", 
915                   detStr.Data()));
916     if (fStopOnError) return kFALSE;
917   }
918
919   if (!manager->GetListOfTasks()->IsEmpty()) {
920     AliInfo("executing digitization");
921     manager->Exec("");
922   }
923
924   delete manager;
925
926   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
927                stopwatch.RealTime(),stopwatch.CpuTime()));
928   
929   return kTRUE;
930 }
931
932 //_____________________________________________________________________________
933 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
934 {
935 // run the digitization and produce digits from hits
936
937   TStopwatch stopwatch;
938   stopwatch.Start();
939
940   AliRunLoader* runLoader = LoadRun("READ");
941   if (!runLoader) return kFALSE;
942
943   TString detStr = detectors;
944   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
945   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
946     AliModule* det = (AliModule*) detArray->At(iDet);
947     if (!det || !det->IsActive()) continue;
948     if (IsSelected(det->GetName(), detStr)) {
949       AliInfo(Form("creating digits from hits for %s", det->GetName()));
950       det->Hits2Digits();
951     }
952   }
953
954   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
955     AliError(Form("the following detectors were not found: %s", 
956                   detStr.Data()));
957     if (fStopOnError) return kFALSE;
958   }
959
960   delete runLoader;
961   //PH Temporary fix to avoid interference with the PHOS loder/getter
962   //PH The problem has to be solved in more general way 09/06/05
963
964   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
965                stopwatch.RealTime(),stopwatch.CpuTime()));
966
967   return kTRUE;
968 }
969
970 //_____________________________________________________________________________
971 Bool_t AliSimulation::WriteRawData(const char* detectors, 
972                                    const char* fileName,
973                                    Bool_t deleteIntermediateFiles)
974 {
975 // convert the digits to raw data
976 // First DDL raw data files for the given detectors are created.
977 // If a file name is given, the DDL files are then converted to a DATE file.
978 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
979 // afterwards.
980 // If the file name has the extension ".root", the DATE file is converted
981 // to a root file.
982 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
983
984   TStopwatch stopwatch;
985   stopwatch.Start();
986
987   if (!WriteRawFiles(detectors)) {
988     if (fStopOnError) return kFALSE;
989   }
990
991   TString dateFileName(fileName);
992   if (!dateFileName.IsNull()) {
993     Bool_t rootOutput = dateFileName.EndsWith(".root");
994     if (rootOutput) dateFileName += ".date";
995     if (!ConvertRawFilesToDate(dateFileName)) {
996       if (fStopOnError) return kFALSE;
997     }
998     if (deleteIntermediateFiles) {
999       AliRunLoader* runLoader = LoadRun("READ");
1000       if (runLoader) for (Int_t iEvent = 0; 
1001                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1002         char command[256];
1003         sprintf(command, "rm -r raw%d", iEvent);
1004         gSystem->Exec(command);
1005       }
1006     }
1007
1008     if (rootOutput) {
1009       if (!ConvertDateToRoot(dateFileName, fileName)) {
1010         if (fStopOnError) return kFALSE;
1011       }
1012       if (deleteIntermediateFiles) {
1013         gSystem->Unlink(dateFileName);
1014       }
1015     }
1016   }
1017
1018   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1019                stopwatch.RealTime(),stopwatch.CpuTime()));
1020
1021   return kTRUE;
1022 }
1023
1024 //_____________________________________________________________________________
1025 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1026 {
1027 // convert the digits to raw data DDL files
1028
1029   AliRunLoader* runLoader = LoadRun("READ");
1030   if (!runLoader) return kFALSE;
1031
1032   // write raw data to DDL files
1033   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1034     AliInfo(Form("processing event %d", iEvent));
1035     runLoader->GetEvent(iEvent);
1036     TString baseDir = gSystem->WorkingDirectory();
1037     char dirName[256];
1038     sprintf(dirName, "raw%d", iEvent);
1039     gSystem->MakeDirectory(dirName);
1040     if (!gSystem->ChangeDirectory(dirName)) {
1041       AliError(Form("couldn't change to directory %s", dirName));
1042       if (fStopOnError) return kFALSE; else continue;
1043     }
1044
1045     TString detStr = detectors;
1046     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1047     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1048       AliModule* det = (AliModule*) detArray->At(iDet);
1049       if (!det || !det->IsActive()) continue;
1050       if (IsSelected(det->GetName(), detStr)) {
1051         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1052         det->Digits2Raw();
1053       }
1054     }
1055
1056     if (!WriteTriggerRawData())
1057       if (fStopOnError) return kFALSE;
1058
1059     gSystem->ChangeDirectory(baseDir);
1060     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1061       AliError(Form("the following detectors were not found: %s", 
1062                     detStr.Data()));
1063       if (fStopOnError) return kFALSE;
1064     }
1065   }
1066
1067   delete runLoader;
1068   return kTRUE;
1069 }
1070
1071 //_____________________________________________________________________________
1072 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1073 {
1074 // convert raw data DDL files to a DATE file with the program "dateStream"
1075
1076   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1077   if (!path) {
1078     AliError("the program dateStream was not found");
1079     if (fStopOnError) return kFALSE;
1080   } else {
1081     delete[] path;
1082   }
1083
1084   AliRunLoader* runLoader = LoadRun("READ");
1085   if (!runLoader) return kFALSE;
1086
1087   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1088   char command[256];
1089   sprintf(command, "dateStream -D -o %s -# %d -C", 
1090           dateFileName, runLoader->GetNumberOfEvents());
1091   FILE* pipe = gSystem->OpenPipe(command, "w");
1092
1093   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1094     fprintf(pipe, "GDC\n");
1095     Float_t ldc = 0;
1096     Int_t prevLDC = -1;
1097
1098     // loop over detectors and DDLs
1099     for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1100       for (Int_t iDDL = 0; iDDL < kDetectorDDLs[iDet]; iDDL++) {
1101
1102         Int_t ddlID = 0x100*iDet + iDDL;
1103         Int_t ldcID = Int_t(ldc + 0.0001);
1104         ldc += kDetectorLDCs[iDet] / kDetectorDDLs[iDet];
1105
1106         char rawFileName[256];
1107         sprintf(rawFileName, "raw%d/%s_%d.ddl", 
1108                 iEvent, kDetectors[iDet], ddlID);
1109
1110         // check existence and size of raw data file
1111         FILE* file = fopen(rawFileName, "rb");
1112         if (!file) continue;
1113         fseek(file, 0, SEEK_END);
1114         unsigned long size = ftell(file);
1115         fclose(file);
1116         if (!size) continue;
1117
1118         if (ldcID != prevLDC) {
1119           fprintf(pipe, " LDC Id %d\n", ldcID);
1120           prevLDC = ldcID;
1121         }
1122         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1123       }
1124     }
1125   }
1126
1127   Int_t result = gSystem->ClosePipe(pipe);
1128
1129   delete runLoader;
1130   return (result == 0);
1131 }
1132
1133 //_____________________________________________________________________________
1134 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1135                                         const char* rootFileName)
1136 {
1137 // convert a DATE file to a root file with the program "alimdc"
1138
1139   // ALIMDC setup
1140   const Int_t kDBSize = 1000000000;
1141   const Int_t kTagDBSize = 1000000000;
1142   const Bool_t kFilter = kFALSE;
1143   const Int_t kCompression = 1;
1144
1145   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1146   if (!path) {
1147     AliError("the program alimdc was not found");
1148     if (fStopOnError) return kFALSE;
1149   } else {
1150     delete[] path;
1151   }
1152
1153   AliInfo(Form("converting DATE file %s to root file %s", 
1154                dateFileName, rootFileName));
1155
1156   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1157   const char* tagDBFS    = "/tmp/mdc1/tags";
1158   const char* runDBFS    = "/tmp/mdc1/meta";
1159
1160   // User defined file system locations
1161   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1162     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1163   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1164     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1165   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1166     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1167   if (gSystem->Getenv("ALIMDC_RUNDB")) 
1168     runDBFS = gSystem->Getenv("ALIMDC_RUNDB");
1169
1170   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1171   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1172   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1173   gSystem->Exec(Form("rm -rf %s",runDBFS));
1174
1175   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1176   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1177   gSystem->Exec(Form("mkdir %s",tagDBFS));
1178   gSystem->Exec(Form("mkdir %s",runDBFS));
1179
1180   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1181                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1182   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1183
1184   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1185   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1186   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1187   gSystem->Exec(Form("rm -rf %s",runDBFS));
1188
1189   return (result == 0);
1190 }
1191
1192
1193 //_____________________________________________________________________________
1194 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1195 {
1196 // delete existing run loaders, open a new one and load gAlice
1197
1198   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1199   AliRunLoader* runLoader = 
1200     AliRunLoader::Open(fGAliceFileName.Data(), 
1201                        AliConfig::GetDefaultEventFolderName(), mode);
1202   if (!runLoader) {
1203     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1204     return NULL;
1205   }
1206   runLoader->LoadgAlice();
1207   gAlice = runLoader->GetAliRun();
1208   if (!gAlice) {
1209     AliError(Form("no gAlice object found in file %s", 
1210                   fGAliceFileName.Data()));
1211     return NULL;
1212   }
1213   return runLoader;
1214 }
1215
1216 //_____________________________________________________________________________
1217 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1218 {
1219 // get or calculate the number of signal events per background event
1220
1221   if (!fBkgrdFileNames) return 1;
1222   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1223   if (nBkgrdFiles == 0) return 1;
1224
1225   // get the number of signal events
1226   if (nEvents <= 0) {
1227     AliRunLoader* runLoader = 
1228       AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1229     if (!runLoader) return 1;
1230     nEvents = runLoader->GetNumberOfEvents();
1231     delete runLoader;
1232   }
1233
1234   Int_t result = 0;
1235   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1236     // get the number of background events
1237     const char* fileName = ((TObjString*)
1238                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1239     AliRunLoader* runLoader = 
1240       AliRunLoader::Open(fileName, "BKGRD");
1241     if (!runLoader) continue;
1242     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1243     delete runLoader;
1244
1245     // get or calculate the number of signal per background events
1246     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1247     if (nSignalPerBkgrd <= 0) {
1248       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1249     } else if (result && (result != nSignalPerBkgrd)) {
1250       AliInfo(Form("the number of signal events per background event "
1251                    "will be changed from %d to %d for stream %d", 
1252                    nSignalPerBkgrd, result, iBkgrdFile+1));
1253       nSignalPerBkgrd = result;
1254     }
1255
1256     if (!result) result = nSignalPerBkgrd;
1257     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1258       AliWarning(Form("not enough background events (%d) for %d signal events "
1259                       "using %d signal per background events for stream %d",
1260                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1261     }
1262   }
1263
1264   return result;
1265 }
1266
1267 //_____________________________________________________________________________
1268 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1269 {
1270 // check whether detName is contained in detectors
1271 // if yes, it is removed from detectors
1272
1273   // check if all detectors are selected
1274   if ((detectors.CompareTo("ALL") == 0) ||
1275       detectors.BeginsWith("ALL ") ||
1276       detectors.EndsWith(" ALL") ||
1277       detectors.Contains(" ALL ")) {
1278     detectors = "ALL";
1279     return kTRUE;
1280   }
1281
1282   // search for the given detector
1283   Bool_t result = kFALSE;
1284   if ((detectors.CompareTo(detName) == 0) ||
1285       detectors.BeginsWith(detName+" ") ||
1286       detectors.EndsWith(" "+detName) ||
1287       detectors.Contains(" "+detName+" ")) {
1288     detectors.ReplaceAll(detName, "");
1289     result = kTRUE;
1290   }
1291
1292   // clean up the detectors string
1293   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1294   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1295   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1296
1297   return result;
1298 }