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