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