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