]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/STEER/AliSimulation.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / STEER / 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: AliSimulation.cxx 63204 2013-06-26 13:33:28Z rgrosso $ */
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 <TFile.h>
109 #include <TGeoGlobalMagField.h>
110 #include <TGeoManager.h>
111 #include <TObjString.h>
112 #include <TROOT.h>
113 #include <TSystem.h>
114 #include <TVirtualMC.h>
115 #include <TVirtualMCApplication.h>
116 #include <TDatime.h>
117 #include <TInterpreter.h>
118
119 #include "AliAlignObj.h"
120 #include "AliCDBEntry.h"
121 #include "AliCDBManager.h"
122 #include "AliGRPManager.h"
123 #include "AliCDBStorage.h"
124 #include "AliCTPRawData.h"
125 #include "AliCentralTrigger.h"
126 #include "AliCentralTrigger.h"
127 #include "AliCodeTimer.h"
128 #include "AliDAQ.h"
129 #include "AliDigitizer.h"
130 #include "AliESDEvent.h"
131 #include "AliGRPObject.h"
132 #include "AliGenEventHeader.h"
133 #include "AliGenerator.h"
134 #include "AliGeomManager.h"
135 #include "AliHLTSimulation.h"
136 #include "AliHeader.h"
137 #include "AliLego.h"
138 #include "AliLegoGenerator.h"
139 #include "AliLog.h"
140 #include "AliMC.h"
141 #include "AliMagF.h"
142 #include "AliModule.h"
143 #include "AliPDG.h"
144 #include "AliRawReaderDate.h"
145 #include "AliRawReaderFile.h"
146 #include "AliRawReaderRoot.h"
147 #include "AliRun.h"
148 #include "AliDigitizationInput.h"
149 #include "AliRunLoader.h"
150 #include "AliStack.h"
151 #include "AliSimulation.h"
152 #include "AliSysInfo.h"
153 #include "AliVertexGenFile.h"
154
155 using std::ofstream;
156 ClassImp(AliSimulation)
157
158 AliSimulation *AliSimulation::fgInstance = 0;
159 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE","AD","MFT","HLT"};
160
161 //_____________________________________________________________________________
162 AliSimulation::AliSimulation(const char* configFileName,
163                              const char* name, const char* title) :
164   TNamed(name, title),
165
166   fRunGeneratorOnly(kFALSE),
167   fRunGeneration(kTRUE),
168   fRunSimulation(kTRUE),
169   fLoadAlignFromCDB(kTRUE),
170   fLoadAlObjsListOfDets("ALL"),
171   fMakeSDigits("ALL"),
172   fMakeDigits("ALL"),
173   fTriggerConfig(""),
174   fMakeDigitsFromHits(""),
175   fWriteRawData(""),
176   fRawDataFileName(""),
177   fDeleteIntermediateFiles(kFALSE),
178   fWriteSelRawData(kFALSE),
179   fStopOnError(kFALSE),
180   fUseMonitoring(kFALSE),
181   fNEvents(1),
182   fConfigFileName(configFileName),
183   fGAliceFileName("galice.root"),
184   fEventsPerFile(),
185   fBkgrdFileNames(NULL),
186   fAlignObjArray(NULL),
187   fUseBkgrdVertex(kTRUE),
188   fRegionOfInterest(kFALSE),
189   fCDBUri(""),
190   fQARefUri(""), 
191   fSpecCDBUri(),
192   fRun(-1),
193   fSeed(0),
194   fInitCDBCalled(kFALSE),
195   fInitRunNumberCalled(kFALSE),
196   fSetRunNumberFromDataCalled(kFALSE),
197   fEmbeddingFlag(kFALSE),
198   fLego(NULL),
199   fKey(0),
200   fUseVertexFromCDB(0),
201   fUseMagFieldFromGRP(0),
202   fGRPWriteLocation(Form("local://%s", gSystem->pwd())),
203   fUseTimeStampFromCDB(0),
204   fTimeStart(0),
205   fTimeEnd(0),
206   fQADetectors("ALL"),                  
207   fQATasks("ALL"),      
208   fRunQA(kTRUE), 
209   fEventSpecie(AliRecoParam::kDefault),
210   fWriteQAExpertData(kTRUE), 
211   fGeometryFile(),
212   fRunHLT("default"),
213   fpHLT(NULL),
214   fWriteGRPEntry(kTRUE)
215 {
216 // create simulation object with default parameters
217   fgInstance = this;
218   SetGAliceFile("galice.root");
219   
220 // for QA
221         AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ; 
222         qam->SetActiveDetectors(fQADetectors) ; 
223         fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ; 
224         qam->SetTasks(fQATasks) ;       
225 }
226
227 //_____________________________________________________________________________
228 AliSimulation::~AliSimulation()
229 {
230 // clean up
231
232   fEventsPerFile.Delete();
233 //  if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
234 //  delete fAlignObjArray; fAlignObjArray=0;
235
236   if (fBkgrdFileNames) {
237     fBkgrdFileNames->Delete();
238     delete fBkgrdFileNames;
239   }
240
241   fSpecCDBUri.Delete();
242   if (fgInstance==this) fgInstance = 0;
243
244   AliQAManager::QAManager()->ShowQA() ; 
245   AliQAManager::Destroy() ;     
246   AliCodeTimer::Instance()->Print();
247 }
248
249
250 //_____________________________________________________________________________
251 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
252 {
253 // set the number of events for one run
254
255   fNEvents = nEvents;
256 }
257
258 //_____________________________________________________________________________
259 void AliSimulation::InitQA()
260 {
261   // activate a default CDB storage
262   // First check if we have any CDB storage set, because it is used 
263   // to retrieve the calibration and alignment constants
264   
265   if (fInitCDBCalled) return;
266   fInitCDBCalled = kTRUE;
267
268   AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ; 
269   qam->SetActiveDetectors(fQADetectors) ; 
270   fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ; 
271   qam->SetTasks(fQATasks) ;
272         if (fWriteQAExpertData)
273     qam->SetWriteExpert() ; 
274   
275   if (qam->IsDefaultStorageSet()) {
276     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
277     AliWarning("Default QA reference storage has been already set !");
278     AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fQARefUri.Data()));
279     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
280     fQARefUri = qam->GetDefaultStorage()->GetURI();
281   } else {
282       if (fQARefUri.Length() > 0) {
283         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
284         AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
285         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
286       } else {
287         fQARefUri="local://$ALICE_ROOT/QARef";
288         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
289         AliWarning("Default QA reference storage not yet set !!!!");
290         AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
291         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
292       }
293     qam->SetDefaultStorage(fQARefUri);
294   }
295 }
296
297 //_____________________________________________________________________________
298 void AliSimulation::InitCDB()
299 {
300 // activate a default CDB storage
301 // First check if we have any CDB storage set, because it is used 
302 // to retrieve the calibration and alignment constants
303
304   if (fInitCDBCalled) return;
305   fInitCDBCalled = kTRUE;
306
307   AliCDBManager* man = AliCDBManager::Instance();
308   if (man->IsDefaultStorageSet())
309   {
310     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
311     AliWarning("Default CDB storage has been already set !");
312     AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
313     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
314     fCDBUri = man->GetDefaultStorage()->GetURI();
315   }
316   else {
317     if (fCDBUri.Length() > 0) 
318     {
319         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
320         AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
321         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
322     } else {
323         fCDBUri="local://$ALICE_ROOT/OCDB";
324         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
325         AliWarning("Default CDB storage not yet set !!!!");
326         AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
327         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
328                 
329     }
330     man->SetDefaultStorage(fCDBUri);
331   }
332
333   // Now activate the detector specific CDB storage locations
334   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
335     TObject* obj = fSpecCDBUri[i];
336     if (!obj) continue;
337     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
338     AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
339     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
340     man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
341   }
342       
343 }
344
345 //_____________________________________________________________________________
346 void AliSimulation::InitRunNumber(){
347 // check run number. If not set, set it to 0 !!!!
348   
349   if (fInitRunNumberCalled) return;
350   fInitRunNumberCalled = kTRUE;
351   
352   AliCDBManager* man = AliCDBManager::Instance();
353   if (man->GetRun() >= 0)
354   {
355         AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
356                         "Use external variable DC_RUN or AliSimulation::SetRun()!"));
357   }
358     
359   if(fRun >= 0) {
360         AliDebug(2, Form("Setting CDB run number to: %d",fRun));
361   } else {
362         fRun=0;
363         AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
364                         fRun));
365   }
366   man->SetRun(fRun);
367
368   man->Print();
369
370 }
371
372 //_____________________________________________________________________________
373 void AliSimulation::SetCDBLock() {
374   // Set CDB lock: from now on it is forbidden to reset the run number
375   // or the default storage or to activate any further storage!
376   
377   ULong64_t key = AliCDBManager::Instance()->SetLock(1);
378   if (key) fKey = key;
379 }
380
381 //_____________________________________________________________________________
382 void AliSimulation::SetDefaultStorage(const char* uri) {
383   // Store the desired default CDB storage location
384   // Activate it later within the Run() method
385   
386   fCDBUri = uri;
387   
388 }
389
390 //_____________________________________________________________________________
391 void AliSimulation::SetQARefDefaultStorage(const char* uri) {
392   // Store the desired default CDB storage location
393   // Activate it later within the Run() method
394   
395   fQARefUri = uri;
396   AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
397 }
398
399 //_____________________________________________________________________________
400 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
401 // Store a detector-specific CDB storage location
402 // Activate it later within the Run() method
403
404   AliCDBPath aPath(calibType);
405   if(!aPath.IsValid()){
406         AliError(Form("Not a valid path: %s", calibType));
407         return;
408   }
409
410   TObject* obj = fSpecCDBUri.FindObject(calibType);
411   if (obj) fSpecCDBUri.Remove(obj);
412   fSpecCDBUri.Add(new TNamed(calibType, uri));
413
414 }
415
416 //_____________________________________________________________________________
417 void AliSimulation::SetRunNumber(Int_t run)
418 {
419 // sets run number
420 // Activate it later within the Run() method
421
422         fRun = run;
423 }
424
425 //_____________________________________________________________________________
426 void AliSimulation::SetSeed(Int_t seed)
427 {
428 // sets seed number
429 // Activate it later within the Run() method
430
431         fSeed = seed;
432 }
433
434 //_____________________________________________________________________________
435 Bool_t AliSimulation::SetRunNumberFromData()
436 {
437   // Set the CDB manager run number
438   // The run number is retrieved from gAlice
439
440     if (fSetRunNumberFromDataCalled) return kTRUE;
441     fSetRunNumberFromDataCalled = kTRUE;    
442   
443     AliCDBManager* man = AliCDBManager::Instance();
444     Int_t runData = -1, runCDB = -1;
445   
446     AliRunLoader* runLoader = LoadRun("READ");
447     if (!runLoader) return kFALSE;
448     else {
449         runData = runLoader->GetHeader()->GetRun();
450         delete runLoader;
451     }
452   
453     runCDB = man->GetRun();
454     if(runCDB >= 0) {
455         if (runCDB != runData) {
456                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
457                 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
458                 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
459                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");    
460         }
461         
462     }
463       
464     man->SetRun(runData);
465     fRun = runData;
466     
467     if(man->GetRun() < 0) {
468         AliError("Run number not properly initalized!");
469         return kFALSE;
470     }
471   
472     man->Print();
473     
474     return kTRUE;
475 }
476
477 //_____________________________________________________________________________
478 void AliSimulation::SetConfigFile(const char* fileName)
479 {
480 // set the name of the config file
481
482   fConfigFileName = fileName;
483 }
484
485 //_____________________________________________________________________________
486 void AliSimulation::SetGAliceFile(const char* fileName)
487 {
488 // set the name of the galice file
489 // the path is converted to an absolute one if it is relative
490
491   fGAliceFileName = fileName;
492   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
493     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
494                                                 fGAliceFileName);
495     fGAliceFileName = absFileName;
496     delete[] absFileName;
497   }
498
499   AliDebug(2, Form("galice file name set to %s", fileName));
500 }
501
502 //_____________________________________________________________________________
503 void AliSimulation::SetEventsPerFile(const char* detector, const char* type, 
504                                      Int_t nEvents)
505 {
506 // set the number of events per file for the given detector and data type
507 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
508
509   TNamed* obj = new TNamed(detector, type);
510   obj->SetUniqueID(nEvents);
511   fEventsPerFile.Add(obj);
512 }
513
514 //_____________________________________________________________________________
515 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
516 {
517   // Read the alignment objects from CDB.
518   // Each detector is supposed to have the
519   // alignment objects in DET/Align/Data CDB path.
520   // All the detector objects are then collected,
521   // sorted by geometry level (starting from ALIC) and
522   // then applied to the TGeo geometry.
523   // Finally an overlaps check is performed.
524
525   if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
526     AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
527     return kFALSE;
528   }  
529   
530   // initialize CDB storage, run number, set CDB lock
531   InitCDB();
532 //  if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
533   SetCDBLock();
534     
535   Bool_t delRunLoader = kFALSE;
536   if (!runLoader) {
537     runLoader = LoadRun("READ");
538     if (!runLoader) return kFALSE;
539     delRunLoader = kTRUE;
540   }
541   
542   // Export ideal geometry 
543   if(!IsGeometryFromFile()) AliGeomManager::GetGeometry()->Export("geometry.root");
544
545   // Load alignment data from CDB and apply to geometry through AliGeomManager
546   if(fLoadAlignFromCDB){
547     
548     TString detStr = fLoadAlObjsListOfDets;
549     TString loadAlObjsListOfDets = "";
550     
551     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
552     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
553       AliModule* det = (AliModule*) detArray->At(iDet);
554       if (!det || !det->IsActive()) continue;
555       if (IsSelected(det->GetName(), detStr)) {
556         //add det to list of dets to be aligned from CDB
557         loadAlObjsListOfDets += det->GetName();
558         loadAlObjsListOfDets += " ";
559       }
560     } // end loop over detectors
561     loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
562     AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
563   }else{
564     // Check if the array with alignment objects was
565     // provided by the user. If yes, apply the objects
566     // to the present TGeo geometry
567     if (fAlignObjArray) {
568       if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
569         AliError("The misalignment of one or more volumes failed!"
570                  "Compare the list of simulated detectors and the list of detector alignment data!");
571         if (delRunLoader) delete runLoader;
572         return kFALSE;
573       }
574     }
575   }
576
577   // Update the internal geometry of modules (ITS needs it)
578   TString detStr = fLoadAlObjsListOfDets;
579   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
580   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
581
582     AliModule* det = (AliModule*) detArray->At(iDet);
583     if (!det || !det->IsActive()) continue;
584     if (IsSelected(det->GetName(), detStr)) {
585       det->UpdateInternalGeometry();
586     }
587   } // end loop over detectors
588
589
590   if (delRunLoader) delete runLoader;
591
592   return kTRUE;
593 }
594
595 //_____________________________________________________________________________
596 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
597 {
598 // add a file with background events for merging
599
600   TObjString* fileNameStr = new TObjString(fileName);
601   fileNameStr->SetUniqueID(nSignalPerBkgrd);
602   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
603   fBkgrdFileNames->Add(fileNameStr);
604 }
605
606 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
607 {
608 // add a file with background events for embeddin
609   MergeWith(fileName, nSignalPerBkgrd);
610   fEmbeddingFlag = kTRUE;
611 }
612
613 //_____________________________________________________________________________
614 Bool_t AliSimulation::Run(Int_t nEvents)
615 {
616 // run the generation, simulation and digitization
617
618  
619   AliCodeTimerAuto("",0)
620   AliSysInfo::AddStamp("Start_Run");
621   
622   // Load run number and seed from environmental vars
623   ProcessEnvironmentVars();
624   AliSysInfo::AddStamp("ProcessEnvironmentVars");
625
626   gRandom->SetSeed(fSeed);
627    
628   if (nEvents > 0) fNEvents = nEvents;
629
630   // Run generator-only code on demand
631   if (fRunGeneratorOnly)
632   {
633     if(!RunGeneratorOnly())
634     {
635       if (fStopOnError) return kFALSE;
636     }
637     else
638       return kTRUE;
639   }
640
641   // create and setup the HLT instance
642   if (!fRunHLT.IsNull() && !CreateHLT()) {
643     if (fStopOnError) return kFALSE;
644     // disable HLT
645     fRunHLT="";
646   }
647   
648   // generation and simulation -> hits
649   if (fRunGeneration) {
650     if (!RunSimulation()) if (fStopOnError) return kFALSE;
651   }
652   AliSysInfo::AddStamp("RunSimulation");
653            
654   // initialize CDB storage from external environment
655   // (either CDB manager or AliSimulation setters),
656   // if not already done in RunSimulation()
657   InitCDB();
658   AliSysInfo::AddStamp("InitCDB");
659   
660   // Set run number in CDBManager from data 
661   // From this point on the run number must be always loaded from data!
662   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
663   
664   // Set CDB lock: from now on it is forbidden to reset the run number
665   // or the default storage or to activate any further storage!
666   SetCDBLock();
667
668   // If RunSimulation was not called, load the geometry and misalign it
669   if (!AliGeomManager::GetGeometry()) {
670     // Initialize the geometry manager
671     AliGeomManager::LoadGeometry("geometry.root");
672     AliSysInfo::AddStamp("GetGeometry");
673 //    // Check that the consistency of symbolic names for the activated subdetectors
674 //    // in the geometry loaded by AliGeomManager
675 //    AliRunLoader* runLoader = LoadRun("READ");
676 //    if (!runLoader) return kFALSE;
677 //
678 //    TString detsToBeChecked = "";
679 //    TObjArray* detArray = runLoader->GetAliRun()->Detectors();
680 //    for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
681 //      AliModule* det = (AliModule*) detArray->At(iDet);
682 //      if (!det || !det->IsActive()) continue;
683 //      detsToBeChecked += det->GetName();
684 //      detsToBeChecked += " ";
685 //    } // end loop over detectors
686 //    if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
687     if(!AliGeomManager::CheckSymNamesLUT("ALL"))
688         AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
689         
690     if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
691     // Misalign geometry
692     if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
693   }
694   AliSysInfo::AddStamp("MissalignGeometry");
695
696
697   // hits -> summable digits
698   AliSysInfo::AddStamp("Start_sdigitization");
699   if (!fMakeSDigits.IsNull()) {
700     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
701  
702   }
703   AliSysInfo::AddStamp("Stop_sdigitization");
704   
705   AliSysInfo::AddStamp("Start_digitization");  
706   // summable digits -> digits  
707   if (!fMakeDigits.IsNull()) {
708     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
709       if (fStopOnError) return kFALSE;
710     }
711    }
712   AliSysInfo::AddStamp("Stop_digitization");
713
714   
715   
716   // hits -> digits
717   if (!fMakeDigitsFromHits.IsNull()) {
718     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
719       AliWarning(Form("Merging and direct creation of digits from hits " 
720                  "was selected for some detectors. "
721                  "No merging will be done for the following detectors: %s",
722                  fMakeDigitsFromHits.Data()));
723     }
724     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
725       if (fStopOnError) return kFALSE;
726     }
727   }
728
729   AliSysInfo::AddStamp("Hits2Digits");
730   
731   
732   // digits -> trigger
733   if (!fTriggerConfig.IsNull() && !RunTrigger(fTriggerConfig,fMakeDigits)) {
734     if (fStopOnError) return kFALSE;
735   }
736
737   AliSysInfo::AddStamp("RunTrigger");
738   
739   
740   // digits -> raw data
741   if (!fWriteRawData.IsNull()) {
742     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
743                       fDeleteIntermediateFiles,fWriteSelRawData)) {
744       if (fStopOnError) return kFALSE;
745     }
746   }
747
748   AliSysInfo::AddStamp("WriteRaw");
749   
750   // run HLT simulation on simulated digit data if raw data is not
751   // simulated, otherwise its called as part of WriteRawData
752   if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
753     if (!RunHLT()) {
754       if (fStopOnError) return kFALSE;
755     }
756   }
757
758   AliSysInfo::AddStamp("RunHLT");
759   
760   //QA
761   if (fRunQA) {
762       Bool_t rv = RunQA() ; 
763       if (!rv)
764           if (fStopOnError) 
765               return kFALSE ;
766   }
767
768   AliSysInfo::AddStamp("RunQA");
769
770   TString snapshotFileOut("");
771   if(TString(gSystem->Getenv("OCDB_SNAPSHOT_CREATE")) == TString("kTRUE")){ 
772       AliInfo(" ******** Creating the snapshot! *********");
773       TString snapshotFile(gSystem->Getenv("OCDB_SNAPSHOT_FILENAME")); 
774       if(!(snapshotFile.IsNull() || snapshotFile.IsWhitespace())) 
775           snapshotFileOut = snapshotFile;
776       else 
777           snapshotFileOut="OCDB.root"; 
778       AliCDBManager::Instance()->DumpToSnapshotFile(snapshotFileOut.Data(),kFALSE); 
779   }
780
781   // Cleanup of CDB manager: cache and active storages!
782   AliCDBManager::Instance()->ClearCache();
783
784   return kTRUE;
785 }
786
787 //_______________________________________________________________________
788 Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
789                      Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
790                      Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
791 {
792   //
793   // Generates lego plots of:
794   //    - radiation length map phi vs theta
795   //    - radiation length map phi vs eta
796   //    - interaction length map
797   //    - g/cm2 length map
798   //
799   //  ntheta    bins in theta, eta
800   //  themin    minimum angle in theta (degrees)
801   //  themax    maximum angle in theta (degrees)
802   //  nphi      bins in phi
803   //  phimin    minimum angle in phi (degrees)
804   //  phimax    maximum angle in phi (degrees)
805   //  rmin      minimum radius
806   //  rmax      maximum radius
807   //  
808   //
809   //  The number of events generated = ntheta*nphi
810   //  run input parameters in macro setup (default="Config.C")
811   //
812   //  Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
813   //Begin_Html
814   /*
815     <img src="picts/AliRunLego1.gif">
816   */
817   //End_Html
818   //Begin_Html
819   /*
820     <img src="picts/AliRunLego2.gif">
821   */
822   //End_Html
823   //Begin_Html
824   /*
825     <img src="picts/AliRunLego3.gif">
826   */
827   //End_Html
828   //
829
830 // run the generation and simulation
831
832   AliCodeTimerAuto("",0)
833
834   // initialize CDB storage and run number from external environment
835   // (either CDB manager or AliSimulation setters)
836   InitCDB();
837   InitRunNumber();
838   SetCDBLock();
839   
840   if (!gAlice) {
841     AliError("no gAlice object. Restart aliroot and try again.");
842     return kFALSE;
843   }
844   if (gAlice->Modules()->GetEntries() > 0) {
845     AliError("gAlice was already run. Restart aliroot and try again.");
846     return kFALSE;
847   }
848
849   AliInfo(Form("initializing gAlice with config file %s",
850           fConfigFileName.Data()));
851
852   // Number of events 
853     if (nev == -1) nev  = nc1 * nc2;
854     
855   // check if initialisation has been done
856   // If runloader has been initialized, set the number of events per file to nc1 * nc2
857     
858   // Set new generator
859   if (!gener) gener  = new AliLegoGenerator();
860   //
861   // Configure Generator
862
863   gener->SetRadiusRange(rmin, rmax);
864   gener->SetZMax(zmax);
865   gener->SetCoor1Range(nc1, c1min, c1max);
866   gener->SetCoor2Range(nc2, c2min, c2max);
867   
868   
869   //Create Lego object  
870   fLego = new AliLego("lego",gener);
871
872   //__________________________________________________________________________
873
874   gAlice->Announce();
875
876   gROOT->LoadMacro(setup);
877   gInterpreter->ProcessLine(gAlice->GetConfigFunction());
878
879   if(AliCDBManager::Instance()->GetRun() >= 0) { 
880     SetRunNumber(AliCDBManager::Instance()->GetRun());
881   } else {
882     AliWarning("Run number not initialized!!");
883   }
884
885   AliRunLoader::Instance()->CdGAFile();
886   
887   AliPDG::AddParticlesToPdgDataBase();  
888   
889   gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
890
891   gAlice->GetMCApp()->Init();
892   
893   
894   //Must be here because some MCs (G4) adds detectors here and not in Config.C
895   gAlice->InitLoaders();
896   AliRunLoader::Instance()->MakeTree("E");
897   
898   //
899   // Save stuff at the beginning of the file to avoid file corruption
900   AliRunLoader::Instance()->CdGAFile();
901   gAlice->Write();
902
903   //Save current generator
904   AliGenerator *gen=gAlice->GetMCApp()->Generator();
905   gAlice->GetMCApp()->ResetGenerator(gener);
906   //Prepare MC for Lego Run
907   gMC->InitLego();
908   
909   //Run Lego Object
910   
911   
912   AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
913   gMC->ProcessRun(nev);
914   
915   // End of this run, close files
916   FinishRun();
917   // Restore current generator
918   gAlice->GetMCApp()->ResetGenerator(gen);
919   // Delete Lego Object
920   delete fLego;
921
922   return kTRUE;
923 }
924
925 //_____________________________________________________________________________
926 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
927 {
928   // run the trigger
929
930   AliCodeTimerAuto("",0)
931
932   // initialize CDB storage from external environment
933   // (either CDB manager or AliSimulation setters),
934   // if not already done in RunSimulation()
935   InitCDB();
936   
937   // Set run number in CDBManager from data 
938   // From this point on the run number must be always loaded from data!
939   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
940   
941   // Set CDB lock: from now on it is forbidden to reset the run number
942   // or the default storage or to activate any further storage!
943   SetCDBLock();
944    
945    AliRunLoader* runLoader = LoadRun("READ");
946    if (!runLoader) return kFALSE;
947    TString trconfiguration = config;
948
949    if (trconfiguration.IsNull()) {
950      if(!fTriggerConfig.IsNull()) {
951        trconfiguration = fTriggerConfig;
952      }
953      else
954        AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
955    }
956
957    runLoader->MakeTree( "GG" );
958    AliCentralTrigger* aCTP = runLoader->GetTrigger();
959    // Load Configuration
960    if (!aCTP->LoadConfiguration( trconfiguration ))
961      return kFALSE;
962
963    // digits -> trigger
964    if( !aCTP->RunTrigger( runLoader , detectors ) ) {
965       if (fStopOnError) {
966         //  delete aCTP;
967         return kFALSE;
968       }
969    }
970
971    delete runLoader;
972
973    return kTRUE;
974 }
975
976 //_____________________________________________________________________________
977 Bool_t AliSimulation::WriteTriggerRawData()
978 {
979   // Writes the CTP (trigger) DDL raw data
980   // Details of the format are given in the
981   // trigger TDR - pages 134 and 135.
982   AliCTPRawData writer;
983   writer.RawData();
984
985   return kTRUE;
986 }
987
988 //_____________________________________________________________________________
989 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
990 {
991 // run the generation and simulation
992
993   AliCodeTimerAuto("",0)
994
995   // initialize CDB storage and run number from external environment
996   // (either CDB manager or AliSimulation setters)
997   AliSysInfo::AddStamp("RunSimulation_Begin");
998   InitCDB();
999   AliSysInfo::AddStamp("RunSimulation_InitCDB");
1000   InitRunNumber();
1001
1002   SetCDBLock();
1003   AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
1004   
1005   if (!gAlice) {
1006     AliError("no gAlice object. Restart aliroot and try again.");
1007     return kFALSE;
1008   }
1009   if (gAlice->Modules()->GetEntries() > 0) {
1010     AliError("gAlice was already run. Restart aliroot and try again.");
1011     return kFALSE;
1012   }
1013   
1014   // Setup monitoring if requested
1015   gAlice->GetMCApp()->SetUseMonitoring(fUseMonitoring);
1016
1017   AliInfo(Form("initializing gAlice with config file %s",
1018           fConfigFileName.Data()));
1019
1020   //
1021   // Initialize ALICE Simulation run
1022   //
1023   gAlice->Announce();
1024
1025   //
1026   // If requested set the mag. field from the GRP entry.
1027   // After this the field is loccked and cannot be changed by Config.C
1028   if (fUseMagFieldFromGRP) {
1029       AliGRPManager grpM;
1030       grpM.ReadGRPEntry();
1031       grpM.SetMagField();
1032       AliInfo("Field is locked now. It cannot be changed in Config.C");
1033   }
1034 //
1035 // Execute Config.C
1036   TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
1037   gROOT->LoadMacro(fConfigFileName.Data());
1038   Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
1039   if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
1040     AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
1041   }
1042   AliSysInfo::AddStamp("RunSimulation_Config");
1043
1044 //
1045 // If requested obtain the vertex position and vertex sigma_z from the CDB
1046 // This overwrites the settings from the Config.C  
1047   if (fUseVertexFromCDB) {
1048       Double_t vtxPos[3] = {0., 0., 0.}; 
1049       Double_t vtxSig[3] = {0., 0., 0.};
1050       AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1051       if (entry) {
1052           AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1053           if (vertex) {
1054               if(vertex->GetXRes()>2.8) { // > pipe radius --> it's a dummy object, don't use it 
1055                   entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1056                   if (entry) vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1057               }
1058           }
1059           if (vertex) {
1060               vertex->GetXYZ(vtxPos);
1061               vertex->GetSigmaXYZ(vtxSig);
1062               AliInfo("Overwriting Config.C vertex settings !");
1063               AliInfo(Form("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
1064                            vtxPos[0], vtxPos[1], vtxPos[2], vtxSig[2]));
1065               
1066               AliGenerator *gen = gAlice->GetMCApp()->Generator();
1067               gen->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]);   // vertex position
1068               gen->SetSigmaZ(vtxSig[2]);
1069           }
1070       }
1071   }
1072
1073   // If requested we take the SOR and EOR time-stamps from the GRP and use them
1074   // in order to generate the event time-stamps
1075   if (fUseTimeStampFromCDB) {
1076     AliGRPManager grpM;
1077     grpM.ReadGRPEntry();
1078     const AliGRPObject *grpObj = grpM.GetGRPData();
1079     if (!grpObj || (grpObj->GetTimeEnd() <= grpObj->GetTimeStart())) {
1080       AliError("Missing GRP or bad SOR/EOR time-stamps! Switching off the time-stamp generation from GRP!");
1081       fTimeStart = fTimeEnd = 0;
1082       fUseTimeStampFromCDB = kFALSE;
1083     }
1084     else {
1085       fTimeStart = grpObj->GetTimeStart();
1086       fTimeEnd = grpObj->GetTimeEnd();
1087     }
1088   }
1089   
1090   if(AliCDBManager::Instance()->GetRun() >= 0) { 
1091     AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
1092     AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
1093   } else {
1094         AliWarning("Run number not initialized!!");
1095   }
1096   
1097    AliRunLoader::Instance()->CdGAFile();
1098    
1099
1100    AliPDG::AddParticlesToPdgDataBase();  
1101
1102    gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1103    AliSysInfo::AddStamp("RunSimulation_GetField");
1104    gAlice->GetMCApp()->Init();
1105    AliSysInfo::AddStamp("RunSimulation_InitMCApp");
1106
1107    //Must be here because some MCs (G4) adds detectors here and not in Config.C
1108    gAlice->InitLoaders();
1109    AliRunLoader::Instance()->MakeTree("E");
1110    AliRunLoader::Instance()->LoadKinematics("RECREATE");
1111    AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1112    AliRunLoader::Instance()->LoadHits("all","RECREATE");
1113    //
1114    // Save stuff at the beginning of the file to avoid file corruption
1115    AliRunLoader::Instance()->CdGAFile();
1116    gAlice->Write();
1117    gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1118    AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1119   //___________________________________________________________________________________________
1120   
1121   AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
1122
1123   // Set run number in CDBManager
1124   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
1125
1126   AliRunLoader* runLoader = AliRunLoader::Instance();
1127   if (!runLoader) {
1128              AliError(Form("gAlice has no run loader object. "
1129                              "Check your config file: %s", fConfigFileName.Data()));
1130              return kFALSE;
1131   }
1132   SetGAliceFile(runLoader->GetFileName());
1133       
1134   // Misalign geometry
1135 #if ROOT_VERSION_CODE < 331527
1136   AliGeomManager::SetGeometry(gGeoManager);
1137   
1138   // Check that the consistency of symbolic names for the activated subdetectors
1139   // in the geometry loaded by AliGeomManager
1140   TString detsToBeChecked = "";
1141   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1142   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1143     AliModule* det = (AliModule*) detArray->At(iDet);
1144     if (!det || !det->IsActive()) continue;
1145     detsToBeChecked += det->GetName();
1146     detsToBeChecked += " ";
1147   } // end loop over detectors
1148   if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
1149     AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
1150   MisalignGeometry(runLoader);
1151   AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
1152 #endif
1153
1154 //   AliRunLoader* runLoader = AliRunLoader::Instance();
1155 //   if (!runLoader) {
1156 //     AliError(Form("gAlice has no run loader object. "
1157 //                   "Check your config file: %s", fConfigFileName.Data()));
1158 //     return kFALSE;
1159 //   }
1160 //   SetGAliceFile(runLoader->GetFileName());
1161
1162   if (!gAlice->GetMCApp()->Generator()) {
1163     AliError(Form("gAlice has no generator object. "
1164                   "Check your config file: %s", fConfigFileName.Data()));
1165     return kFALSE;
1166   }
1167
1168   // Write GRP entry corresponding to the setting found in Cofig.C
1169   if (fWriteGRPEntry)
1170     WriteGRPEntry();
1171   AliSysInfo::AddStamp("RunSimulation_WriteGRP");
1172
1173   if (nEvents <= 0) nEvents = fNEvents;
1174
1175   // get vertex from background file in case of merging
1176   if (fUseBkgrdVertex &&
1177       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
1178     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
1179     const char* fileName = ((TObjString*)
1180                             (fBkgrdFileNames->At(0)))->GetName();
1181     AliInfo(Form("The vertex will be taken from the background "
1182                  "file %s with nSignalPerBackground = %d", 
1183                  fileName, signalPerBkgrd));
1184     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
1185     gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
1186   }
1187
1188   if (!fRunSimulation) {
1189     gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
1190   }
1191
1192   // set the number of events per file for given detectors and data types
1193   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
1194     if (!fEventsPerFile[i]) continue;
1195     const char* detName = fEventsPerFile[i]->GetName();
1196     const char* typeName = fEventsPerFile[i]->GetTitle();
1197     TString loaderName(detName);
1198     loaderName += "Loader";
1199     AliLoader* loader = runLoader->GetLoader(loaderName);
1200     if (!loader) {
1201       AliError(Form("RunSimulation no loader for %s found\n Number of events per file not set for %s %s", 
1202                     detName, typeName, detName));
1203       continue;
1204     }
1205     AliDataLoader* dataLoader = 
1206       loader->GetDataLoader(typeName);
1207     if (!dataLoader) {
1208       AliError(Form("no data loader for %s found\n"
1209                     "Number of events per file not set for %s %s", 
1210                     typeName, detName, typeName));
1211       continue;
1212     }
1213     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
1214     AliDebug(1, Form("number of events per file set to %d for %s %s",
1215                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
1216   }
1217
1218   AliInfo("running gAlice");
1219   AliSysInfo::AddStamp("Start_ProcessRun");
1220
1221   // Create the Root Tree with one branch per detector
1222   //Hits moved to begin event -> now we are crating separate tree for each event
1223   gMC->ProcessRun(nEvents);
1224
1225   // End of this run, close files
1226   if(nEvents>0) FinishRun();
1227
1228   AliSysInfo::AddStamp("Stop_ProcessRun");
1229   delete runLoader;
1230
1231   return kTRUE;
1232 }
1233
1234 //_____________________________________________________________________________
1235 Bool_t AliSimulation::RunGeneratorOnly()
1236 {
1237   // Execute Config.C
1238   TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
1239   gROOT->LoadMacro(fConfigFileName.Data());
1240   Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
1241   if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
1242     AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
1243   }
1244
1245   // Setup the runloader and generator, check if everything is OK
1246   AliRunLoader* runLoader = AliRunLoader::Instance();
1247   AliGenerator* generator = gAlice->GetMCApp()->Generator();
1248   if (!runLoader) {
1249     AliError(Form("gAlice has no run loader object. "
1250                   "Check your config file: %s", fConfigFileName.Data()));
1251     return kFALSE;
1252   }
1253   if (!generator) {
1254     AliError(Form("gAlice has no generator object. "
1255                   "Check your config file: %s", fConfigFileName.Data()));
1256     return kFALSE;
1257   }
1258
1259   runLoader->LoadKinematics("RECREATE");
1260   runLoader->MakeTree("E");
1261
1262   // Create stack and header
1263   runLoader->MakeStack();
1264   AliStack*  stack      = runLoader->Stack();
1265   AliHeader* header     = runLoader->GetHeader();
1266
1267   // Intialize generator
1268   generator->Init();
1269   generator->SetStack(stack);
1270
1271   // Run main generator loop
1272
1273   for (Int_t iev=0; iev<fNEvents; iev++)
1274   {
1275     // Initialize event
1276     header->Reset(0,iev);
1277     runLoader->SetEventNumber(iev);
1278     stack->Reset();
1279     runLoader->MakeTree("K");
1280
1281     // Generate event
1282     generator->Generate();
1283
1284     // Finish event
1285     header->SetNprimary(stack->GetNprimary());
1286     header->SetNtrack(stack->GetNtrack());
1287     stack->FinishEvent();
1288     header->SetStack(stack);
1289     runLoader->TreeE()->Fill();
1290     runLoader->WriteKinematics("OVERWRITE");
1291   }
1292
1293   // Finalize
1294   generator->FinishRun();
1295   // Write file
1296   runLoader->WriteHeader("OVERWRITE");
1297   generator->Write();
1298   runLoader->Write();
1299
1300   return kTRUE;
1301 }
1302
1303 //_____________________________________________________________________________
1304 Bool_t AliSimulation::RunSDigitization(const char* detectors)
1305 {
1306 // run the digitization and produce summable digits
1307   static Int_t eventNr=0;
1308   AliCodeTimerAuto("",0) ;
1309
1310   // initialize CDB storage, run number, set CDB lock
1311   InitCDB();
1312   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1313   SetCDBLock();
1314   
1315   AliRunLoader* runLoader = LoadRun();
1316   if (!runLoader) return kFALSE;
1317
1318   TString detStr = detectors;
1319   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1320   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1321     AliModule* det = (AliModule*) detArray->At(iDet);
1322     if (!det || !det->IsActive()) continue;
1323     if (IsSelected(det->GetName(), detStr)) {
1324       AliInfo(Form("creating summable digits for %s", det->GetName()));
1325       AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
1326       det->Hits2SDigits();
1327       AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
1328       AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
1329     }
1330   }
1331
1332   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1333     AliError(Form("the following detectors were not found: %s",
1334                   detStr.Data()));
1335     if (fStopOnError) return kFALSE;
1336   }
1337   eventNr++;
1338   delete runLoader;
1339
1340   return kTRUE;
1341 }
1342
1343
1344 //_____________________________________________________________________________
1345 Bool_t AliSimulation::RunDigitization(const char* detectors, 
1346                                       const char* excludeDetectors)
1347 {
1348 // run the digitization and produce digits from sdigits
1349
1350   AliCodeTimerAuto("",0)
1351
1352   // initialize CDB storage, run number, set CDB lock
1353   InitCDB();
1354   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1355   SetCDBLock();
1356   
1357   delete AliRunLoader::Instance();
1358   delete gAlice;
1359   gAlice = NULL;
1360
1361   Int_t nStreams = 1;
1362   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1363   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1364   AliDigitizationInput digInp(nStreams, signalPerBkgrd);
1365   // digInp.SetEmbeddingFlag(fEmbeddingFlag);
1366   digInp.SetRegionOfInterest(fRegionOfInterest);
1367   digInp.SetInputStream(0, fGAliceFileName.Data());
1368   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1369     const char* fileName = ((TObjString*)(fBkgrdFileNames->At(iStream-1)))->GetName();
1370     digInp.SetInputStream(iStream, fileName);
1371   }
1372   TObjArray detArr;
1373   detArr.SetOwner(kTRUE);
1374   TString detStr = detectors;
1375   TString detExcl = excludeDetectors;
1376   if (!static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice()) {
1377     AliError("Error occured while getting gAlice from Input 0");
1378     return kFALSE;
1379   }
1380   AliRunLoader* runLoader = AliRunLoader::GetRunLoader(digInp.GetInputStream(0)->GetFolderName());
1381   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1382   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1383     AliModule* det = (AliModule*) detArray->At(iDet);
1384     if (!det || !det->IsActive()) continue;
1385     if (!IsSelected(det->GetName(), detStr) || IsSelected(det->GetName(), detExcl)) continue;
1386     AliDigitizer* digitizer = det->CreateDigitizer(&digInp);
1387     if (!digitizer || !digitizer->Init()) {
1388       AliError(Form("no digitizer for %s", det->GetName()));
1389       if (fStopOnError) return kFALSE;
1390       else continue;
1391     }
1392     detArr.AddLast(digitizer);    
1393     AliInfo(Form("Created digitizer from SDigits -> Digits for %s", det->GetName()));    
1394   }
1395   //
1396   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1397     AliError(Form("the following detectors were not found: %s", detStr.Data()));
1398     if (fStopOnError) return kFALSE;
1399   }
1400   //
1401   Int_t ndigs = detArr.GetEntriesFast();
1402   Int_t eventsCreated = 0;
1403   AliRunLoader* outRl =  digInp.GetOutRunLoader();
1404   while ((eventsCreated++ < fNEvents) || (fNEvents < 0)) {
1405     if (!digInp.ConnectInputTrees()) break;
1406     digInp.InitEvent(); //this must be after call of Connect Input tress.
1407     if (outRl) outRl->SetEventNumber(eventsCreated-1);
1408     static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice(); // use gAlice of the first input stream
1409     for (int id=0;id<ndigs;id++) ((AliDigitizer*)detArr[id])->Digitize("");
1410     digInp.FinishEvent();
1411   };
1412   digInp.FinishGlobal();
1413   //
1414   return kTRUE;
1415 }
1416
1417 //_____________________________________________________________________________
1418 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1419 {
1420 // run the digitization and produce digits from hits
1421
1422   AliCodeTimerAuto("",0)
1423
1424   // initialize CDB storage, run number, set CDB lock
1425   InitCDB();
1426   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1427   SetCDBLock();
1428   
1429   AliRunLoader* runLoader = LoadRun("READ");
1430   if (!runLoader) return kFALSE;
1431
1432   TString detStr = detectors;
1433   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1434   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1435     AliModule* det = (AliModule*) detArray->At(iDet);
1436     if (!det || !det->IsActive()) continue;
1437     if (IsSelected(det->GetName(), detStr)) {
1438       AliInfo(Form("creating digits from hits for %s", det->GetName()));
1439       det->Hits2Digits();
1440     }
1441   }
1442
1443   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1444     AliError(Form("the following detectors were not found: %s", 
1445                   detStr.Data()));
1446     if (fStopOnError) return kFALSE;
1447   }
1448
1449   return kTRUE;
1450 }
1451
1452 //_____________________________________________________________________________
1453 Bool_t AliSimulation::WriteRawData(const char* detectors, 
1454                                    const char* fileName,
1455                                    Bool_t deleteIntermediateFiles,
1456                                    Bool_t selrawdata)
1457 {
1458 // convert the digits to raw data
1459 // First DDL raw data files for the given detectors are created.
1460 // If a file name is given, the DDL files are then converted to a DATE file.
1461 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
1462 // afterwards.
1463 // If the file name has the extension ".root", the DATE file is converted
1464 // to a root file.
1465 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1466 // 'selrawdata' flag can be used to enable writing of detectors raw data
1467 // accoring to the trigger cluster.
1468
1469   AliCodeTimerAuto("",0)
1470   AliSysInfo::AddStamp("WriteRawData_Start");
1471   
1472   TString detStr = detectors;
1473   if (!WriteRawFiles(detStr.Data())) {
1474     if (fStopOnError) return kFALSE;
1475   }
1476   AliSysInfo::AddStamp("WriteRawFiles");
1477
1478   // run HLT simulation on simulated DDL raw files
1479   // and produce HLT ddl raw files to be included in date/root file
1480   // bugfix 2009-06-26: the decision whether to write HLT raw data
1481   // is taken in RunHLT. Here HLT always needs to be run in order to
1482   // create HLT digits, unless its switched off. This is due to the
1483   // special placement of the HLT between the generation of DDL files
1484   // and conversion to DATE/Root file.
1485   detStr.ReplaceAll("HLT", "");
1486   if (!fRunHLT.IsNull()) {
1487     if (!RunHLT()) {
1488       if (fStopOnError) return kFALSE;
1489     }
1490   }
1491   AliSysInfo::AddStamp("WriteRawData_RunHLT");
1492
1493   TString dateFileName(fileName);
1494   if (!dateFileName.IsNull()) {
1495     Bool_t rootOutput = dateFileName.EndsWith(".root");
1496     if (rootOutput) dateFileName += ".date";
1497     TString selDateFileName;
1498     if (selrawdata) {
1499       selDateFileName = "selected.";
1500       selDateFileName+= dateFileName;
1501     }
1502     if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1503       if (fStopOnError) return kFALSE;
1504     }
1505     AliSysInfo::AddStamp("ConvertRawFilesToDate");
1506     if (deleteIntermediateFiles) {
1507       AliRunLoader* runLoader = LoadRun("READ");
1508       if (runLoader) for (Int_t iEvent = 0; 
1509                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1510         char command[256];
1511         snprintf(command, 256, "rm -r raw%d", iEvent);
1512         gSystem->Exec(command);
1513       }
1514       delete runLoader;
1515     }
1516
1517     if (rootOutput) {
1518       if (!ConvertDateToRoot(dateFileName, fileName)) {
1519         if (fStopOnError) return kFALSE;
1520       }
1521       AliSysInfo::AddStamp("ConvertDateToRoot");
1522       if (deleteIntermediateFiles) {
1523         gSystem->Unlink(dateFileName);
1524       }
1525       if (selrawdata) {
1526         TString selFileName = "selected.";
1527         selFileName        += fileName;
1528         if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1529           if (fStopOnError) return kFALSE;
1530         }
1531         if (deleteIntermediateFiles) {
1532           gSystem->Unlink(selDateFileName);
1533         }
1534       }
1535     }
1536   }
1537
1538   return kTRUE;
1539 }
1540
1541 //_____________________________________________________________________________
1542 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1543 {
1544 // convert the digits to raw data DDL files
1545
1546   AliCodeTimerAuto("",0)
1547   
1548   AliRunLoader* runLoader = LoadRun("READ");
1549   if (!runLoader) return kFALSE;
1550
1551   // write raw data to DDL files
1552   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1553     AliInfo(Form("processing event %d", iEvent));
1554     runLoader->GetEvent(iEvent);
1555     TString baseDir = gSystem->WorkingDirectory();
1556     char dirName[256];
1557     snprintf(dirName, 256, "raw%d", iEvent);
1558     gSystem->MakeDirectory(dirName);
1559     if (!gSystem->ChangeDirectory(dirName)) {
1560       AliError(Form("couldn't change to directory %s", dirName));
1561       if (fStopOnError) return kFALSE; else continue;
1562     }
1563
1564     ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1565     runNbFile.close();
1566
1567     TString detStr = detectors;
1568     if (IsSelected("HLT", detStr)) {
1569       // Do nothing. "HLT" will be removed from detStr and HLT raw
1570       // data files are generated in RunHLT.
1571     }
1572
1573     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1574     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1575       AliModule* det = (AliModule*) detArray->At(iDet);
1576       if (!det || !det->IsActive()) continue;
1577       if (IsSelected(det->GetName(), detStr)) {
1578         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1579         det->Digits2Raw();
1580       }
1581     }
1582
1583     if (!WriteTriggerRawData())
1584       if (fStopOnError) return kFALSE;
1585
1586     gSystem->ChangeDirectory(baseDir);
1587     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1588       AliError(Form("the following detectors were not found: %s", 
1589                     detStr.Data()));
1590       if (fStopOnError) return kFALSE;
1591     }
1592   }
1593
1594   delete runLoader;
1595   
1596   return kTRUE;
1597 }
1598
1599 //_____________________________________________________________________________
1600 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1601                                             const char* selDateFileName)
1602 {
1603 // convert raw data DDL files to a DATE file with the program "dateStream"
1604 // The second argument is not empty when the user decides to write
1605 // the detectors raw data according to the trigger cluster.
1606
1607   AliCodeTimerAuto("",0)
1608   
1609   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1610   if (!path) {
1611     AliError("the program dateStream was not found");
1612     if (fStopOnError) return kFALSE;
1613   } else {
1614     delete[] path;
1615   }
1616
1617   AliRunLoader* runLoader = LoadRun("READ");
1618   if (!runLoader) return kFALSE;
1619
1620   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1621   Bool_t selrawdata = kFALSE;
1622   if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1623
1624   char command[256];
1625   // Note the option -s. It is used in order to avoid
1626   // the generation of SOR/EOR events.
1627   snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1628           dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1629   FILE* pipe = gSystem->OpenPipe(command, "w");
1630
1631   if (!pipe) {
1632     AliError(Form("Cannot execute command: %s",command));
1633     return kFALSE;
1634   }
1635
1636   Int_t selEvents = 0;
1637   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1638
1639     UInt_t detectorPattern = 0;
1640     runLoader->GetEvent(iEvent);
1641     if (!runLoader->LoadTrigger()) {
1642       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1643       detectorPattern = aCTP->GetClusterMask();
1644       // Check if the event was triggered by CTP
1645       if (selrawdata) {
1646         if (aCTP->GetClassMask()) selEvents++;
1647       }
1648     }
1649     else {
1650       AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1651       if (selrawdata) {
1652         AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1653         selrawdata = kFALSE;
1654       }
1655     }
1656
1657     fprintf(pipe, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
1658     Float_t ldc = 0;
1659     Int_t prevLDC = -1;
1660
1661     // loop over detectors and DDLs
1662     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1663       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1664
1665         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1666         Int_t ldcID = Int_t(ldc + 0.0001);
1667         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1668
1669         char rawFileName[256];
1670         snprintf(rawFileName, 256, "raw%d/%s", 
1671                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1672
1673         // check existence and size of raw data file
1674         FILE* file = fopen(rawFileName, "rb");
1675         if (!file) continue;
1676         fseek(file, 0, SEEK_END);
1677         unsigned long size = ftell(file);
1678         fclose(file);
1679         if (!size) continue;
1680
1681         if (ldcID != prevLDC) {
1682           fprintf(pipe, " LDC Id %d\n", ldcID);
1683           prevLDC = ldcID;
1684         }
1685         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1686       }
1687     }
1688   }
1689
1690   Int_t result = gSystem->ClosePipe(pipe);
1691
1692   if (!(selrawdata && selEvents > 0)) {
1693     delete runLoader;
1694     return (result == 0);
1695   }
1696
1697   AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1698   
1699   snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1700           selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1701   FILE* pipe2 = gSystem->OpenPipe(command, "w");
1702
1703   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1704
1705     // Get the trigger decision and cluster
1706     UInt_t detectorPattern = 0;
1707     TString detClust;
1708     runLoader->GetEvent(iEvent);
1709     if (!runLoader->LoadTrigger()) {
1710       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1711       if (aCTP->GetClassMask() == 0) continue;
1712       detectorPattern = aCTP->GetClusterMask();
1713       detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1714       AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1715     }
1716
1717     fprintf(pipe2, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
1718     Float_t ldc = 0;
1719     Int_t prevLDC = -1;
1720
1721     // loop over detectors and DDLs
1722     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1723       // Write only raw data from detectors that
1724       // are contained in the trigger cluster(s)
1725       if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1726
1727       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1728
1729         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1730         Int_t ldcID = Int_t(ldc + 0.0001);
1731         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1732
1733         char rawFileName[256];
1734         snprintf(rawFileName, 256, "raw%d/%s", 
1735                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1736
1737         // check existence and size of raw data file
1738         FILE* file = fopen(rawFileName, "rb");
1739         if (!file) continue;
1740         fseek(file, 0, SEEK_END);
1741         unsigned long size = ftell(file);
1742         fclose(file);
1743         if (!size) continue;
1744
1745         if (ldcID != prevLDC) {
1746           fprintf(pipe2, " LDC Id %d\n", ldcID);
1747           prevLDC = ldcID;
1748         }
1749         fprintf(pipe2, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1750       }
1751     }
1752   }
1753
1754   Int_t result2 = gSystem->ClosePipe(pipe2);
1755
1756   delete runLoader;
1757   return ((result == 0) && (result2 == 0));
1758 }
1759
1760 //_____________________________________________________________________________
1761 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1762                                         const char* rootFileName)
1763 {
1764 // convert a DATE file to a root file with the program "alimdc"
1765
1766   // ALIMDC setup
1767   const Int_t kDBSize = 2000000000;
1768   const Int_t kTagDBSize = 1000000000;
1769   const Bool_t kFilter = kFALSE;
1770   const Int_t kCompression = 1;
1771
1772   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1773   if (!path) {
1774     AliError("the program alimdc was not found");
1775     if (fStopOnError) return kFALSE;
1776   } else {
1777     delete[] path;
1778   }
1779
1780   AliInfo(Form("converting DATE file %s to root file %s", 
1781                dateFileName, rootFileName));
1782
1783   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1784   const char* tagDBFS    = "/tmp/mdc1/tags";
1785
1786   // User defined file system locations
1787   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1788     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1789   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1790     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1791   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1792     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1793
1794   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1795   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1796   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1797
1798   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1799   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1800   gSystem->Exec(Form("mkdir %s",tagDBFS));
1801
1802   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1803                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1804   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1805
1806   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1807   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1808   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1809
1810   return (result == 0);
1811 }
1812
1813
1814 //_____________________________________________________________________________
1815 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1816 {
1817 // delete existing run loaders, open a new one and load gAlice
1818
1819   delete AliRunLoader::Instance();
1820   AliRunLoader* runLoader = 
1821     AliRunLoader::Open(fGAliceFileName.Data(), 
1822                        AliConfig::GetDefaultEventFolderName(), mode);
1823   if (!runLoader) {
1824     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1825     return NULL;
1826   }
1827   runLoader->LoadgAlice();
1828   runLoader->LoadHeader();
1829   gAlice = runLoader->GetAliRun();
1830   if (!gAlice) {
1831     AliError(Form("no gAlice object found in file %s", 
1832                   fGAliceFileName.Data()));
1833     return NULL;
1834   }
1835   return runLoader;
1836 }
1837
1838 //_____________________________________________________________________________
1839 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1840 {
1841 // get or calculate the number of signal events per background event
1842
1843   if (!fBkgrdFileNames) return 1;
1844   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1845   if (nBkgrdFiles == 0) return 1;
1846
1847   // get the number of signal events
1848   if (nEvents <= 0) {
1849     AliRunLoader* runLoader = 
1850         AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1851     if (!runLoader) return 1;
1852     
1853     nEvents = runLoader->GetNumberOfEvents();
1854     delete runLoader;
1855   }
1856
1857   Int_t result = 0;
1858   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1859     // get the number of background events
1860     const char* fileName = ((TObjString*)
1861                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1862     AliRunLoader* runLoader =
1863       AliRunLoader::Open(fileName, "BKGRD");
1864     if (!runLoader) continue;
1865     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1866     delete runLoader;
1867   
1868     // get or calculate the number of signal per background events
1869     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1870     if (nSignalPerBkgrd <= 0) {
1871       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1872     } else if (result && (result != nSignalPerBkgrd)) {
1873       AliInfo(Form("the number of signal events per background event "
1874                    "will be changed from %d to %d for stream %d", 
1875                    nSignalPerBkgrd, result, iBkgrdFile+1));
1876       nSignalPerBkgrd = result;
1877     }
1878
1879     if (!result) result = nSignalPerBkgrd;
1880     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1881       AliWarning(Form("not enough background events (%d) for %d signal events "
1882                       "using %d signal per background events for stream %d",
1883                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1884     }
1885   }
1886
1887   return result;
1888 }
1889
1890 //_____________________________________________________________________________
1891 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1892 {
1893 // check whether detName is contained in detectors
1894 // if yes, it is removed from detectors
1895
1896   // check if all detectors are selected
1897   if ((detectors.CompareTo("ALL") == 0) ||
1898       detectors.BeginsWith("ALL ") ||
1899       detectors.EndsWith(" ALL") ||
1900       detectors.Contains(" ALL ")) {
1901     detectors = "ALL";
1902     return kTRUE;
1903   }
1904
1905   // search for the given detector
1906   Bool_t result = kFALSE;
1907   if ((detectors.CompareTo(detName) == 0) ||
1908       detectors.BeginsWith(detName+" ") ||
1909       detectors.EndsWith(" "+detName) ||
1910       detectors.Contains(" "+detName+" ")) {
1911     detectors.ReplaceAll(detName, "");
1912     result = kTRUE;
1913   }
1914
1915   // clean up the detectors string
1916   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1917   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1918   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1919
1920   return result;
1921 }
1922
1923 //_____________________________________________________________________________
1924 Int_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName, Int_t N, Int_t nSkip)
1925 {
1926 //
1927 // Steering routine  to convert raw data in directory rawDirectory/ to fake SDigits. 
1928 // These can be used for embedding of MC tracks into RAW data using the standard 
1929 // merging procedure.
1930 //
1931 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1932 //
1933   if (!gAlice) {
1934     AliError("no gAlice object. Restart aliroot and try again.");
1935     return kFALSE;
1936   }
1937   if (gAlice->Modules()->GetEntries() > 0) {
1938     AliError("gAlice was already run. Restart aliroot and try again.");
1939     return kFALSE;
1940   }
1941   
1942   AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1943   
1944   gAlice->Announce();
1945   
1946   gROOT->LoadMacro(fConfigFileName.Data());
1947   gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1948   
1949   if(AliCDBManager::Instance()->GetRun() >= 0) { 
1950         SetRunNumber(AliCDBManager::Instance()->GetRun());
1951   } else {
1952         AliWarning("Run number not initialized!!");
1953   }
1954   
1955    AliRunLoader::Instance()->CdGAFile();
1956     
1957    AliPDG::AddParticlesToPdgDataBase();  
1958
1959    gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1960    
1961    gAlice->GetMCApp()->Init();
1962
1963    //Must be here because some MCs (G4) adds detectors here and not in Config.C
1964    gAlice->InitLoaders();
1965    AliRunLoader::Instance()->MakeTree("E");
1966    AliRunLoader::Instance()->LoadKinematics("RECREATE");
1967    AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1968    AliRunLoader::Instance()->LoadHits("all","RECREATE");
1969
1970    //
1971    // Save stuff at the beginning of the file to avoid file corruption
1972    AliRunLoader::Instance()->CdGAFile();
1973    gAlice->Write();
1974 //
1975 //  Initialize CDB     
1976     InitCDB();
1977     //AliCDBManager* man = AliCDBManager::Instance();
1978     //man->SetRun(0); // Should this come from rawdata header ?
1979     
1980     Int_t iDet;
1981     //
1982     // Get the runloader
1983     AliRunLoader* runLoader = AliRunLoader::Instance();
1984     //
1985     // Open esd file if available
1986     TFile* esdFile = 0;
1987     TTree* treeESD = 0;
1988     AliESDEvent* esd = 0;
1989     if (esdFileName && (strlen(esdFileName)>0)) {
1990       esdFile = TFile::Open(esdFileName);
1991       if (esdFile) {
1992         esd = new AliESDEvent();
1993         esdFile->GetObject("esdTree", treeESD);
1994                   if (treeESD) {
1995                           esd->ReadFromTree(treeESD);
1996                           if (nSkip>0) {
1997                                   AliInfo(Form("Asking to skip first %d ESDs events",nSkip));
1998                           } else {
1999                                   nSkip=0;
2000                           }
2001                   }
2002       }
2003     }
2004
2005     //
2006     // Create the RawReader
2007     TString fileName(rawDirectory);
2008     AliRawReader* rawReader = AliRawReader::Create(fileName.Data());
2009     if (!rawReader) return (kFALSE);
2010     
2011 //     if (!fEquipIdMap.IsNull() && fRawReader)
2012 //       fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
2013     //
2014     // Get list of detectors
2015     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
2016     //
2017     // Get Header
2018     AliHeader* header = runLoader->GetHeader();
2019     // Event loop
2020     Int_t nev = 0;
2021     while(kTRUE) {
2022         if (!(rawReader->NextEvent())) break;
2023         runLoader->SetEventNumber(nev);
2024         runLoader->GetHeader()->Reset(rawReader->GetRunNumber(), 
2025                                       nev, nev);
2026         runLoader->GetEvent(nev);
2027         AliInfo(Form("We are at event %d",nev));
2028         //
2029         // Detector loop
2030         TString detStr = fMakeSDigits;
2031         for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
2032             AliModule* det = (AliModule*) detArray->At(iDet);
2033             if (!det || !det->IsActive()) continue;
2034             if (IsSelected(det->GetName(), detStr)) {
2035               AliInfo(Form("Calling Raw2SDigits for %s", det->GetName()));
2036               det->Raw2SDigits(rawReader);
2037               rawReader->Reset();
2038             }
2039         } // detectors
2040         
2041
2042         //
2043         //  If ESD information available obtain reconstructed vertex and store in header.
2044         if (treeESD) {
2045                 AliInfo(Form("Selected event %d correspond to event %d in raw and to %d in esd",nev,rawReader->GetEventIndex(),nSkip+rawReader->GetEventIndex()));
2046             treeESD->GetEvent(nSkip+rawReader->GetEventIndex());
2047             const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
2048             Double_t position[3];
2049             esdVertex->GetXYZ(position);
2050             AliGenEventHeader* mcHeader = new  AliGenEventHeader("ESD");
2051             TArrayF mcV;
2052             mcV.Set(3);
2053             for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
2054             mcHeader->SetPrimaryVertex(mcV);
2055             header->Reset(0,nev);
2056             header->SetGenEventHeader(mcHeader);
2057             AliInfo(Form("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]));
2058         }
2059 //
2060 //      Finish the event
2061         runLoader->TreeE()->Fill();
2062         AliInfo(Form("Finished event %d",nev));
2063         nev++;
2064         if (N>0&&nev>=N)
2065           break;
2066     } // events
2067  
2068     delete rawReader;
2069 //
2070 //  Finish the run 
2071     runLoader->CdGAFile();
2072     runLoader->WriteHeader("OVERWRITE");
2073     runLoader->WriteRunLoader();
2074
2075     return nev;
2076 }
2077
2078 //_____________________________________________________________________________
2079 void AliSimulation::FinishRun()
2080 {
2081   //
2082   // Called at the end of the run.
2083   //
2084
2085   if(IsLegoRun()) 
2086    {
2087     AliDebug(1, "Finish Lego");
2088     AliRunLoader::Instance()->CdGAFile();
2089     fLego->FinishRun();
2090    }
2091   
2092   // Clean detector information
2093   TIter next(gAlice->Modules());
2094   AliModule *detector;
2095   while((detector = dynamic_cast<AliModule*>(next()))) {
2096     AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
2097     detector->FinishRun();
2098   }
2099   
2100   AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
2101   AliRunLoader::Instance()->WriteHeader("OVERWRITE");
2102
2103   // Write AliRun info and all detectors parameters
2104   AliRunLoader::Instance()->CdGAFile();
2105   gAlice->Write(0,TObject::kOverwrite);//write AliRun
2106   AliRunLoader::Instance()->Write(0,TObject::kOverwrite);//write RunLoader itself
2107   
2108   if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();  
2109   AliRunLoader::Instance()->Synchronize();
2110 }
2111
2112 //_____________________________________________________________________________
2113 Int_t AliSimulation::GetDetIndex(const char* detector)
2114 {
2115   // return the detector index corresponding to detector
2116   Int_t index = -1 ; 
2117   for (index = 0; index < fgkNDetectors ; index++) {
2118     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2119           break ; 
2120   }     
2121   return index ; 
2122 }
2123
2124 //_____________________________________________________________________________
2125 Bool_t AliSimulation::CreateHLT()
2126 {
2127   // Init the HLT simulation.
2128   // The function  loads the library and creates the instance of AliHLTSimulation.
2129   // the main reason for the decoupled creation is to set the transient OCDB
2130   // objects before the OCDB is locked
2131
2132   // load the library dynamically
2133   gSystem->Load(ALIHLTSIMULATION_LIBRARY);
2134
2135   // check for the library version
2136   AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
2137   if (!fctVersion) {
2138     AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
2139     return kFALSE;
2140   }
2141   if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
2142     AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
2143   }
2144
2145   // print compile info
2146   typedef void (*CompileInfo)( const char*& date, const char*& time);
2147   CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
2148   if (fctInfo) {
2149     const char* date="";
2150     const char* time="";
2151     (*fctInfo)(date, time);
2152     if (!date) date="unknown";
2153     if (!time) time="unknown";
2154     AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
2155   } else {
2156     AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
2157   }
2158
2159   // create instance of the HLT simulation
2160   AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
2161   if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
2162     AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
2163     return kFALSE;    
2164   }
2165
2166   TString specObjects;
2167   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
2168     if (specObjects.Length()>0) specObjects+=" ";
2169     specObjects+=fSpecCDBUri[i]->GetName();
2170   }
2171
2172   AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
2173   if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
2174     AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
2175   }
2176
2177   return kTRUE;
2178 }
2179
2180 //_____________________________________________________________________________
2181 Bool_t AliSimulation::RunHLT()
2182 {
2183   // Run the HLT simulation
2184   // HLT simulation is implemented in HLT/sim/AliHLTSimulation
2185   // Disabled if fRunHLT is empty, default vaule is "default".
2186   // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
2187   // The default simulation depends on the HLT component libraries and their
2188   // corresponding agents which define components and chains to run. See
2189   // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
2190   // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
2191   //
2192   // The libraries to be loaded can be specified as an option.
2193   // <pre>
2194   // AliSimulation sim;
2195   // sim.SetRunHLT("libAliHLTSample.so");
2196   // </pre>
2197   // will only load <tt>libAliHLTSample.so</tt>
2198
2199   // Other available options:
2200   // \li loglevel=<i>level</i> <br>
2201   //     logging level for this processing
2202   // \li alilog=off
2203   //     disable redirection of log messages to AliLog class
2204   // \li config=<i>macro</i>
2205   //     configuration macro
2206   // \li chains=<i>configuration</i>
2207   //     comma separated list of configurations to be run during simulation
2208   // \li rawfile=<i>file</i>
2209   //     source for the RawReader to be created, the default is <i>./</i> if
2210   //     raw data is simulated
2211
2212   int iResult=0;
2213
2214   if (!fpHLT && !CreateHLT()) {
2215     return kFALSE;
2216   }
2217   AliHLTSimulation* pHLT=fpHLT;
2218
2219   AliRunLoader* pRunLoader = LoadRun("READ");
2220   if (!pRunLoader) return kFALSE;
2221
2222   // initialize CDB storage, run number, set CDB lock
2223   // thats for the case of running HLT simulation without all the other steps
2224   // multiple calls are handled by the function, so we can just call
2225   InitCDB();
2226   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
2227   SetCDBLock();
2228   
2229   // init the HLT simulation
2230   TString options;
2231   if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2232   TString detStr = fWriteRawData;
2233   if (!IsSelected("HLT", detStr)) {
2234     options+=" writerawfiles=";
2235   } else {
2236     options+=" writerawfiles=HLT";
2237   }
2238
2239   if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2240     // as a matter of fact, HLT will run reconstruction and needs the RawReader
2241     // in order to get detector data. By default, RawReaderFile is used to read
2242     // the already simulated ddl files. Date and Root files from the raw data
2243     // are generated after the HLT simulation.
2244     options+=" rawfile=./";
2245   }
2246
2247   AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2248   if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2249     AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2250   } else {
2251     // run the HLT simulation
2252     AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2253     if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2254       AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2255     }
2256   }
2257
2258   // delete the instance
2259   AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2260   if (fctDelete==NULL || fctDelete(pHLT)<0) {
2261     AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2262   }
2263   pHLT=NULL;
2264
2265   return iResult>=0?kTRUE:kFALSE;
2266 }
2267
2268 //_____________________________________________________________________________
2269 Bool_t AliSimulation::RunQA()
2270 {
2271         // run the QA on summable hits, digits or digits
2272         
2273     //if(!gAlice) return kFALSE;
2274         AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2275
2276         TString detectorsw("") ;  
2277         Bool_t rv = kTRUE ; 
2278   AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2279         detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ; 
2280         if ( detectorsw.IsNull() ) 
2281                 rv = kFALSE ; 
2282         return rv ; 
2283 }
2284
2285 //_____________________________________________________________________________
2286 Bool_t AliSimulation::SetRunQA(TString detAndAction) 
2287 {
2288         // Allows to run QA for a selected set of detectors
2289         // and a selected set of tasks among HITS, SDIGITS and DIGITS
2290         // all selected detectors run the same selected tasks
2291         
2292         if (!detAndAction.Contains(":")) {
2293                 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2294                 fRunQA = kFALSE ;
2295                 return kFALSE ;                 
2296         }
2297         Int_t colon = detAndAction.Index(":") ; 
2298         fQADetectors = detAndAction(0, colon) ; 
2299         if (fQADetectors.Contains("ALL") ){
2300     TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ; 
2301     Int_t minus = fQADetectors.Last('-') ; 
2302     TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ; 
2303     TString toRemove("") ;
2304     while (minus >= 0) {
2305       toRemove = fQADetectors(minus+1, fQADetectors.Length()) ; 
2306       toRemove = toRemove.Strip() ; 
2307       toKeep.ReplaceAll(toRemove, "") ; 
2308       fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ; 
2309       minus = fQADetectors.Last('-') ; 
2310     }
2311     fQADetectors = toKeep ; 
2312   }
2313   fQATasks   = detAndAction(colon+1, detAndAction.Sizeof() ) ; 
2314         if (fQATasks.Contains("ALL") ) {
2315                 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ; 
2316         } else {
2317                 fQATasks.ToUpper() ; 
2318                 TString tempo("") ; 
2319                 if ( fQATasks.Contains("HIT") ) 
2320                         tempo = Form("%d ", AliQAv1::kHITS) ; 
2321                 if ( fQATasks.Contains("SDIGIT") ) 
2322                         tempo += Form("%d ", AliQAv1::kSDIGITS) ; 
2323                 if ( fQATasks.Contains("DIGIT") ) 
2324                         tempo += Form("%d ", AliQAv1::kDIGITS) ; 
2325                 fQATasks = tempo ; 
2326                 if (fQATasks.IsNull()) {
2327                         AliInfo("No QA requested\n")  ;
2328                         fRunQA = kFALSE ;
2329                         return kTRUE ; 
2330                 }
2331         }       
2332         TString tempo(fQATasks) ; 
2333     tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS))  ;
2334     tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;  
2335     tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;    
2336         AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;  
2337         fRunQA = kTRUE ;
2338         AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ; 
2339         AliQAManager::QAManager()->SetTasks(fQATasks) ; 
2340   for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) 
2341     AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2342   
2343         return kTRUE; 
2344
2345
2346 //_____________________________________________________________________________
2347 void AliSimulation::ProcessEnvironmentVars()
2348 {
2349 // Extract run number and random generator seed from env variables
2350
2351     AliInfo("Processing environment variables");
2352     
2353     // Random Number seed
2354     
2355     // first check that seed is not already set
2356     if (fSeed == 0) {
2357         if (gSystem->Getenv("CONFIG_SEED")) {
2358                 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2359         }
2360     } else {
2361         if (gSystem->Getenv("CONFIG_SEED")) {
2362                 AliInfo(Form("Seed for random number generation already set (%d)"
2363                              ": CONFIG_SEED variable ignored!", fSeed));
2364         }
2365     }
2366    
2367     AliInfo(Form("Seed for random number generation = %d ", fSeed)); 
2368
2369     // Run Number
2370     
2371     // first check that run number is not already set
2372     if(fRun < 0) {    
2373         if (gSystem->Getenv("DC_RUN")) {
2374                 fRun = atoi(gSystem->Getenv("DC_RUN"));
2375         }
2376     } else {
2377         if (gSystem->Getenv("DC_RUN")) {
2378                 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2379         }
2380     }
2381     
2382     AliInfo(Form("Run number = %d", fRun)); 
2383 }
2384
2385 //---------------------------------------------------------------------
2386 void AliSimulation::WriteGRPEntry()
2387 {
2388   // Get the necessary information from galice (generator, trigger etc) and
2389   // write a GRP entry corresponding to the settings in the Config.C used
2390   // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2391
2392
2393   AliInfo("Writing global run parameters entry into the OCDB");
2394
2395   AliGRPObject* grpObj = new AliGRPObject();
2396
2397   grpObj->SetRunType("PHYSICS");
2398   grpObj->SetTimeStart(fTimeStart);
2399   grpObj->SetTimeEnd(fTimeEnd); 
2400   grpObj->SetBeamEnergyIsSqrtSHalfGeV(); // new format of GRP: store sqrt(s)/2 in GeV
2401
2402   const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2403   Int_t a = 0;
2404   Int_t z = 0;
2405
2406   if (gen) {
2407     TString projectile;
2408     gen->GetProjectile(projectile,a,z);
2409     TString target;
2410     gen->GetTarget(target,a,z);
2411     TString beamType = projectile + "-" + target;
2412     beamType.ReplaceAll(" ","");
2413     if (!beamType.CompareTo("-")) {
2414       grpObj->SetBeamType("UNKNOWN");
2415       grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
2416     }
2417     else {
2418       grpObj->SetBeamType(beamType);
2419       if (z != 0) {
2420           grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 * a / z);
2421       } else {
2422           grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 );
2423       }
2424       // Heavy ion run, the event specie is set to kHighMult
2425       fEventSpecie = AliRecoParam::kHighMult;
2426       if ((strcmp(beamType,"p-p") == 0) ||
2427           (strcmp(beamType,"p-")  == 0) ||
2428           (strcmp(beamType,"-p")  == 0) ||
2429           (strcmp(beamType,"P-P") == 0) ||
2430           (strcmp(beamType,"P-")  == 0) ||
2431           (strcmp(beamType,"-P")  == 0)) {
2432         // Proton run, the event specie is set to kLowMult
2433         fEventSpecie = AliRecoParam::kLowMult;
2434       } 
2435     }
2436   } else {
2437     AliWarning("Unknown beam type and energy! Setting energy to 0");
2438     grpObj->SetBeamEnergy(0);
2439     grpObj->SetBeamType("UNKNOWN");
2440   }
2441
2442   UInt_t detectorPattern  = 0;
2443   Int_t nDets = 0;
2444   TObjArray *detArray = gAlice->Detectors();
2445   for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2446     if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2447       AliDebug(1, Form("Detector #%d found: %s", iDet, AliDAQ::OfflineModuleName(iDet)));
2448       detectorPattern |= (1 << iDet);
2449       nDets++;
2450     }
2451   }
2452   // CTP
2453   if (!fTriggerConfig.IsNull())
2454     detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2455
2456   // HLT
2457   if (!fRunHLT.IsNull())
2458     detectorPattern |= (1 << AliDAQ::kHLTId);
2459
2460   grpObj->SetNumberOfDetectors((Char_t)nDets);
2461   grpObj->SetDetectorMask((Int_t)detectorPattern);
2462   grpObj->SetLHCPeriod("LHC08c");
2463   grpObj->SetLHCState("STABLE_BEAMS");
2464   //
2465   AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2466   Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2467
2468   Float_t factorSol     = field ? field->GetFactorSol() : 0;
2469   Float_t currentSol    = TMath::Abs(factorSol)>1E-6 ? 
2470     TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
2471   //
2472   Float_t factorDip     = field ? field->GetFactorDip() : 0;
2473   Float_t currentDip    = 6000.*TMath::Abs(factorDip);
2474   //
2475   grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
2476   grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);  
2477   grpObj->SetL3Polarity(factorSol>0 ? 0:1);  
2478   grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
2479   if (field) grpObj->SetUniformBMap(field->IsUniform()); // for special MC with k5kGUniform map
2480   grpObj->SetPolarityConventionLHC();                    // LHC convention +/+ current -> -/- field main components
2481   //
2482   grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2483   
2484   //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2485
2486   // Now store the entry in OCDB
2487   AliCDBManager* man = AliCDBManager::Instance();
2488   
2489   man->SetLock(0, fKey);
2490   
2491   AliCDBStorage* sto = man->GetStorage(fGRPWriteLocation.Data());
2492   
2493
2494   AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2495   AliCDBMetaData *metadata= new AliCDBMetaData();
2496
2497   metadata->SetResponsible("alice-off@cern.ch");
2498   metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2499  
2500   sto->Put(grpObj,id,metadata);
2501   man->SetLock(1, fKey);
2502 }
2503
2504 //_____________________________________________________________________________
2505 time_t AliSimulation::GenerateTimeStamp() const
2506 {
2507   // Generate event time-stamp according to
2508   // SOR/EOR time from GRP
2509   if (fUseTimeStampFromCDB)
2510     return fTimeStart + gRandom->Integer(fTimeEnd-fTimeStart);
2511   else
2512     return 0;
2513 }