]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
Geometry builder classes moved from base to sim.
[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->GetTopNode()->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)) {
736 //     delete runLoader;
737 //     return kFALSE;
738 //   }
739   MisalignGeometry(runLoader);
740
741   // Temporary fix by A.Gheata
742   // Could be removed with the next Root version (>5.11)
743   if (gGeoManager) {
744     TIter next(gGeoManager->GetListOfVolumes());
745     TGeoVolume *vol;
746     while ((vol = (TGeoVolume *)next())) {
747       if (vol->GetVoxels()) {
748         if (vol->GetVoxels()->NeedRebuild()) {
749           vol->GetVoxels()->Voxelize();
750           vol->FindOverlaps();
751         }
752       }
753     }
754   }
755
756   // Export (mis)aligned geometry 
757   if (gGeoManager) gGeoManager->Export("misaligned_geometry.root");
758
759 //   AliRunLoader* runLoader = gAlice->GetRunLoader();
760 //   if (!runLoader) {
761 //     AliError(Form("gAlice has no run loader object. "
762 //                   "Check your config file: %s", fConfigFileName.Data()));
763 //     return kFALSE;
764 //   }
765 //   SetGAliceFile(runLoader->GetFileName());
766
767   if (!gAlice->Generator()) {
768     AliError(Form("gAlice has no generator object. "
769                   "Check your config file: %s", fConfigFileName.Data()));
770     return kFALSE;
771   }
772   if (nEvents <= 0) nEvents = fNEvents;
773
774   // get vertex from background file in case of merging
775   if (fUseBkgrdVertex &&
776       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
777     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
778     const char* fileName = ((TObjString*)
779                             (fBkgrdFileNames->At(0)))->GetName();
780     AliInfo(Form("The vertex will be taken from the background "
781                  "file %s with nSignalPerBackground = %d", 
782                  fileName, signalPerBkgrd));
783     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
784     gAlice->Generator()->SetVertexGenerator(vtxGen);
785   }
786
787   if (!fRunSimulation) {
788     gAlice->Generator()->SetTrackingFlag(0);
789   }
790
791   // set the number of events per file for given detectors and data types
792   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
793     if (!fEventsPerFile[i]) continue;
794     const char* detName = fEventsPerFile[i]->GetName();
795     const char* typeName = fEventsPerFile[i]->GetTitle();
796     TString loaderName(detName);
797     loaderName += "Loader";
798     AliLoader* loader = runLoader->GetLoader(loaderName);
799     if (!loader) {
800       AliError(Form("RunSimulation", "no loader for %s found\n"
801                     "Number of events per file not set for %s %s", 
802                     detName, typeName, detName));
803       continue;
804     }
805     AliDataLoader* dataLoader = 
806       loader->GetDataLoader(typeName);
807     if (!dataLoader) {
808       AliError(Form("no data loader for %s found\n"
809                     "Number of events per file not set for %s %s", 
810                     typeName, detName, typeName));
811       continue;
812     }
813     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
814     AliDebug(1, Form("number of events per file set to %d for %s %s",
815                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
816   }
817
818   AliInfo("running gAlice");
819   StdoutToAliInfo(StderrToAliError(
820     gAlice->Run(nEvents);
821   ););
822
823   delete runLoader;
824
825   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
826                stopwatch.RealTime(),stopwatch.CpuTime()));
827
828   return kTRUE;
829 }
830
831 //_____________________________________________________________________________
832 Bool_t AliSimulation::RunSDigitization(const char* detectors)
833 {
834 // run the digitization and produce summable digits
835
836   TStopwatch stopwatch;
837   stopwatch.Start();
838
839   AliRunLoader* runLoader = LoadRun();
840   if (!runLoader) return kFALSE;
841
842   TString detStr = detectors;
843   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
844   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
845     AliModule* det = (AliModule*) detArray->At(iDet);
846     if (!det || !det->IsActive()) continue;
847     if (IsSelected(det->GetName(), detStr)) {
848       AliInfo(Form("creating summable digits for %s", det->GetName()));
849       TStopwatch stopwatchDet;
850       stopwatchDet.Start();
851       det->Hits2SDigits();
852       AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
853            det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
854     }
855   }
856
857   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
858     AliError(Form("the following detectors were not found: %s",
859                   detStr.Data()));
860     if (fStopOnError) return kFALSE;
861   }
862
863   delete runLoader;
864
865   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
866            stopwatch.RealTime(),stopwatch.CpuTime()));
867
868   return kTRUE;
869 }
870
871
872 //_____________________________________________________________________________
873 Bool_t AliSimulation::RunDigitization(const char* detectors, 
874                                       const char* excludeDetectors)
875 {
876 // run the digitization and produce digits from sdigits
877
878   TStopwatch stopwatch;
879   stopwatch.Start();
880
881   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
882   if (gAlice) delete gAlice;
883   gAlice = NULL;
884
885   Int_t nStreams = 1;
886   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
887   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
888   AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
889   manager->SetInputStream(0, fGAliceFileName.Data());
890   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
891     const char* fileName = ((TObjString*)
892                             (fBkgrdFileNames->At(iStream-1)))->GetName();
893     manager->SetInputStream(iStream, fileName);
894   }
895
896   TString detStr = detectors;
897   TString detExcl = excludeDetectors;
898   manager->GetInputStream(0)->ImportgAlice();
899   AliRunLoader* runLoader = 
900     AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
901   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
902   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
903     AliModule* det = (AliModule*) detArray->At(iDet);
904     if (!det || !det->IsActive()) continue;
905     if (IsSelected(det->GetName(), detStr) && 
906         !IsSelected(det->GetName(), detExcl)) {
907       AliDigitizer* digitizer = det->CreateDigitizer(manager);
908       if (!digitizer) {
909         AliError(Form("no digitizer for %s", det->GetName()));
910         if (fStopOnError) return kFALSE;
911       } else {
912         digitizer->SetRegionOfInterest(fRegionOfInterest);
913       }
914     }
915   }
916
917   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
918     AliError(Form("the following detectors were not found: %s", 
919                   detStr.Data()));
920     if (fStopOnError) return kFALSE;
921   }
922
923   if (!manager->GetListOfTasks()->IsEmpty()) {
924     AliInfo("executing digitization");
925     manager->Exec("");
926   }
927
928   delete manager;
929
930   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
931                stopwatch.RealTime(),stopwatch.CpuTime()));
932   
933   return kTRUE;
934 }
935
936 //_____________________________________________________________________________
937 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
938 {
939 // run the digitization and produce digits from hits
940
941   TStopwatch stopwatch;
942   stopwatch.Start();
943
944   AliRunLoader* runLoader = LoadRun("READ");
945   if (!runLoader) return kFALSE;
946
947   TString detStr = detectors;
948   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
949   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
950     AliModule* det = (AliModule*) detArray->At(iDet);
951     if (!det || !det->IsActive()) continue;
952     if (IsSelected(det->GetName(), detStr)) {
953       AliInfo(Form("creating digits from hits for %s", det->GetName()));
954       det->Hits2Digits();
955     }
956   }
957
958   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
959     AliError(Form("the following detectors were not found: %s", 
960                   detStr.Data()));
961     if (fStopOnError) return kFALSE;
962   }
963
964   delete runLoader;
965   //PH Temporary fix to avoid interference with the PHOS loder/getter
966   //PH The problem has to be solved in more general way 09/06/05
967
968   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
969                stopwatch.RealTime(),stopwatch.CpuTime()));
970
971   return kTRUE;
972 }
973
974 //_____________________________________________________________________________
975 Bool_t AliSimulation::WriteRawData(const char* detectors, 
976                                    const char* fileName,
977                                    Bool_t deleteIntermediateFiles)
978 {
979 // convert the digits to raw data
980 // First DDL raw data files for the given detectors are created.
981 // If a file name is given, the DDL files are then converted to a DATE file.
982 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
983 // afterwards.
984 // If the file name has the extension ".root", the DATE file is converted
985 // to a root file.
986 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
987
988   TStopwatch stopwatch;
989   stopwatch.Start();
990
991   if (!WriteRawFiles(detectors)) {
992     if (fStopOnError) return kFALSE;
993   }
994
995   TString dateFileName(fileName);
996   if (!dateFileName.IsNull()) {
997     Bool_t rootOutput = dateFileName.EndsWith(".root");
998     if (rootOutput) dateFileName += ".date";
999     if (!ConvertRawFilesToDate(dateFileName)) {
1000       if (fStopOnError) return kFALSE;
1001     }
1002     if (deleteIntermediateFiles) {
1003       AliRunLoader* runLoader = LoadRun("READ");
1004       if (runLoader) for (Int_t iEvent = 0; 
1005                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1006         char command[256];
1007         sprintf(command, "rm -r raw%d", iEvent);
1008         gSystem->Exec(command);
1009       }
1010     }
1011
1012     if (rootOutput) {
1013       if (!ConvertDateToRoot(dateFileName, fileName)) {
1014         if (fStopOnError) return kFALSE;
1015       }
1016       if (deleteIntermediateFiles) {
1017         gSystem->Unlink(dateFileName);
1018       }
1019     }
1020   }
1021
1022   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1023                stopwatch.RealTime(),stopwatch.CpuTime()));
1024
1025   return kTRUE;
1026 }
1027
1028 //_____________________________________________________________________________
1029 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1030 {
1031 // convert the digits to raw data DDL files
1032
1033   AliRunLoader* runLoader = LoadRun("READ");
1034   if (!runLoader) return kFALSE;
1035
1036   // write raw data to DDL files
1037   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1038     AliInfo(Form("processing event %d", iEvent));
1039     runLoader->GetEvent(iEvent);
1040     TString baseDir = gSystem->WorkingDirectory();
1041     char dirName[256];
1042     sprintf(dirName, "raw%d", iEvent);
1043     gSystem->MakeDirectory(dirName);
1044     if (!gSystem->ChangeDirectory(dirName)) {
1045       AliError(Form("couldn't change to directory %s", dirName));
1046       if (fStopOnError) return kFALSE; else continue;
1047     }
1048
1049     TString detStr = detectors;
1050     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1051     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1052       AliModule* det = (AliModule*) detArray->At(iDet);
1053       if (!det || !det->IsActive()) continue;
1054       if (IsSelected(det->GetName(), detStr)) {
1055         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1056         det->Digits2Raw();
1057       }
1058     }
1059
1060     if (!WriteTriggerRawData())
1061       if (fStopOnError) return kFALSE;
1062
1063     gSystem->ChangeDirectory(baseDir);
1064     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1065       AliError(Form("the following detectors were not found: %s", 
1066                     detStr.Data()));
1067       if (fStopOnError) return kFALSE;
1068     }
1069   }
1070
1071   delete runLoader;
1072   return kTRUE;
1073 }
1074
1075 //_____________________________________________________________________________
1076 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1077 {
1078 // convert raw data DDL files to a DATE file with the program "dateStream"
1079
1080   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1081   if (!path) {
1082     AliError("the program dateStream was not found");
1083     if (fStopOnError) return kFALSE;
1084   } else {
1085     delete[] path;
1086   }
1087
1088   AliRunLoader* runLoader = LoadRun("READ");
1089   if (!runLoader) return kFALSE;
1090
1091   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1092   char command[256];
1093   sprintf(command, "dateStream -D -o %s -# %d -C", 
1094           dateFileName, runLoader->GetNumberOfEvents());
1095   FILE* pipe = gSystem->OpenPipe(command, "w");
1096
1097   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1098     fprintf(pipe, "GDC\n");
1099     Float_t ldc = 0;
1100     Int_t prevLDC = -1;
1101
1102     // loop over detectors and DDLs
1103     for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1104       for (Int_t iDDL = 0; iDDL < kDetectorDDLs[iDet]; iDDL++) {
1105
1106         Int_t ddlID = 0x100*iDet + iDDL;
1107         Int_t ldcID = Int_t(ldc + 0.0001);
1108         ldc += kDetectorLDCs[iDet] / kDetectorDDLs[iDet];
1109
1110         char rawFileName[256];
1111         sprintf(rawFileName, "raw%d/%s_%d.ddl", 
1112                 iEvent, kDetectors[iDet], ddlID);
1113
1114         // check existence and size of raw data file
1115         FILE* file = fopen(rawFileName, "rb");
1116         if (!file) continue;
1117         fseek(file, 0, SEEK_END);
1118         unsigned long size = ftell(file);
1119         fclose(file);
1120         if (!size) continue;
1121
1122         if (ldcID != prevLDC) {
1123           fprintf(pipe, " LDC Id %d\n", ldcID);
1124           prevLDC = ldcID;
1125         }
1126         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1127       }
1128     }
1129   }
1130
1131   Int_t result = gSystem->ClosePipe(pipe);
1132
1133   delete runLoader;
1134   return (result == 0);
1135 }
1136
1137 //_____________________________________________________________________________
1138 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1139                                         const char* rootFileName)
1140 {
1141 // convert a DATE file to a root file with the program "alimdc"
1142
1143   // ALIMDC setup
1144   const Int_t kDBSize = 1000000000;
1145   const Int_t kTagDBSize = 1000000000;
1146   const Bool_t kFilter = kFALSE;
1147   const Int_t kCompression = 1;
1148
1149   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1150   if (!path) {
1151     AliError("the program alimdc was not found");
1152     if (fStopOnError) return kFALSE;
1153   } else {
1154     delete[] path;
1155   }
1156
1157   AliInfo(Form("converting DATE file %s to root file %s", 
1158                dateFileName, rootFileName));
1159
1160   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1161   const char* tagDBFS    = "/tmp/mdc1/tags";
1162   const char* runDBFS    = "/tmp/mdc1/meta";
1163
1164   // User defined file system locations
1165   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1166     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1167   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1168     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1169   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1170     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1171   if (gSystem->Getenv("ALIMDC_RUNDB")) 
1172     runDBFS = gSystem->Getenv("ALIMDC_RUNDB");
1173
1174   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1175   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1176   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1177   gSystem->Exec(Form("rm -rf %s",runDBFS));
1178
1179   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1180   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1181   gSystem->Exec(Form("mkdir %s",tagDBFS));
1182   gSystem->Exec(Form("mkdir %s",runDBFS));
1183
1184   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1185                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1186   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1187
1188   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1189   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1190   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1191   gSystem->Exec(Form("rm -rf %s",runDBFS));
1192
1193   return (result == 0);
1194 }
1195
1196
1197 //_____________________________________________________________________________
1198 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1199 {
1200 // delete existing run loaders, open a new one and load gAlice
1201
1202   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1203   AliRunLoader* runLoader = 
1204     AliRunLoader::Open(fGAliceFileName.Data(), 
1205                        AliConfig::GetDefaultEventFolderName(), mode);
1206   if (!runLoader) {
1207     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1208     return NULL;
1209   }
1210   runLoader->LoadgAlice();
1211   gAlice = runLoader->GetAliRun();
1212   if (!gAlice) {
1213     AliError(Form("no gAlice object found in file %s", 
1214                   fGAliceFileName.Data()));
1215     return NULL;
1216   }
1217   return runLoader;
1218 }
1219
1220 //_____________________________________________________________________________
1221 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1222 {
1223 // get or calculate the number of signal events per background event
1224
1225   if (!fBkgrdFileNames) return 1;
1226   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1227   if (nBkgrdFiles == 0) return 1;
1228
1229   // get the number of signal events
1230   if (nEvents <= 0) {
1231     AliRunLoader* runLoader = 
1232       AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1233     if (!runLoader) return 1;
1234     nEvents = runLoader->GetNumberOfEvents();
1235     delete runLoader;
1236   }
1237
1238   Int_t result = 0;
1239   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1240     // get the number of background events
1241     const char* fileName = ((TObjString*)
1242                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1243     AliRunLoader* runLoader = 
1244       AliRunLoader::Open(fileName, "BKGRD");
1245     if (!runLoader) continue;
1246     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1247     delete runLoader;
1248
1249     // get or calculate the number of signal per background events
1250     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1251     if (nSignalPerBkgrd <= 0) {
1252       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1253     } else if (result && (result != nSignalPerBkgrd)) {
1254       AliInfo(Form("the number of signal events per background event "
1255                    "will be changed from %d to %d for stream %d", 
1256                    nSignalPerBkgrd, result, iBkgrdFile+1));
1257       nSignalPerBkgrd = result;
1258     }
1259
1260     if (!result) result = nSignalPerBkgrd;
1261     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1262       AliWarning(Form("not enough background events (%d) for %d signal events "
1263                       "using %d signal per background events for stream %d",
1264                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1265     }
1266   }
1267
1268   return result;
1269 }
1270
1271 //_____________________________________________________________________________
1272 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1273 {
1274 // check whether detName is contained in detectors
1275 // if yes, it is removed from detectors
1276
1277   // check if all detectors are selected
1278   if ((detectors.CompareTo("ALL") == 0) ||
1279       detectors.BeginsWith("ALL ") ||
1280       detectors.EndsWith(" ALL") ||
1281       detectors.Contains(" ALL ")) {
1282     detectors = "ALL";
1283     return kTRUE;
1284   }
1285
1286   // search for the given detector
1287   Bool_t result = kFALSE;
1288   if ((detectors.CompareTo(detName) == 0) ||
1289       detectors.BeginsWith(detName+" ") ||
1290       detectors.EndsWith(" "+detName) ||
1291       detectors.Contains(" "+detName+" ")) {
1292     detectors.ReplaceAll(detName, "");
1293     result = kTRUE;
1294   }
1295
1296   // clean up the detectors string
1297   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1298   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1299   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1300
1301   return result;
1302 }