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