]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
reorganisation of the creation of the QADataMaker solving mem leaks
[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 <TVirtualMCApplication.h>
109 #include <TGeoManager.h>
110 #include <TObjString.h>
111 #include <TSystem.h>
112 #include <TFile.h>
113
114 #include "AliCodeTimer.h"
115 #include "AliCDBStorage.h"
116 #include "AliCDBEntry.h"
117 #include "AliCDBManager.h"
118 #include "AliGeomManager.h"
119 #include "AliAlignObj.h"
120 #include "AliCentralTrigger.h"
121 #include "AliDAQ.h"
122 #include "AliDigitizer.h"
123 #include "AliGenerator.h"
124 #include "AliLog.h"
125 #include "AliModule.h"
126 #include "AliRun.h"
127 #include "AliRunDigitizer.h"
128 #include "AliRunLoader.h"
129 #include "AliSimulation.h"
130 #include "AliVertexGenFile.h"
131 #include "AliCentralTrigger.h"
132 #include "AliCTPRawData.h"
133 #include "AliRawReaderFile.h"
134 #include "AliRawReaderRoot.h"
135 #include "AliRawReaderDate.h"
136 #include "AliESD.h"
137 #include "AliHeader.h"
138 #include "AliGenEventHeader.h"
139 #include "AliMC.h"
140 #include "AliHLTSimulation.h"
141 #include "AliQADataMakerSteer.h"
142 #include "AliSysInfo.h"
143
144 ClassImp(AliSimulation)
145
146 AliSimulation *AliSimulation::fgInstance = 0;
147 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
148
149 //_____________________________________________________________________________
150 AliSimulation::AliSimulation(const char* configFileName,
151                              const char* name, const char* title) :
152   TNamed(name, title),
153
154   fRunGeneration(kTRUE),
155   fRunSimulation(kTRUE),
156   fLoadAlignFromCDB(kTRUE),
157   fLoadAlObjsListOfDets("ALL"),
158   fMakeSDigits("ALL"),
159   fMakeDigits("ALL"),
160   fMakeTrigger(""),
161   fMakeDigitsFromHits(""),
162   fWriteRawData(""),
163   fRawDataFileName(""),
164   fDeleteIntermediateFiles(kFALSE),
165   fWriteSelRawData(kFALSE),
166   fStopOnError(kFALSE),
167
168   fNEvents(1),
169   fConfigFileName(configFileName),
170   fGAliceFileName("galice.root"),
171   fEventsPerFile(),
172   fBkgrdFileNames(NULL),
173   fAlignObjArray(NULL),
174   fUseBkgrdVertex(kTRUE),
175   fRegionOfInterest(kFALSE),
176   fCDBUri(""),
177   fSpecCDBUri(),
178   fRun(-1),
179   fSeed(0),
180   fInitCDBCalled(kFALSE),
181   fInitRunNumberCalled(kFALSE),
182   fSetRunNumberFromDataCalled(kFALSE),
183   fEmbeddingFlag(kFALSE),
184   fQADetectors("ALL"),                  
185   fQATasks("ALL"),      
186   fRunQA(kTRUE), 
187   fRunHLT("default")
188 {
189 // create simulation object with default parameters
190   fgInstance = this;
191   SetGAliceFile("galice.root");
192   
193 // for QA
194    for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) 
195         fQACycles[iDet] = 999999;
196 }
197
198 //_____________________________________________________________________________
199 AliSimulation::AliSimulation(const AliSimulation& sim) :
200   TNamed(sim),
201
202   fRunGeneration(sim.fRunGeneration),
203   fRunSimulation(sim.fRunSimulation),
204   fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
205   fLoadAlObjsListOfDets(sim.fLoadAlObjsListOfDets),
206   fMakeSDigits(sim.fMakeSDigits),
207   fMakeDigits(sim.fMakeDigits),
208   fMakeTrigger(sim.fMakeTrigger),
209   fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
210   fWriteRawData(sim.fWriteRawData),
211   fRawDataFileName(""),
212   fDeleteIntermediateFiles(kFALSE),
213   fWriteSelRawData(kFALSE),
214   fStopOnError(sim.fStopOnError),
215
216   fNEvents(sim.fNEvents),
217   fConfigFileName(sim.fConfigFileName),
218   fGAliceFileName(sim.fGAliceFileName),
219   fEventsPerFile(),
220   fBkgrdFileNames(NULL),
221   fAlignObjArray(NULL),
222   fUseBkgrdVertex(sim.fUseBkgrdVertex),
223   fRegionOfInterest(sim.fRegionOfInterest),
224   fCDBUri(sim.fCDBUri),
225   fSpecCDBUri(),
226   fRun(-1),
227   fSeed(0),
228   fInitCDBCalled(sim.fInitCDBCalled),
229   fInitRunNumberCalled(sim.fInitRunNumberCalled),
230   fSetRunNumberFromDataCalled(sim.fSetRunNumberFromDataCalled),
231   fEmbeddingFlag(sim.fEmbeddingFlag),
232   fQADetectors(sim.fQADetectors),                  
233   fQATasks(sim.fQATasks),       
234   fRunQA(sim.fRunQA), 
235   fRunHLT(sim.fRunHLT)
236 {
237 // copy constructor
238
239   for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
240     if (!sim.fEventsPerFile[i]) continue;
241     fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
242   }
243
244   fBkgrdFileNames = new TObjArray;
245   for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
246     if (!sim.fBkgrdFileNames->At(i)) continue;
247     fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
248   }
249
250   for (Int_t i = 0; i < sim.fSpecCDBUri.GetEntriesFast(); i++) {
251     if (sim.fSpecCDBUri[i]) fSpecCDBUri.Add(sim.fSpecCDBUri[i]->Clone());
252   }
253   fgInstance = this;
254
255 // for QA
256    for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) 
257         fQACycles[iDet] = sim.fQACycles[iDet];
258 }
259
260 //_____________________________________________________________________________
261 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
262 {
263 // assignment operator
264
265   this->~AliSimulation();
266   new(this) AliSimulation(sim);
267   return *this;
268 }
269
270 //_____________________________________________________________________________
271 AliSimulation::~AliSimulation()
272 {
273 // clean up
274
275   fEventsPerFile.Delete();
276 //  if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
277 //  delete fAlignObjArray; fAlignObjArray=0;
278
279   if (fBkgrdFileNames) {
280     fBkgrdFileNames->Delete();
281     delete fBkgrdFileNames;
282   }
283
284   fSpecCDBUri.Delete();
285   if (fgInstance==this) fgInstance = 0;
286
287   AliCodeTimer::Instance()->Print();
288 }
289
290
291 //_____________________________________________________________________________
292 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
293 {
294 // set the number of events for one run
295
296   fNEvents = nEvents;
297 }
298
299 //_____________________________________________________________________________
300 void AliSimulation::InitCDB()
301 {
302 // activate a default CDB storage
303 // First check if we have any CDB storage set, because it is used 
304 // to retrieve the calibration and alignment constants
305
306   if (fInitCDBCalled) return;
307   fInitCDBCalled = kTRUE;
308
309   AliCDBManager* man = AliCDBManager::Instance();
310   if (man->IsDefaultStorageSet())
311   {
312     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
313     AliWarning("Default CDB storage has been already set !");
314     AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
315     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
316     fCDBUri = man->GetDefaultStorage()->GetURI();
317   }
318   else {
319     if (fCDBUri.Length() > 0) 
320     {
321         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
322         AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
323         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
324     } else {
325         fCDBUri="local://$ALICE_ROOT";
326         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
327         AliWarning("Default CDB storage not yet set !!!!");
328         AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
329         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
330                 
331     }
332     man->SetDefaultStorage(fCDBUri);
333   }
334
335   // Now activate the detector specific CDB storage locations
336   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
337     TObject* obj = fSpecCDBUri[i];
338     if (!obj) continue;
339     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
340     AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
341     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
342     man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
343   }
344       
345 }
346
347 //_____________________________________________________________________________
348 void AliSimulation::InitRunNumber(){
349 // check run number. If not set, set it to 0 !!!!
350   
351   if (fInitRunNumberCalled) return;
352   fInitRunNumberCalled = kTRUE;
353   
354   AliCDBManager* man = AliCDBManager::Instance();
355   if (man->GetRun() >= 0)
356   {
357         AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
358                         "Use external variable DC_RUN or AliSimulation::SetRun()!"));
359   }
360     
361   if(fRun >= 0) {
362         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
363         AliDebug(2, Form("Setting CDB run number to: %d",fRun));
364         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
365   } else {
366         fRun=0;
367         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
368         AliWarning("Run number not yet set !!!!");
369         AliWarning(Form("Setting it now to: %d", fRun));
370         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
371         
372   }
373   man->SetRun(fRun);
374
375   man->Print();
376
377 }
378
379 //_____________________________________________________________________________
380 void AliSimulation::SetCDBLock() {
381   // Set CDB lock: from now on it is forbidden to reset the run number
382   // or the default storage or to activate any further storage!
383   
384   AliCDBManager::Instance()->SetLock(1);
385 }
386
387 //_____________________________________________________________________________
388 void AliSimulation::SetDefaultStorage(const char* uri) {
389 // Store the desired default CDB storage location
390 // Activate it later within the Run() method
391
392   fCDBUri = uri;
393
394 }
395
396 //_____________________________________________________________________________
397 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
398 // Store a detector-specific CDB storage location
399 // Activate it later within the Run() method
400
401   AliCDBPath aPath(calibType);
402   if(!aPath.IsValid()){
403         AliError(Form("Not a valid path: %s", calibType));
404         return;
405   }
406
407   TObject* obj = fSpecCDBUri.FindObject(calibType);
408   if (obj) fSpecCDBUri.Remove(obj);
409   fSpecCDBUri.Add(new TNamed(calibType, uri));
410
411 }
412
413 //_____________________________________________________________________________
414 void AliSimulation::SetRunNumber(Int_t run)
415 {
416 // sets run number
417 // Activate it later within the Run() method
418
419         fRun = run;
420 }
421
422 //_____________________________________________________________________________
423 void AliSimulation::SetSeed(Int_t seed)
424 {
425 // sets seed number
426 // Activate it later within the Run() method
427
428         fSeed = seed;
429 }
430
431 //_____________________________________________________________________________
432 Bool_t AliSimulation::SetRunNumberFromData()
433 {
434   // Set the CDB manager run number
435   // The run number is retrieved from gAlice
436
437     if (fSetRunNumberFromDataCalled) return kTRUE;
438     fSetRunNumberFromDataCalled = kTRUE;    
439   
440     AliCDBManager* man = AliCDBManager::Instance();
441     Int_t runData = -1, runCDB = -1;
442   
443     AliRunLoader* runLoader = LoadRun("READ");
444     if (!runLoader) return kFALSE;
445     else {
446         runData = runLoader->GetAliRun()->GetHeader()->GetRun();
447         delete runLoader;
448     }
449   
450     runCDB = man->GetRun();
451     if(runCDB >= 0) {
452         if (runCDB != runData) {
453                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
454                 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
455                 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
456                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");    
457         }
458         
459     }
460       
461     man->SetRun(runData);
462     fRun = runData;
463     
464     if(man->GetRun() < 0) {
465         AliError("Run number not properly initalized!");
466         return kFALSE;
467     }
468   
469     man->Print();
470     
471     return kTRUE;
472 }
473
474 //_____________________________________________________________________________
475 void AliSimulation::SetConfigFile(const char* fileName)
476 {
477 // set the name of the config file
478
479   fConfigFileName = fileName;
480 }
481
482 //_____________________________________________________________________________
483 void AliSimulation::SetGAliceFile(const char* fileName)
484 {
485 // set the name of the galice file
486 // the path is converted to an absolute one if it is relative
487
488   fGAliceFileName = fileName;
489   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
490     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
491                                                 fGAliceFileName);
492     fGAliceFileName = absFileName;
493     delete[] absFileName;
494   }
495
496   AliDebug(2, Form("galice file name set to %s", fileName));
497 }
498
499 //_____________________________________________________________________________
500 void AliSimulation::SetEventsPerFile(const char* detector, const char* type, 
501                                      Int_t nEvents)
502 {
503 // set the number of events per file for the given detector and data type
504 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
505
506   TNamed* obj = new TNamed(detector, type);
507   obj->SetUniqueID(nEvents);
508   fEventsPerFile.Add(obj);
509 }
510
511 //_____________________________________________________________________________
512 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
513 {
514   // Read the alignment objects from CDB.
515   // Each detector is supposed to have the
516   // alignment objects in DET/Align/Data CDB path.
517   // All the detector objects are then collected,
518   // sorted by geometry level (starting from ALIC) and
519   // then applied to the TGeo geometry.
520   // Finally an overlaps check is performed.
521
522   if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
523     AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
524     return kFALSE;
525   }  
526   
527   // initialize CDB storage, run number, set CDB lock
528   InitCDB();
529 //  if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
530   SetCDBLock();
531     
532   Bool_t delRunLoader = kFALSE;
533   if (!runLoader) {
534     runLoader = LoadRun("READ");
535     if (!runLoader) return kFALSE;
536     delRunLoader = kTRUE;
537   }
538   
539   // Export ideal geometry 
540   if(!gAlice->IsRootGeometry()) AliGeomManager::GetGeometry()->Export("geometry.root");
541
542   // Load alignment data from CDB and apply to geometry through AliGeomManager
543   if(fLoadAlignFromCDB){
544     
545     TString detStr = fLoadAlObjsListOfDets;
546     TString loadAlObjsListOfDets = "";
547     
548     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
549     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
550       AliModule* det = (AliModule*) detArray->At(iDet);
551       if (!det || !det->IsActive()) continue;
552       if (IsSelected(det->GetName(), detStr)) {
553         //add det to list of dets to be aligned from CDB
554         loadAlObjsListOfDets += det->GetName();
555         loadAlObjsListOfDets += " ";
556       }
557     } // end loop over detectors
558     loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
559     AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
560   }else{
561     // Check if the array with alignment objects was
562     // provided by the user. If yes, apply the objects
563     // to the present TGeo geometry
564     if (fAlignObjArray) {
565       if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
566         AliError("The misalignment of one or more volumes failed!"
567                  "Compare the list of simulated detectors and the list of detector alignment data!");
568         if (delRunLoader) delete runLoader;
569         return kFALSE;
570       }
571     }
572   }
573
574   // Update the internal geometry of modules (ITS needs it)
575   TString detStr = fLoadAlObjsListOfDets;
576   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
577   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
578
579     AliModule* det = (AliModule*) detArray->At(iDet);
580     if (!det || !det->IsActive()) continue;
581     if (IsSelected(det->GetName(), detStr)) {
582       det->UpdateInternalGeometry();
583     }
584   } // end loop over detectors
585
586
587   if (delRunLoader) delete runLoader;
588
589   return kTRUE;
590 }
591
592 //_____________________________________________________________________________
593 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
594 {
595 // add a file with background events for merging
596
597   TObjString* fileNameStr = new TObjString(fileName);
598   fileNameStr->SetUniqueID(nSignalPerBkgrd);
599   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
600   fBkgrdFileNames->Add(fileNameStr);
601 }
602
603 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
604 {
605 // add a file with background events for embeddin
606   MergeWith(fileName, nSignalPerBkgrd);
607   fEmbeddingFlag = kTRUE;
608 }
609
610 //_____________________________________________________________________________
611 Bool_t AliSimulation::Run(Int_t nEvents)
612 {
613 // run the generation, simulation and digitization
614
615  
616   AliCodeTimerAuto("")
617   
618   // Load run number and seed from environmental vars
619   ProcessEnvironmentVars();
620
621   gRandom->SetSeed(fSeed);
622    
623   if (nEvents > 0) fNEvents = nEvents;
624
625   // generation and simulation -> hits
626   if (fRunGeneration) {
627     if (!RunSimulation()) if (fStopOnError) return kFALSE;
628   }
629            
630   // initialize CDB storage from external environment
631   // (either CDB manager or AliSimulation setters),
632   // if not already done in RunSimulation()
633   InitCDB();
634   
635   // Set run number in CDBManager from data 
636   // From this point on the run number must be always loaded from data!
637   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
638   
639   // Set CDB lock: from now on it is forbidden to reset the run number
640   // or the default storage or to activate any further storage!
641   SetCDBLock();
642
643   // If RunSimulation was not called, load the geometry and misalign it
644   if (!AliGeomManager::GetGeometry()) {
645     // Initialize the geometry manager
646     AliGeomManager::LoadGeometry("geometry.root");
647     
648 //    // Check that the consistency of symbolic names for the activated subdetectors
649 //    // in the geometry loaded by AliGeomManager
650 //    AliRunLoader* runLoader = LoadRun("READ");
651 //    if (!runLoader) return kFALSE;
652 //
653 //    TString detsToBeChecked = "";
654 //    TObjArray* detArray = runLoader->GetAliRun()->Detectors();
655 //    for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
656 //      AliModule* det = (AliModule*) detArray->At(iDet);
657 //      if (!det || !det->IsActive()) continue;
658 //      detsToBeChecked += det->GetName();
659 //      detsToBeChecked += " ";
660 //    } // end loop over detectors
661 //    if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
662     if(!AliGeomManager::CheckSymNamesLUT("ALL"))
663         AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
664         
665     if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
666     // Misalign geometry
667     if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
668   }
669
670
671   // hits -> summable digits
672   AliSysInfo::AddStamp("Start_sdigitization");
673   if (!fMakeSDigits.IsNull()) {
674     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
675  
676   }
677   AliSysInfo::AddStamp("Stop_sdigitization");
678   
679   AliSysInfo::AddStamp("Start_digitization");  
680   // summable digits -> digits  
681   if (!fMakeDigits.IsNull()) {
682     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
683       if (fStopOnError) return kFALSE;
684     }
685    }
686   AliSysInfo::AddStamp("Stop_digitization");
687
688   
689   
690   // hits -> digits
691   if (!fMakeDigitsFromHits.IsNull()) {
692     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
693       AliWarning(Form("Merging and direct creation of digits from hits " 
694                  "was selected for some detectors. "
695                  "No merging will be done for the following detectors: %s",
696                  fMakeDigitsFromHits.Data()));
697     }
698     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
699       if (fStopOnError) return kFALSE;
700     }
701   }
702
703   
704   
705   // digits -> trigger
706   if (!RunTrigger(fMakeTrigger,fMakeDigits)) {
707     if (fStopOnError) return kFALSE;
708   }
709
710   
711   
712   // digits -> raw data
713   if (!fWriteRawData.IsNull()) {
714     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
715                       fDeleteIntermediateFiles,fWriteSelRawData)) {
716       if (fStopOnError) return kFALSE;
717     }
718   }
719
720   
721   
722   // run HLT simulation
723   if (!fRunHLT.IsNull()) {
724     if (!RunHLT()) {
725       if (fStopOnError) return kFALSE;
726     }
727   }
728   
729   //QA
730         if (fRunQA) {
731                 Bool_t rv = RunQA() ; 
732                 if (!rv)
733                         if (fStopOnError) 
734                                 return kFALSE ;         
735         }
736
737   // Cleanup of CDB manager: cache and active storages!
738   AliCDBManager::Instance()->ClearCache();
739
740   return kTRUE;
741 }
742
743 //_____________________________________________________________________________
744 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
745 {
746   // run the trigger
747
748   AliCodeTimerAuto("")
749
750   // initialize CDB storage from external environment
751   // (either CDB manager or AliSimulation setters),
752   // if not already done in RunSimulation()
753   InitCDB();
754   
755   // Set run number in CDBManager from data 
756   // From this point on the run number must be always loaded from data!
757   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
758   
759   // Set CDB lock: from now on it is forbidden to reset the run number
760   // or the default storage or to activate any further storage!
761   SetCDBLock();
762    
763    AliRunLoader* runLoader = LoadRun("READ");
764    if (!runLoader) return kFALSE;
765    TString trconfiguration = config;
766
767    if (trconfiguration.IsNull()) {
768      if (gAlice->GetTriggerDescriptor() != "") {
769        trconfiguration = gAlice->GetTriggerDescriptor();
770      }
771      else
772        AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
773    }
774
775    runLoader->MakeTree( "GG" );
776    AliCentralTrigger* aCTP = runLoader->GetTrigger();
777    // Load Configuration
778    if (!aCTP->LoadConfiguration( trconfiguration ))
779      return kFALSE;
780
781    // digits -> trigger
782    if( !aCTP->RunTrigger( runLoader , detectors ) ) {
783       if (fStopOnError) {
784         //  delete aCTP;
785         return kFALSE;
786       }
787    }
788
789    delete runLoader;
790
791    return kTRUE;
792 }
793
794 //_____________________________________________________________________________
795 Bool_t AliSimulation::WriteTriggerRawData()
796 {
797   // Writes the CTP (trigger) DDL raw data
798   // Details of the format are given in the
799   // trigger TDR - pages 134 and 135.
800   AliCTPRawData writer;
801   writer.RawData();
802
803   return kTRUE;
804 }
805
806 //_____________________________________________________________________________
807 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
808 {
809 // run the generation and simulation
810
811   AliCodeTimerAuto("")
812
813   // initialize CDB storage and run number from external environment
814   // (either CDB manager or AliSimulation setters)
815   InitCDB();
816   InitRunNumber();
817   SetCDBLock();
818   
819   if (!gAlice) {
820     AliError("no gAlice object. Restart aliroot and try again.");
821     return kFALSE;
822   }
823   if (gAlice->Modules()->GetEntries() > 0) {
824     AliError("gAlice was already run. Restart aliroot and try again.");
825     return kFALSE;
826   }
827
828   AliInfo(Form("initializing gAlice with config file %s",
829           fConfigFileName.Data()));
830   StdoutToAliInfo(StderrToAliError(
831     gAlice->Init(fConfigFileName.Data());
832   ););
833   
834   // Get the trigger descriptor string
835   // Either from AliSimulation or from
836   // gAlice
837   if (fMakeTrigger.IsNull()) {
838     if (gAlice->GetTriggerDescriptor() != "")
839       fMakeTrigger = gAlice->GetTriggerDescriptor();
840   }
841   else
842     gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
843
844   // Set run number in CDBManager
845   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
846
847   AliRunLoader* runLoader = gAlice->GetRunLoader();
848   if (!runLoader) {
849              AliError(Form("gAlice has no run loader object. "
850                              "Check your config file: %s", fConfigFileName.Data()));
851              return kFALSE;
852   }
853   SetGAliceFile(runLoader->GetFileName());
854       
855   // Misalign geometry
856 #if ROOT_VERSION_CODE < 331527
857   AliGeomManager::SetGeometry(gGeoManager);
858   
859   // Check that the consistency of symbolic names for the activated subdetectors
860   // in the geometry loaded by AliGeomManager
861   TString detsToBeChecked = "";
862   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
863   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
864     AliModule* det = (AliModule*) detArray->At(iDet);
865     if (!det || !det->IsActive()) continue;
866     detsToBeChecked += det->GetName();
867     detsToBeChecked += " ";
868   } // end loop over detectors
869   if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
870     AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
871   MisalignGeometry(runLoader);
872 #endif
873
874 //   AliRunLoader* runLoader = gAlice->GetRunLoader();
875 //   if (!runLoader) {
876 //     AliError(Form("gAlice has no run loader object. "
877 //                   "Check your config file: %s", fConfigFileName.Data()));
878 //     return kFALSE;
879 //   }
880 //   SetGAliceFile(runLoader->GetFileName());
881
882   if (!gAlice->Generator()) {
883     AliError(Form("gAlice has no generator object. "
884                   "Check your config file: %s", fConfigFileName.Data()));
885     return kFALSE;
886   }
887   if (nEvents <= 0) nEvents = fNEvents;
888
889   // get vertex from background file in case of merging
890   if (fUseBkgrdVertex &&
891       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
892     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
893     const char* fileName = ((TObjString*)
894                             (fBkgrdFileNames->At(0)))->GetName();
895     AliInfo(Form("The vertex will be taken from the background "
896                  "file %s with nSignalPerBackground = %d", 
897                  fileName, signalPerBkgrd));
898     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
899     gAlice->Generator()->SetVertexGenerator(vtxGen);
900   }
901
902   if (!fRunSimulation) {
903     gAlice->Generator()->SetTrackingFlag(0);
904   }
905
906   // set the number of events per file for given detectors and data types
907   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
908     if (!fEventsPerFile[i]) continue;
909     const char* detName = fEventsPerFile[i]->GetName();
910     const char* typeName = fEventsPerFile[i]->GetTitle();
911     TString loaderName(detName);
912     loaderName += "Loader";
913     AliLoader* loader = runLoader->GetLoader(loaderName);
914     if (!loader) {
915       AliError(Form("RunSimulation", "no loader for %s found\n"
916                     "Number of events per file not set for %s %s", 
917                     detName, typeName, detName));
918       continue;
919     }
920     AliDataLoader* dataLoader = 
921       loader->GetDataLoader(typeName);
922     if (!dataLoader) {
923       AliError(Form("no data loader for %s found\n"
924                     "Number of events per file not set for %s %s", 
925                     typeName, detName, typeName));
926       continue;
927     }
928     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
929     AliDebug(1, Form("number of events per file set to %d for %s %s",
930                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
931   }
932
933   AliInfo("running gAlice");
934   AliSysInfo::AddStamp("Start_simulation");
935   StdoutToAliInfo(StderrToAliError(
936     gAlice->Run(nEvents);
937   ););
938   AliSysInfo::AddStamp("Stop_simulation");
939   delete runLoader;
940
941   return kTRUE;
942 }
943
944 //_____________________________________________________________________________
945 Bool_t AliSimulation::RunSDigitization(const char* detectors)
946 {
947 // run the digitization and produce summable digits
948   static Int_t eventNr=0;
949   AliCodeTimerAuto("")
950
951   // initialize CDB storage, run number, set CDB lock
952   InitCDB();
953   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
954   SetCDBLock();
955   
956   AliRunLoader* runLoader = LoadRun();
957   if (!runLoader) return kFALSE;
958
959   TString detStr = detectors;
960   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
961   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
962     AliModule* det = (AliModule*) detArray->At(iDet);
963     if (!det || !det->IsActive()) continue;
964     if (IsSelected(det->GetName(), detStr)) {
965       AliInfo(Form("creating summable digits for %s", det->GetName()));
966       AliCodeTimerAuto(Form("creating summable digits for %s", det->GetName()));
967       det->Hits2SDigits();
968       AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
969     }
970   }
971
972   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
973     AliError(Form("the following detectors were not found: %s",
974                   detStr.Data()));
975     if (fStopOnError) return kFALSE;
976   }
977   eventNr++;
978   delete runLoader;
979
980   return kTRUE;
981 }
982
983
984 //_____________________________________________________________________________
985 Bool_t AliSimulation::RunDigitization(const char* detectors, 
986                                       const char* excludeDetectors)
987 {
988 // run the digitization and produce digits from sdigits
989
990   AliCodeTimerAuto("")
991
992   // initialize CDB storage, run number, set CDB lock
993   InitCDB();
994   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
995   SetCDBLock();
996   
997   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
998   if (gAlice) delete gAlice;
999   gAlice = NULL;
1000
1001   Int_t nStreams = 1;
1002   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1003   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1004   AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
1005   // manager->SetEmbeddingFlag(fEmbeddingFlag);
1006   manager->SetInputStream(0, fGAliceFileName.Data());
1007   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1008     const char* fileName = ((TObjString*)
1009                             (fBkgrdFileNames->At(iStream-1)))->GetName();
1010     manager->SetInputStream(iStream, fileName);
1011   }
1012
1013   TString detStr = detectors;
1014   TString detExcl = excludeDetectors;
1015   manager->GetInputStream(0)->ImportgAlice();
1016   AliRunLoader* runLoader = 
1017     AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
1018   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1019   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1020     AliModule* det = (AliModule*) detArray->At(iDet);
1021     if (!det || !det->IsActive()) continue;
1022     if (IsSelected(det->GetName(), detStr) && 
1023         !IsSelected(det->GetName(), detExcl)) {
1024       AliDigitizer* digitizer = det->CreateDigitizer(manager);
1025       
1026       if (!digitizer) {
1027         AliError(Form("no digitizer for %s", det->GetName()));
1028         if (fStopOnError) return kFALSE;
1029       } else {
1030         digitizer->SetRegionOfInterest(fRegionOfInterest);
1031       }
1032     }
1033   }
1034
1035   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1036     AliError(Form("the following detectors were not found: %s", 
1037                   detStr.Data()));
1038     if (fStopOnError) return kFALSE;
1039   }
1040
1041   if (!manager->GetListOfTasks()->IsEmpty()) {
1042     AliInfo("executing digitization");
1043     manager->Exec("");
1044   }
1045
1046   delete manager;
1047
1048   return kTRUE;
1049 }
1050
1051 //_____________________________________________________________________________
1052 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1053 {
1054 // run the digitization and produce digits from hits
1055
1056   AliCodeTimerAuto("")
1057
1058   // initialize CDB storage, run number, set CDB lock
1059   InitCDB();
1060   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1061   SetCDBLock();
1062   
1063   AliRunLoader* runLoader = LoadRun("READ");
1064   if (!runLoader) return kFALSE;
1065
1066   TString detStr = detectors;
1067   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1068   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1069     AliModule* det = (AliModule*) detArray->At(iDet);
1070     if (!det || !det->IsActive()) continue;
1071     if (IsSelected(det->GetName(), detStr)) {
1072       AliInfo(Form("creating digits from hits for %s", det->GetName()));
1073       det->Hits2Digits();
1074     }
1075   }
1076
1077   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1078     AliError(Form("the following detectors were not found: %s", 
1079                   detStr.Data()));
1080     if (fStopOnError) return kFALSE;
1081   }
1082
1083   delete runLoader;
1084   //PH Temporary fix to avoid interference with the PHOS loder/getter
1085   //PH The problem has to be solved in more general way 09/06/05
1086
1087   return kTRUE;
1088 }
1089
1090 //_____________________________________________________________________________
1091 Bool_t AliSimulation::WriteRawData(const char* detectors, 
1092                                    const char* fileName,
1093                                    Bool_t deleteIntermediateFiles,
1094                                    Bool_t selrawdata)
1095 {
1096 // convert the digits to raw data
1097 // First DDL raw data files for the given detectors are created.
1098 // If a file name is given, the DDL files are then converted to a DATE file.
1099 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
1100 // afterwards.
1101 // If the file name has the extension ".root", the DATE file is converted
1102 // to a root file.
1103 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1104 // 'selrawdata' flag can be used to enable writing of detectors raw data
1105 // accoring to the trigger cluster.
1106
1107   AliCodeTimerAuto("")
1108   
1109   TString detStr = detectors;
1110   if (IsSelected("HLT", detStr))
1111   {
1112         // Do nothing. "HLT" will be removed from detStr because the HLT raw
1113         // data files are generated in RunHLT.
1114   }
1115
1116   if (!WriteRawFiles(detStr.Data())) {
1117     if (fStopOnError) return kFALSE;
1118   }
1119
1120   TString dateFileName(fileName);
1121   if (!dateFileName.IsNull()) {
1122     Bool_t rootOutput = dateFileName.EndsWith(".root");
1123     if (rootOutput) dateFileName += ".date";
1124     TString selDateFileName;
1125     if (selrawdata) {
1126       selDateFileName = "selected.";
1127       selDateFileName+= dateFileName;
1128     }
1129     if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1130       if (fStopOnError) return kFALSE;
1131     }
1132     if (deleteIntermediateFiles) {
1133       AliRunLoader* runLoader = LoadRun("READ");
1134       if (runLoader) for (Int_t iEvent = 0; 
1135                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1136         char command[256];
1137         sprintf(command, "rm -r raw%d", iEvent);
1138         gSystem->Exec(command);
1139       }
1140     }
1141
1142     if (rootOutput) {
1143       if (!ConvertDateToRoot(dateFileName, fileName)) {
1144         if (fStopOnError) return kFALSE;
1145       }
1146       if (deleteIntermediateFiles) {
1147         gSystem->Unlink(dateFileName);
1148       }
1149       if (selrawdata) {
1150         TString selFileName = "selected.";
1151         selFileName        += fileName;
1152         if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1153           if (fStopOnError) return kFALSE;
1154         }
1155         if (deleteIntermediateFiles) {
1156           gSystem->Unlink(selDateFileName);
1157         }
1158       }
1159     }
1160   }
1161
1162   return kTRUE;
1163 }
1164
1165 //_____________________________________________________________________________
1166 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1167 {
1168 // convert the digits to raw data DDL files
1169
1170   AliCodeTimerAuto("")
1171   
1172   AliRunLoader* runLoader = LoadRun("READ");
1173   if (!runLoader) return kFALSE;
1174
1175   // write raw data to DDL files
1176   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1177     AliInfo(Form("processing event %d", iEvent));
1178     runLoader->GetEvent(iEvent);
1179     TString baseDir = gSystem->WorkingDirectory();
1180     char dirName[256];
1181     sprintf(dirName, "raw%d", iEvent);
1182     gSystem->MakeDirectory(dirName);
1183     if (!gSystem->ChangeDirectory(dirName)) {
1184       AliError(Form("couldn't change to directory %s", dirName));
1185       if (fStopOnError) return kFALSE; else continue;
1186     }
1187
1188     TString detStr = detectors;
1189     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1190     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1191       AliModule* det = (AliModule*) detArray->At(iDet);
1192       if (!det || !det->IsActive()) continue;
1193       if (IsSelected(det->GetName(), detStr)) {
1194         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1195         det->Digits2Raw();
1196       }
1197     }
1198
1199     if (!WriteTriggerRawData())
1200       if (fStopOnError) return kFALSE;
1201
1202     gSystem->ChangeDirectory(baseDir);
1203     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1204       AliError(Form("the following detectors were not found: %s", 
1205                     detStr.Data()));
1206       if (fStopOnError) return kFALSE;
1207     }
1208   }
1209
1210   delete runLoader;
1211   
1212   return kTRUE;
1213 }
1214
1215 //_____________________________________________________________________________
1216 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1217                                             const char* selDateFileName)
1218 {
1219 // convert raw data DDL files to a DATE file with the program "dateStream"
1220 // The second argument is not empty when the user decides to write
1221 // the detectors raw data according to the trigger cluster.
1222
1223   AliCodeTimerAuto("")
1224   
1225   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1226   if (!path) {
1227     AliError("the program dateStream was not found");
1228     if (fStopOnError) return kFALSE;
1229   } else {
1230     delete[] path;
1231   }
1232
1233   AliRunLoader* runLoader = LoadRun("READ");
1234   if (!runLoader) return kFALSE;
1235
1236   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1237   Bool_t selrawdata = kFALSE;
1238   if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1239
1240   char command[256];
1241   // Note the option -s. It is used in order to avoid
1242   // the generation of SOR/EOR events.
1243   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1244           dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1245   FILE* pipe = gSystem->OpenPipe(command, "w");
1246
1247   Int_t selEvents = 0;
1248   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1249     fprintf(pipe, "GDC\n");
1250     Float_t ldc = 0;
1251     Int_t prevLDC = -1;
1252
1253     if (selrawdata) {
1254       // Check if the event was triggered by CTP
1255       runLoader->GetEvent(iEvent);
1256       if (!runLoader->LoadTrigger()) {
1257         AliCentralTrigger *aCTP = runLoader->GetTrigger();
1258         if (aCTP->GetClassMask()) selEvents++;
1259       }
1260       else {
1261         AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1262         selrawdata = kFALSE;
1263       }
1264     }
1265
1266     // loop over detectors and DDLs
1267     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1268       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1269
1270         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1271         Int_t ldcID = Int_t(ldc + 0.0001);
1272         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1273
1274         char rawFileName[256];
1275         sprintf(rawFileName, "raw%d/%s", 
1276                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1277
1278         // check existence and size of raw data file
1279         FILE* file = fopen(rawFileName, "rb");
1280         if (!file) continue;
1281         fseek(file, 0, SEEK_END);
1282         unsigned long size = ftell(file);
1283         fclose(file);
1284         if (!size) continue;
1285
1286         if (ldcID != prevLDC) {
1287           fprintf(pipe, " LDC Id %d\n", ldcID);
1288           prevLDC = ldcID;
1289         }
1290         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1291       }
1292     }
1293   }
1294
1295   Int_t result = gSystem->ClosePipe(pipe);
1296
1297   if (!(selrawdata && selEvents > 0)) {
1298     delete runLoader;
1299     return (result == 0);
1300   }
1301
1302   AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1303   
1304   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1305           selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1306   FILE* pipe2 = gSystem->OpenPipe(command, "w");
1307
1308   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1309
1310     // Get the trigger decision and cluster
1311     TString detClust;
1312     runLoader->GetEvent(iEvent);
1313     if (!runLoader->LoadTrigger()) {
1314       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1315       if (aCTP->GetClassMask() == 0) continue;
1316       detClust = aCTP->GetTriggeredDetectors();
1317       AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1318     }
1319
1320     fprintf(pipe2, "GDC\n");
1321     Float_t ldc = 0;
1322     Int_t prevLDC = -1;
1323
1324     // loop over detectors and DDLs
1325     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1326       // Write only raw data from detectors that
1327       // are contained in the trigger cluster(s)
1328       if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1329
1330       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1331
1332         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1333         Int_t ldcID = Int_t(ldc + 0.0001);
1334         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1335
1336         char rawFileName[256];
1337         sprintf(rawFileName, "raw%d/%s", 
1338                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1339
1340         // check existence and size of raw data file
1341         FILE* file = fopen(rawFileName, "rb");
1342         if (!file) continue;
1343         fseek(file, 0, SEEK_END);
1344         unsigned long size = ftell(file);
1345         fclose(file);
1346         if (!size) continue;
1347
1348         if (ldcID != prevLDC) {
1349           fprintf(pipe2, " LDC Id %d\n", ldcID);
1350           prevLDC = ldcID;
1351         }
1352         fprintf(pipe2, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1353       }
1354     }
1355   }
1356
1357   Int_t result2 = gSystem->ClosePipe(pipe2);
1358
1359   delete runLoader;
1360   return ((result == 0) && (result2 == 0));
1361 }
1362
1363 //_____________________________________________________________________________
1364 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1365                                         const char* rootFileName)
1366 {
1367 // convert a DATE file to a root file with the program "alimdc"
1368
1369   // ALIMDC setup
1370   const Int_t kDBSize = 2000000000;
1371   const Int_t kTagDBSize = 1000000000;
1372   const Bool_t kFilter = kFALSE;
1373   const Int_t kCompression = 1;
1374
1375   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1376   if (!path) {
1377     AliError("the program alimdc was not found");
1378     if (fStopOnError) return kFALSE;
1379   } else {
1380     delete[] path;
1381   }
1382
1383   AliInfo(Form("converting DATE file %s to root file %s", 
1384                dateFileName, rootFileName));
1385
1386   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1387   const char* tagDBFS    = "/tmp/mdc1/tags";
1388
1389   // User defined file system locations
1390   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1391     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1392   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1393     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1394   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1395     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1396
1397   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1398   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1399   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1400
1401   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1402   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1403   gSystem->Exec(Form("mkdir %s",tagDBFS));
1404
1405   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1406                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1407   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1408
1409   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1410   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1411   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1412
1413   return (result == 0);
1414 }
1415
1416
1417 //_____________________________________________________________________________
1418 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1419 {
1420 // delete existing run loaders, open a new one and load gAlice
1421
1422   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1423   AliRunLoader* runLoader = 
1424     AliRunLoader::Open(fGAliceFileName.Data(), 
1425                        AliConfig::GetDefaultEventFolderName(), mode);
1426   if (!runLoader) {
1427     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1428     return NULL;
1429   }
1430   runLoader->LoadgAlice();
1431   runLoader->LoadHeader();
1432   gAlice = runLoader->GetAliRun();
1433   if (!gAlice) {
1434     AliError(Form("no gAlice object found in file %s", 
1435                   fGAliceFileName.Data()));
1436     return NULL;
1437   }
1438   return runLoader;
1439 }
1440
1441 //_____________________________________________________________________________
1442 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1443 {
1444 // get or calculate the number of signal events per background event
1445
1446   if (!fBkgrdFileNames) return 1;
1447   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1448   if (nBkgrdFiles == 0) return 1;
1449
1450   // get the number of signal events
1451   if (nEvents <= 0) {
1452     AliRunLoader* runLoader = 
1453         AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1454     if (!runLoader) return 1;
1455     
1456     nEvents = runLoader->GetNumberOfEvents();
1457     delete runLoader;
1458   }
1459
1460   Int_t result = 0;
1461   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1462     // get the number of background events
1463     const char* fileName = ((TObjString*)
1464                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1465     AliRunLoader* runLoader =
1466       AliRunLoader::Open(fileName, "BKGRD");
1467     if (!runLoader) continue;
1468     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1469     delete runLoader;
1470   
1471     // get or calculate the number of signal per background events
1472     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1473     if (nSignalPerBkgrd <= 0) {
1474       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1475     } else if (result && (result != nSignalPerBkgrd)) {
1476       AliInfo(Form("the number of signal events per background event "
1477                    "will be changed from %d to %d for stream %d", 
1478                    nSignalPerBkgrd, result, iBkgrdFile+1));
1479       nSignalPerBkgrd = result;
1480     }
1481
1482     if (!result) result = nSignalPerBkgrd;
1483     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1484       AliWarning(Form("not enough background events (%d) for %d signal events "
1485                       "using %d signal per background events for stream %d",
1486                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1487     }
1488   }
1489
1490   return result;
1491 }
1492
1493 //_____________________________________________________________________________
1494 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1495 {
1496 // check whether detName is contained in detectors
1497 // if yes, it is removed from detectors
1498
1499   // check if all detectors are selected
1500   if ((detectors.CompareTo("ALL") == 0) ||
1501       detectors.BeginsWith("ALL ") ||
1502       detectors.EndsWith(" ALL") ||
1503       detectors.Contains(" ALL ")) {
1504     detectors = "ALL";
1505     return kTRUE;
1506   }
1507
1508   // search for the given detector
1509   Bool_t result = kFALSE;
1510   if ((detectors.CompareTo(detName) == 0) ||
1511       detectors.BeginsWith(detName+" ") ||
1512       detectors.EndsWith(" "+detName) ||
1513       detectors.Contains(" "+detName+" ")) {
1514     detectors.ReplaceAll(detName, "");
1515     result = kTRUE;
1516   }
1517
1518   // clean up the detectors string
1519   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1520   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1521   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1522
1523   return result;
1524 }
1525
1526 //_____________________________________________________________________________
1527 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName) 
1528 {
1529 //
1530 // Steering routine  to convert raw data in directory rawDirectory/ to fake SDigits. 
1531 // These can be used for embedding of MC tracks into RAW data using the standard 
1532 // merging procedure.
1533 //
1534 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1535 //
1536     if (!gAlice) {
1537         AliError("no gAlice object. Restart aliroot and try again.");
1538         return kFALSE;
1539     }
1540     if (gAlice->Modules()->GetEntries() > 0) {
1541         AliError("gAlice was already run. Restart aliroot and try again.");
1542         return kFALSE;
1543     }
1544     
1545     AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1546     StdoutToAliInfo(StderrToAliError(gAlice->Init(fConfigFileName.Data());););
1547 //
1548 //  Initialize CDB     
1549     InitCDB();
1550     //AliCDBManager* man = AliCDBManager::Instance();
1551     //man->SetRun(0); // Should this come from rawdata header ?
1552     
1553     Int_t iDet;
1554     //
1555     // Get the runloader
1556     AliRunLoader* runLoader = gAlice->GetRunLoader();
1557     //
1558     // Open esd file if available
1559     TFile* esdFile = TFile::Open(esdFileName);
1560     Bool_t esdOK = (esdFile != 0);
1561     AliESD* esd = new AliESD;
1562     TTree* treeESD = 0;
1563     if (esdOK) {
1564         treeESD = (TTree*) esdFile->Get("esdTree");
1565         if (!treeESD) {
1566             AliWarning("No ESD tree found");
1567             esdOK = kFALSE;
1568         } else {
1569             treeESD->SetBranchAddress("ESD", &esd);
1570         }
1571     }
1572     //
1573     // Create the RawReader
1574     TString fileName(rawDirectory);
1575     AliRawReader* rawReader = 0x0;
1576     if (fileName.EndsWith("/")) {
1577       rawReader = new AliRawReaderFile(fileName);
1578     } else if (fileName.EndsWith(".root")) {
1579       rawReader = new AliRawReaderRoot(fileName);
1580     } else if (!fileName.IsNull()) {
1581       rawReader = new AliRawReaderDate(fileName);
1582     }
1583 //     if (!fEquipIdMap.IsNull() && fRawReader)
1584 //       fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1585     //
1586     // Get list of detectors
1587     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1588     //
1589     // Get Header
1590     AliHeader* header = runLoader->GetHeader();
1591     //
1592     TString detStr = fMakeSDigits;
1593     // Event loop
1594     Int_t nev = 0;
1595     while(kTRUE) {
1596         if (!(rawReader->NextEvent())) break;
1597         //
1598         // Detector loop
1599         for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1600             AliModule* det = (AliModule*) detArray->At(iDet);
1601             if (!det || !det->IsActive()) continue;
1602             if (IsSelected(det->GetName(), detStr)) {
1603               AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1604               det->Raw2SDigits(rawReader);
1605               rawReader->Reset();
1606             }
1607         } // detectors
1608
1609
1610         //
1611         //  If ESD information available obtain reconstructed vertex and store in header.
1612         if (esdOK) {
1613             treeESD->GetEvent(nev);
1614             const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1615             Double_t position[3];
1616             esdVertex->GetXYZ(position);
1617             AliGenEventHeader* mcHeader = new  AliGenEventHeader("ESD");
1618             TArrayF mcV;
1619             mcV.Set(3);
1620             for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1621             mcHeader->SetPrimaryVertex(mcV);
1622             header->Reset(0,nev);
1623             header->SetGenEventHeader(mcHeader);
1624             printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1625         }
1626         nev++;
1627 //
1628 //      Finish the event
1629         runLoader->TreeE()->Fill();
1630         runLoader->SetNextEvent();
1631     } // events
1632  
1633     delete rawReader;
1634 //
1635 //  Finish the run 
1636     runLoader->CdGAFile();
1637     runLoader->WriteHeader("OVERWRITE");
1638     runLoader->WriteRunLoader();
1639
1640     return kTRUE;
1641 }
1642
1643 //_____________________________________________________________________________
1644 Int_t AliSimulation::GetDetIndex(const char* detector)
1645 {
1646   // return the detector index corresponding to detector
1647   Int_t index = -1 ; 
1648   for (index = 0; index < fgkNDetectors ; index++) {
1649     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1650           break ; 
1651   }     
1652   return index ; 
1653 }
1654
1655 //_____________________________________________________________________________
1656 Bool_t AliSimulation::RunHLT()
1657 {
1658   // Run the HLT simulation
1659   // HLT simulation is implemented in HLT/sim/AliHLTSimulation
1660   // Disabled if fRunHLT is empty, default vaule is "default".
1661   // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
1662   // The default simulation depends on the HLT component libraries and their
1663   // corresponding agents which define components and chains to run. See
1664   // http://web.ift.uib.no/~kjeks/doc/alice-hlt/
1665   // http://web.ift.uib.no/~kjeks/doc/alice-hlt/classAliHLTModuleAgent.html
1666   //
1667   // The libraries to be loaded can be specified as an option.
1668   // <pre>
1669   // AliSimulation sim;
1670   // sim.SetRunHLT("libAliHLTSample.so");
1671   // </pre>
1672   // will only load <tt>libAliHLTSample.so</tt>
1673
1674   // Other available options:
1675   // \li loglevel=<i>level</i> <br>
1676   //     logging level for this processing
1677   // \li alilog=off
1678   //     disable redirection of log messages to AliLog class
1679   // \li config=<i>macro</i>
1680   //     configuration macro
1681   // \li localrec=<i>configuration</i>
1682   //     comma separated list of configurations to be run during simulation
1683
1684   int iResult=0;
1685   AliRunLoader* pRunLoader = LoadRun("READ");
1686   if (!pRunLoader) return kFALSE;
1687
1688   // initialize CDB storage, run number, set CDB lock
1689   InitCDB();
1690   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1691   SetCDBLock();
1692   
1693   // load the library dynamically
1694   gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1695
1696   // check for the library version
1697   AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
1698   if (!fctVersion) {
1699     AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
1700     return kFALSE;
1701   }
1702   if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
1703     AliError(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
1704     return kFALSE;
1705   }
1706
1707   // print compile info
1708   typedef void (*CompileInfo)( char*& date, char*& time);
1709   CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
1710   if (fctInfo) {
1711     char* date="";
1712     char* time="";
1713     (*fctInfo)(date, time);
1714     if (!date) date="unknown";
1715     if (!time) time="unknown";
1716     AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
1717   } else {
1718     AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
1719   }
1720
1721   // create instance of the HLT simulation
1722   AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
1723   AliHLTSimulation* pHLT=NULL;
1724   if (fctCreate==NULL || (pHLT=(fctCreate()))==NULL) {
1725     AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
1726     return kFALSE;    
1727   }
1728
1729   // init the HLT simulation
1730   TString options;
1731   if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
1732   if (!IsSelected("HLT", fWriteRawData)) {
1733     options+=" writerawfiles=";
1734   } else {
1735     options+=" writerawfiles=HLT";
1736   }
1737   AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
1738   if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
1739     AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
1740   } else {
1741     // run the HLT simulation
1742     AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
1743     if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
1744       AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
1745     }
1746   }
1747
1748   // delete the instance
1749   AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
1750   if (fctDelete==NULL || fctDelete(pHLT)<0) {
1751     AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
1752   }
1753   pHLT=NULL;
1754
1755   return iResult>=0?kTRUE:kFALSE;
1756 }
1757
1758 //_____________________________________________________________________________
1759 Bool_t AliSimulation::RunQA()
1760 {
1761         // run the QA on summable hits, digits or digits
1762
1763         AliQADataMakerSteer qas("sim") ; 
1764     qas.SetRunLoader(gAlice->GetRunLoader()) ;
1765
1766         TString detectorsw("") ;  
1767         Bool_t rv = kTRUE ; 
1768         if (fQATasks.Contains(Form("%d", AliQA::kHITS))) 
1769                 detectorsw =  qas.Run(fQADetectors.Data(), AliQA::kHITS) ; 
1770 //      qas.Reset() ; 
1771         if (fQATasks.Contains(Form("%d", AliQA::kSDIGITS))) 
1772                 detectorsw += qas.Run(fQADetectors.Data(), AliQA::kSDIGITS) ;   
1773 //      qas.Reset() ; 
1774         if (fQATasks.Contains(Form("%d", AliQA::kDIGITS))) 
1775                 detectorsw += qas.Run(fQADetectors.Data(), AliQA::kDIGITS) ;    
1776
1777         if ( detectorsw.IsNull() ) 
1778                 rv = kFALSE ; 
1779         return rv ; 
1780 }
1781
1782 //_____________________________________________________________________________
1783 Bool_t AliSimulation::SetRunQA(TString detAndAction) 
1784 {
1785         // Allows to run QA for a selected set of detectors
1786         // and a selected set of tasks among HITS, SDIGITS and DIGITS
1787         // all selected detectors run the same selected tasks
1788         
1789         if (!detAndAction.Contains(":")) {
1790                 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
1791                 fRunQA = kFALSE ;
1792                 return kFALSE ;                 
1793         }
1794         Int_t colon = detAndAction.Index(":") ; 
1795         fQADetectors = detAndAction(0, colon) ; 
1796         if (fQADetectors.Contains("ALL") )
1797                 fQADetectors = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ; 
1798         fQATasks   = detAndAction(colon+1, detAndAction.Sizeof() ) ; 
1799         if (fQATasks.Contains("ALL") ) {
1800                 fQATasks = Form("%d %d %d", AliQA::kHITS, AliQA::kSDIGITS, AliQA::kDIGITS) ; 
1801         } else {
1802                 fQATasks.ToUpper() ; 
1803                 TString tempo("") ; 
1804                 if ( fQATasks.Contains("HIT") ) 
1805                         tempo = Form("%d ", AliQA::kHITS) ; 
1806                 if ( fQATasks.Contains("SDIGIT") ) 
1807                         tempo += Form("%d ", AliQA::kSDIGITS) ; 
1808                 if ( fQATasks.Contains("DIGIT") ) 
1809                         tempo += Form("%d ", AliQA::kDIGITS) ; 
1810                 fQATasks = tempo ; 
1811                 if (fQATasks.IsNull()) {
1812                         AliInfo("No QA requested\n")  ;
1813                         fRunQA = kFALSE ;
1814                         return kTRUE ; 
1815                 }
1816         }       
1817         TString tempo(fQATasks) ; 
1818     tempo.ReplaceAll(Form("%d", AliQA::kHITS), AliQA::GetTaskName(AliQA::kHITS))        ;
1819     tempo.ReplaceAll(Form("%d", AliQA::kSDIGITS), AliQA::GetTaskName(AliQA::kSDIGITS)) ;        
1820     tempo.ReplaceAll(Form("%d", AliQA::kDIGITS), AliQA::GetTaskName(AliQA::kDIGITS)) ;  
1821         AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;  
1822         fRunQA = kTRUE ;
1823         return kTRUE; 
1824
1825
1826 //_____________________________________________________________________________
1827 void AliSimulation::ProcessEnvironmentVars()
1828 {
1829 // Extract run number and random generator seed from env variables
1830
1831     AliInfo("Processing environment variables");
1832     
1833     // Random Number seed
1834     
1835     // first check that seed is not already set
1836     if (fSeed == 0) {
1837         if (gSystem->Getenv("CONFIG_SEED")) {
1838                 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
1839         }
1840     } else {
1841         if (gSystem->Getenv("CONFIG_SEED")) {
1842                 AliInfo(Form("Seed for random number generation already set (%d)"
1843                              ": CONFIG_SEED variable ignored!", fSeed));
1844         }
1845     }
1846    
1847     AliInfo(Form("Seed for random number generation = %d ", fSeed)); 
1848
1849     // Run Number
1850     
1851     // first check that run number is not already set
1852     if(fRun < 0) {    
1853         if (gSystem->Getenv("DC_RUN")) {
1854                 fRun = atoi(gSystem->Getenv("DC_RUN"));
1855         }
1856     } else {
1857         if (gSystem->Getenv("DC_RUN")) {
1858                 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
1859         }
1860     }
1861     
1862     AliInfo(Form("Run number = %d", fRun)); 
1863 }
1864
1865