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