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