]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
Added new materials
[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   if (!WriteRawFiles(detectors)) {
1110     if (fStopOnError) return kFALSE;
1111   }
1112
1113   TString dateFileName(fileName);
1114   if (!dateFileName.IsNull()) {
1115     Bool_t rootOutput = dateFileName.EndsWith(".root");
1116     if (rootOutput) dateFileName += ".date";
1117     TString selDateFileName;
1118     if (selrawdata) {
1119       selDateFileName = "selected.";
1120       selDateFileName+= dateFileName;
1121     }
1122     if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1123       if (fStopOnError) return kFALSE;
1124     }
1125     if (deleteIntermediateFiles) {
1126       AliRunLoader* runLoader = LoadRun("READ");
1127       if (runLoader) for (Int_t iEvent = 0; 
1128                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1129         char command[256];
1130         sprintf(command, "rm -r raw%d", iEvent);
1131         gSystem->Exec(command);
1132       }
1133     }
1134
1135     if (rootOutput) {
1136       if (!ConvertDateToRoot(dateFileName, fileName)) {
1137         if (fStopOnError) return kFALSE;
1138       }
1139       if (deleteIntermediateFiles) {
1140         gSystem->Unlink(dateFileName);
1141       }
1142       if (selrawdata) {
1143         TString selFileName = "selected.";
1144         selFileName        += fileName;
1145         if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1146           if (fStopOnError) return kFALSE;
1147         }
1148         if (deleteIntermediateFiles) {
1149           gSystem->Unlink(selDateFileName);
1150         }
1151       }
1152     }
1153   }
1154
1155   return kTRUE;
1156 }
1157
1158 //_____________________________________________________________________________
1159 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1160 {
1161 // convert the digits to raw data DDL files
1162
1163   AliCodeTimerAuto("")
1164   
1165   AliRunLoader* runLoader = LoadRun("READ");
1166   if (!runLoader) return kFALSE;
1167
1168   // write raw data to DDL files
1169   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1170     AliInfo(Form("processing event %d", iEvent));
1171     runLoader->GetEvent(iEvent);
1172     TString baseDir = gSystem->WorkingDirectory();
1173     char dirName[256];
1174     sprintf(dirName, "raw%d", iEvent);
1175     gSystem->MakeDirectory(dirName);
1176     if (!gSystem->ChangeDirectory(dirName)) {
1177       AliError(Form("couldn't change to directory %s", dirName));
1178       if (fStopOnError) return kFALSE; else continue;
1179     }
1180
1181     TString detStr = detectors;
1182     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1183     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1184       AliModule* det = (AliModule*) detArray->At(iDet);
1185       if (!det || !det->IsActive()) continue;
1186       if (IsSelected(det->GetName(), detStr)) {
1187         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1188         det->Digits2Raw();
1189       }
1190     }
1191
1192     if (!WriteTriggerRawData())
1193       if (fStopOnError) return kFALSE;
1194
1195     gSystem->ChangeDirectory(baseDir);
1196     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1197       AliError(Form("the following detectors were not found: %s", 
1198                     detStr.Data()));
1199       if (fStopOnError) return kFALSE;
1200     }
1201   }
1202
1203   delete runLoader;
1204   
1205   return kTRUE;
1206 }
1207
1208 //_____________________________________________________________________________
1209 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1210                                             const char* selDateFileName)
1211 {
1212 // convert raw data DDL files to a DATE file with the program "dateStream"
1213 // The second argument is not empty when the user decides to write
1214 // the detectors raw data according to the trigger cluster.
1215
1216   AliCodeTimerAuto("")
1217   
1218   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1219   if (!path) {
1220     AliError("the program dateStream was not found");
1221     if (fStopOnError) return kFALSE;
1222   } else {
1223     delete[] path;
1224   }
1225
1226   AliRunLoader* runLoader = LoadRun("READ");
1227   if (!runLoader) return kFALSE;
1228
1229   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1230   Bool_t selrawdata = kFALSE;
1231   if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1232
1233   char command[256];
1234   // Note the option -s. It is used in order to avoid
1235   // the generation of SOR/EOR events.
1236   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1237           dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1238   FILE* pipe = gSystem->OpenPipe(command, "w");
1239
1240   Int_t selEvents = 0;
1241   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1242     fprintf(pipe, "GDC\n");
1243     Float_t ldc = 0;
1244     Int_t prevLDC = -1;
1245
1246     if (selrawdata) {
1247       // Check if the event was triggered by CTP
1248       runLoader->GetEvent(iEvent);
1249       if (!runLoader->LoadTrigger()) {
1250         AliCentralTrigger *aCTP = runLoader->GetTrigger();
1251         if (aCTP->GetClassMask()) selEvents++;
1252       }
1253       else {
1254         AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1255         selrawdata = kFALSE;
1256       }
1257     }
1258
1259     // loop over detectors and DDLs
1260     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1261       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1262
1263         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1264         Int_t ldcID = Int_t(ldc + 0.0001);
1265         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1266
1267         char rawFileName[256];
1268         sprintf(rawFileName, "raw%d/%s", 
1269                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1270
1271         // check existence and size of raw data file
1272         FILE* file = fopen(rawFileName, "rb");
1273         if (!file) continue;
1274         fseek(file, 0, SEEK_END);
1275         unsigned long size = ftell(file);
1276         fclose(file);
1277         if (!size) continue;
1278
1279         if (ldcID != prevLDC) {
1280           fprintf(pipe, " LDC Id %d\n", ldcID);
1281           prevLDC = ldcID;
1282         }
1283         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1284       }
1285     }
1286   }
1287
1288   Int_t result = gSystem->ClosePipe(pipe);
1289
1290   if (!(selrawdata && selEvents > 0)) {
1291     delete runLoader;
1292     return (result == 0);
1293   }
1294
1295   AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1296   
1297   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1298           selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1299   FILE* pipe2 = gSystem->OpenPipe(command, "w");
1300
1301   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1302
1303     // Get the trigger decision and cluster
1304     TString detClust;
1305     runLoader->GetEvent(iEvent);
1306     if (!runLoader->LoadTrigger()) {
1307       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1308       if (aCTP->GetClassMask() == 0) continue;
1309       detClust = aCTP->GetTriggeredDetectors();
1310       AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1311     }
1312
1313     fprintf(pipe2, "GDC\n");
1314     Float_t ldc = 0;
1315     Int_t prevLDC = -1;
1316
1317     // loop over detectors and DDLs
1318     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1319       // Write only raw data from detectors that
1320       // are contained in the trigger cluster(s)
1321       if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1322
1323       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1324
1325         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1326         Int_t ldcID = Int_t(ldc + 0.0001);
1327         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1328
1329         char rawFileName[256];
1330         sprintf(rawFileName, "raw%d/%s", 
1331                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1332
1333         // check existence and size of raw data file
1334         FILE* file = fopen(rawFileName, "rb");
1335         if (!file) continue;
1336         fseek(file, 0, SEEK_END);
1337         unsigned long size = ftell(file);
1338         fclose(file);
1339         if (!size) continue;
1340
1341         if (ldcID != prevLDC) {
1342           fprintf(pipe2, " LDC Id %d\n", ldcID);
1343           prevLDC = ldcID;
1344         }
1345         fprintf(pipe2, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1346       }
1347     }
1348   }
1349
1350   Int_t result2 = gSystem->ClosePipe(pipe2);
1351
1352   delete runLoader;
1353   return ((result == 0) && (result2 == 0));
1354 }
1355
1356 //_____________________________________________________________________________
1357 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1358                                         const char* rootFileName)
1359 {
1360 // convert a DATE file to a root file with the program "alimdc"
1361
1362   // ALIMDC setup
1363   const Int_t kDBSize = 2000000000;
1364   const Int_t kTagDBSize = 1000000000;
1365   const Bool_t kFilter = kFALSE;
1366   const Int_t kCompression = 1;
1367
1368   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1369   if (!path) {
1370     AliError("the program alimdc was not found");
1371     if (fStopOnError) return kFALSE;
1372   } else {
1373     delete[] path;
1374   }
1375
1376   AliInfo(Form("converting DATE file %s to root file %s", 
1377                dateFileName, rootFileName));
1378
1379   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1380   const char* tagDBFS    = "/tmp/mdc1/tags";
1381
1382   // User defined file system locations
1383   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1384     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1385   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1386     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1387   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1388     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1389
1390   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1391   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1392   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1393
1394   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1395   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1396   gSystem->Exec(Form("mkdir %s",tagDBFS));
1397
1398   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1399                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1400   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1401
1402   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1403   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1404   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1405
1406   return (result == 0);
1407 }
1408
1409
1410 //_____________________________________________________________________________
1411 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1412 {
1413 // delete existing run loaders, open a new one and load gAlice
1414
1415   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1416   AliRunLoader* runLoader = 
1417     AliRunLoader::Open(fGAliceFileName.Data(), 
1418                        AliConfig::GetDefaultEventFolderName(), mode);
1419   if (!runLoader) {
1420     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1421     return NULL;
1422   }
1423   runLoader->LoadgAlice();
1424   runLoader->LoadHeader();
1425   gAlice = runLoader->GetAliRun();
1426   if (!gAlice) {
1427     AliError(Form("no gAlice object found in file %s", 
1428                   fGAliceFileName.Data()));
1429     return NULL;
1430   }
1431   return runLoader;
1432 }
1433
1434 //_____________________________________________________________________________
1435 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1436 {
1437 // get or calculate the number of signal events per background event
1438
1439   if (!fBkgrdFileNames) return 1;
1440   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1441   if (nBkgrdFiles == 0) return 1;
1442
1443   // get the number of signal events
1444   if (nEvents <= 0) {
1445     AliRunLoader* runLoader = 
1446         AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1447     if (!runLoader) return 1;
1448     
1449     nEvents = runLoader->GetNumberOfEvents();
1450     delete runLoader;
1451   }
1452
1453   Int_t result = 0;
1454   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1455     // get the number of background events
1456     const char* fileName = ((TObjString*)
1457                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1458     AliRunLoader* runLoader =
1459       AliRunLoader::Open(fileName, "BKGRD");
1460     if (!runLoader) continue;
1461     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1462     delete runLoader;
1463   
1464     // get or calculate the number of signal per background events
1465     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1466     if (nSignalPerBkgrd <= 0) {
1467       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1468     } else if (result && (result != nSignalPerBkgrd)) {
1469       AliInfo(Form("the number of signal events per background event "
1470                    "will be changed from %d to %d for stream %d", 
1471                    nSignalPerBkgrd, result, iBkgrdFile+1));
1472       nSignalPerBkgrd = result;
1473     }
1474
1475     if (!result) result = nSignalPerBkgrd;
1476     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1477       AliWarning(Form("not enough background events (%d) for %d signal events "
1478                       "using %d signal per background events for stream %d",
1479                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1480     }
1481   }
1482
1483   return result;
1484 }
1485
1486 //_____________________________________________________________________________
1487 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1488 {
1489 // check whether detName is contained in detectors
1490 // if yes, it is removed from detectors
1491
1492   // check if all detectors are selected
1493   if ((detectors.CompareTo("ALL") == 0) ||
1494       detectors.BeginsWith("ALL ") ||
1495       detectors.EndsWith(" ALL") ||
1496       detectors.Contains(" ALL ")) {
1497     detectors = "ALL";
1498     return kTRUE;
1499   }
1500
1501   // search for the given detector
1502   Bool_t result = kFALSE;
1503   if ((detectors.CompareTo(detName) == 0) ||
1504       detectors.BeginsWith(detName+" ") ||
1505       detectors.EndsWith(" "+detName) ||
1506       detectors.Contains(" "+detName+" ")) {
1507     detectors.ReplaceAll(detName, "");
1508     result = kTRUE;
1509   }
1510
1511   // clean up the detectors string
1512   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1513   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1514   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1515
1516   return result;
1517 }
1518
1519 //_____________________________________________________________________________
1520 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName) 
1521 {
1522 //
1523 // Steering routine  to convert raw data in directory rawDirectory/ to fake SDigits. 
1524 // These can be used for embedding of MC tracks into RAW data using the standard 
1525 // merging procedure.
1526 //
1527 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1528 //
1529     if (!gAlice) {
1530         AliError("no gAlice object. Restart aliroot and try again.");
1531         return kFALSE;
1532     }
1533     if (gAlice->Modules()->GetEntries() > 0) {
1534         AliError("gAlice was already run. Restart aliroot and try again.");
1535         return kFALSE;
1536     }
1537     
1538     AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1539     StdoutToAliInfo(StderrToAliError(gAlice->Init(fConfigFileName.Data());););
1540 //
1541 //  Initialize CDB     
1542     InitCDB();
1543     //AliCDBManager* man = AliCDBManager::Instance();
1544     //man->SetRun(0); // Should this come from rawdata header ?
1545     
1546     Int_t iDet;
1547     //
1548     // Get the runloader
1549     AliRunLoader* runLoader = gAlice->GetRunLoader();
1550     //
1551     // Open esd file if available
1552     TFile* esdFile = TFile::Open(esdFileName);
1553     Bool_t esdOK = (esdFile != 0);
1554     AliESD* esd = new AliESD;
1555     TTree* treeESD = 0;
1556     if (esdOK) {
1557         treeESD = (TTree*) esdFile->Get("esdTree");
1558         if (!treeESD) {
1559             AliWarning("No ESD tree found");
1560             esdOK = kFALSE;
1561         } else {
1562             treeESD->SetBranchAddress("ESD", &esd);
1563         }
1564     }
1565     //
1566     // Create the RawReader
1567     TString fileName(rawDirectory);
1568     AliRawReader* rawReader = 0x0;
1569     if (fileName.EndsWith("/")) {
1570       rawReader = new AliRawReaderFile(fileName);
1571     } else if (fileName.EndsWith(".root")) {
1572       rawReader = new AliRawReaderRoot(fileName);
1573     } else if (!fileName.IsNull()) {
1574       rawReader = new AliRawReaderDate(fileName);
1575       rawReader->SelectEvents(7);
1576     }
1577 //     if (!fEquipIdMap.IsNull() && fRawReader)
1578 //       fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1579     //
1580     // Get list of detectors
1581     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1582     //
1583     // Get Header
1584     AliHeader* header = runLoader->GetHeader();
1585     //
1586     TString detStr = fMakeSDigits;
1587     // Event loop
1588     Int_t nev = 0;
1589     while(kTRUE) {
1590         if (!(rawReader->NextEvent())) break;
1591         //
1592         // Detector loop
1593         for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1594             AliModule* det = (AliModule*) detArray->At(iDet);
1595             if (!det || !det->IsActive()) continue;
1596             if (IsSelected(det->GetName(), detStr)) {
1597               AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1598               det->Raw2SDigits(rawReader);
1599               rawReader->Reset();
1600             }
1601         } // detectors
1602
1603
1604         //
1605         //  If ESD information available obtain reconstructed vertex and store in header.
1606         if (esdOK) {
1607             treeESD->GetEvent(nev);
1608             const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1609             Double_t position[3];
1610             esdVertex->GetXYZ(position);
1611             AliGenEventHeader* mcHeader = new  AliGenEventHeader("ESD");
1612             TArrayF mcV;
1613             mcV.Set(3);
1614             for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1615             mcHeader->SetPrimaryVertex(mcV);
1616             header->Reset(0,nev);
1617             header->SetGenEventHeader(mcHeader);
1618             printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1619         }
1620         nev++;
1621 //
1622 //      Finish the event
1623         runLoader->TreeE()->Fill();
1624         runLoader->SetNextEvent();
1625     } // events
1626  
1627     delete rawReader;
1628 //
1629 //  Finish the run 
1630     runLoader->CdGAFile();
1631     runLoader->WriteHeader("OVERWRITE");
1632     runLoader->WriteRunLoader();
1633
1634     return kTRUE;
1635 }
1636
1637 //_____________________________________________________________________________
1638 Int_t AliSimulation::GetDetIndex(const char* detector)
1639 {
1640   // return the detector index corresponding to detector
1641   Int_t index = -1 ; 
1642   for (index = 0; index < fgkNDetectors ; index++) {
1643     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1644           break ; 
1645   }     
1646   return index ; 
1647 }
1648
1649 //_____________________________________________________________________________
1650 Bool_t AliSimulation::RunHLT()
1651 {
1652   // Run the HLT simulation
1653   // HLT simulation is implemented in HLT/sim/AliHLTSimulation
1654   // Disabled if fRunHLT is empty, default vaule is "default".
1655   // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
1656   // The default simulation depends on the HLT component libraries and their
1657   // corresponding agents which define components and chains to run. See
1658   // http://web.ift.uib.no/~kjeks/doc/alice-hlt/
1659   // http://web.ift.uib.no/~kjeks/doc/alice-hlt/classAliHLTModuleAgent.html
1660   //
1661   // The libraries to be loaded can be specified as an option.
1662   // <pre>
1663   // AliSimulation sim;
1664   // sim.SetRunHLT("libAliHLTSample.so");
1665   // </pre>
1666   // will only load <tt>libAliHLTSample.so</tt>
1667
1668   // Other available options:
1669   // \li loglevel=<i>level</i> <br>
1670   //     logging level for this processing
1671   // \li alilog=off
1672   //     disable redirection of log messages to AliLog class
1673   // \li config=<i>macro</i>
1674   //     configuration macro
1675   // \li localrec=<i>configuration</i>
1676   //     comma separated list of configurations to be run during simulation
1677
1678   int iResult=0;
1679   AliRunLoader* pRunLoader = LoadRun("READ");
1680   if (!pRunLoader) return kFALSE;
1681
1682   // initialize CDB storage, run number, set CDB lock
1683   InitCDB();
1684   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1685   SetCDBLock();
1686   
1687   // load the library dynamically
1688   gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1689
1690   // check for the library version
1691   AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
1692   if (!fctVersion) {
1693     AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
1694     return kFALSE;
1695   }
1696   if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
1697     AliError(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
1698     return kFALSE;
1699   }
1700
1701   // print compile info
1702   typedef void (*CompileInfo)( char*& date, char*& time);
1703   CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
1704   if (fctInfo) {
1705     char* date="";
1706     char* time="";
1707     (*fctInfo)(date, time);
1708     if (!date) date="unknown";
1709     if (!time) time="unknown";
1710     AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
1711   } else {
1712     AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
1713   }
1714
1715   // create instance of the HLT simulation
1716   AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
1717   AliHLTSimulation* pHLT=NULL;
1718   if (fctCreate==NULL || (pHLT=(fctCreate()))==NULL) {
1719     AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
1720     return kFALSE;    
1721   }
1722
1723   // init the HLT simulation
1724   TString options;
1725   if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
1726   if (!IsSelected("HLT", fWriteRawData)) {
1727     options+=" writerawfiles=";
1728   } else {
1729     options+=" writerawfiles=HLT";
1730   }
1731   AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
1732   if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
1733     AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
1734   } else {
1735     // run the HLT simulation
1736     AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
1737     if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
1738       AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
1739     }
1740   }
1741
1742   // delete the instance
1743   AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
1744   if (fctDelete==NULL || fctDelete(pHLT)<0) {
1745     AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
1746   }
1747   pHLT=NULL;
1748
1749   return iResult>=0?kTRUE:kFALSE;
1750 }
1751
1752 //_____________________________________________________________________________
1753 Bool_t AliSimulation::RunQA()
1754 {
1755         // run the QA on summable hits, digits or digits
1756
1757         AliQADataMakerSteer qas ; 
1758     qas.SetRunLoader(gAlice->GetRunLoader()) ;
1759
1760         TString detectorsw("") ;  
1761         Bool_t rv = kTRUE ; 
1762         if (fQATasks.Contains(Form("%d", AliQA::kHITS))) 
1763                 detectorsw =  qas.Run(fQADetectors.Data(), AliQA::kHITS) ; 
1764 //      qas.Reset() ; 
1765         if (fQATasks.Contains(Form("%d", AliQA::kSDIGITS))) 
1766                 detectorsw += qas.Run(fQADetectors.Data(), AliQA::kSDIGITS) ;   
1767 //      qas.Reset() ; 
1768         if (fQATasks.Contains(Form("%d", AliQA::kDIGITS))) 
1769                 detectorsw += qas.Run(fQADetectors.Data(), AliQA::kDIGITS) ;    
1770
1771         if ( detectorsw.IsNull() ) 
1772                 rv = kFALSE ; 
1773         return rv ; 
1774 }
1775
1776 //_____________________________________________________________________________
1777 Bool_t AliSimulation::SetRunQA(TString detAndAction) 
1778 {
1779         // Allows to run QA for a selected set of detectors
1780         // and a selected set of tasks among HITS, SDIGITS and DIGITS
1781         // all selected detectors run the same selected tasks
1782         
1783         if (!detAndAction.Contains(":")) {
1784                 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
1785                 fRunQA = kFALSE ;
1786                 return kFALSE ;                 
1787         }
1788         Int_t colon = detAndAction.Index(":") ; 
1789         fQADetectors = detAndAction(0, colon) ; 
1790         if (fQADetectors.Contains("ALL") )
1791                 fQADetectors = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ; 
1792         fQATasks   = detAndAction(colon+1, detAndAction.Sizeof() ) ; 
1793         if (fQATasks.Contains("ALL") ) {
1794                 fQATasks = Form("%d %d %d", AliQA::kHITS, AliQA::kSDIGITS, AliQA::kDIGITS) ; 
1795         } else {
1796                 fQATasks.ToUpper() ; 
1797                 TString tempo("") ; 
1798                 if ( fQATasks.Contains("HIT") ) 
1799                         tempo = Form("%d ", AliQA::kHITS) ; 
1800                 if ( fQATasks.Contains("SDIGIT") ) 
1801                         tempo += Form("%d ", AliQA::kSDIGITS) ; 
1802                 if ( fQATasks.Contains("DIGIT") ) 
1803                         tempo += Form("%d ", AliQA::kDIGITS) ; 
1804                 fQATasks = tempo ; 
1805                 if (fQATasks.IsNull()) {
1806                         AliInfo("No QA requested\n")  ;
1807                         fRunQA = kFALSE ;
1808                         return kTRUE ; 
1809                 }
1810         }       
1811         TString tempo(fQATasks) ; 
1812     tempo.ReplaceAll(Form("%d", AliQA::kHITS), AliQA::GetTaskName(AliQA::kHITS))        ;
1813     tempo.ReplaceAll(Form("%d", AliQA::kSDIGITS), AliQA::GetTaskName(AliQA::kSDIGITS)) ;        
1814     tempo.ReplaceAll(Form("%d", AliQA::kDIGITS), AliQA::GetTaskName(AliQA::kDIGITS)) ;  
1815         AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;  
1816         fRunQA = kTRUE ;
1817         return kTRUE; 
1818
1819
1820 //_____________________________________________________________________________
1821 void AliSimulation::ProcessEnvironmentVars()
1822 {
1823 // Extract run number and random generator seed from env variables
1824
1825     AliInfo("Processing environment variables");
1826     
1827     // Random Number seed
1828     
1829     // first check that seed is not already set
1830     if (fSeed == 0) {
1831         if (gSystem->Getenv("CONFIG_SEED")) {
1832                 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
1833         }
1834     } else {
1835         if (gSystem->Getenv("CONFIG_SEED")) {
1836                 AliInfo(Form("Seed for random number generation already set (%d)"
1837                              ": CONFIG_SEED variable ignored!", fSeed));
1838         }
1839     }
1840    
1841     AliInfo(Form("Seed for random number generation = %d ", fSeed)); 
1842
1843     // Run Number
1844     
1845     // first check that run number is not already set
1846     if(fRun < 0) {    
1847         if (gSystem->Getenv("DC_RUN")) {
1848                 fRun = atoi(gSystem->Getenv("DC_RUN"));
1849         }
1850     } else {
1851         if (gSystem->Getenv("DC_RUN")) {
1852                 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
1853         }
1854     }
1855     
1856     AliInfo(Form("Run number = %d", fRun)); 
1857 }
1858
1859