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