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