]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
4507d905ad302f8edcdf6b2e708fb6cd2031587d
[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 <TVirtualMCApplication.h>
109 #include <TGeoManager.h>
110 #include <TObjString.h>
111 #include <TSystem.h>
112 #include <TFile.h>
113
114 #include "AliCodeTimer.h"
115 #include "AliCDBStorage.h"
116 #include "AliCDBEntry.h"
117 #include "AliCDBManager.h"
118 #include "AliGeomManager.h"
119 #include "AliAlignObj.h"
120 #include "AliCentralTrigger.h"
121 #include "AliDAQ.h"
122 #include "AliDigitizer.h"
123 #include "AliGenerator.h"
124 #include "AliLog.h"
125 #include "AliModule.h"
126 #include "AliRun.h"
127 #include "AliRunDigitizer.h"
128 #include "AliRunLoader.h"
129 #include "AliSimulation.h"
130 #include "AliVertexGenFile.h"
131 #include "AliCentralTrigger.h"
132 #include "AliCTPRawData.h"
133 #include "AliRawReaderFile.h"
134 #include "AliRawReaderRoot.h"
135 #include "AliRawReaderDate.h"
136 #include "AliESD.h"
137 #include "AliHeader.h"
138 #include "AliGenEventHeader.h"
139 #include "AliMC.h"
140 #include "AliHLTSimulation.h"
141 #include "AliQADataMakerSteer.h"
142
143 ClassImp(AliSimulation)
144
145 AliSimulation *AliSimulation::fgInstance = 0;
146 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
147
148 //_____________________________________________________________________________
149 AliSimulation::AliSimulation(const char* configFileName,
150                              const char* name, const char* title) :
151   TNamed(name, title),
152
153   fRunGeneration(kTRUE),
154   fRunSimulation(kTRUE),
155   fLoadAlignFromCDB(kTRUE),
156   fLoadAlObjsListOfDets("ALL"),
157   fMakeSDigits("ALL"),
158   fMakeDigits("ALL"),
159   fMakeTrigger(""),
160   fMakeDigitsFromHits(""),
161   fWriteRawData(""),
162   fRawDataFileName(""),
163   fDeleteIntermediateFiles(kFALSE),
164   fStopOnError(kFALSE),
165
166   fNEvents(1),
167   fConfigFileName(configFileName),
168   fGAliceFileName("galice.root"),
169   fEventsPerFile(),
170   fBkgrdFileNames(NULL),
171   fAlignObjArray(NULL),
172   fUseBkgrdVertex(kTRUE),
173   fRegionOfInterest(kFALSE),
174   fCDBUri(""),
175   fSpecCDBUri(),
176   fRun(-1),
177   fSeed(0),
178   fInitCDBCalled(kFALSE),
179   fInitRunNumberCalled(kFALSE),
180   fSetRunNumberFromDataCalled(kFALSE),
181   fEmbeddingFlag(kFALSE),
182   fRunQA(kTRUE), 
183   fRunHLT("default")
184 {
185 // create simulation object with default parameters
186   fgInstance = this;
187   SetGAliceFile("galice.root");
188   
189 // for QA
190    for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) 
191         fQACycles[iDet] = 999999;
192 }
193
194 //_____________________________________________________________________________
195 AliSimulation::AliSimulation(const AliSimulation& sim) :
196   TNamed(sim),
197
198   fRunGeneration(sim.fRunGeneration),
199   fRunSimulation(sim.fRunSimulation),
200   fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
201   fLoadAlObjsListOfDets(sim.fLoadAlObjsListOfDets),
202   fMakeSDigits(sim.fMakeSDigits),
203   fMakeDigits(sim.fMakeDigits),
204   fMakeTrigger(sim.fMakeTrigger),
205   fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
206   fWriteRawData(sim.fWriteRawData),
207   fRawDataFileName(""),
208   fDeleteIntermediateFiles(kFALSE),
209   fStopOnError(sim.fStopOnError),
210
211   fNEvents(sim.fNEvents),
212   fConfigFileName(sim.fConfigFileName),
213   fGAliceFileName(sim.fGAliceFileName),
214   fEventsPerFile(),
215   fBkgrdFileNames(NULL),
216   fAlignObjArray(NULL),
217   fUseBkgrdVertex(sim.fUseBkgrdVertex),
218   fRegionOfInterest(sim.fRegionOfInterest),
219   fCDBUri(sim.fCDBUri),
220   fSpecCDBUri(),
221   fRun(-1),
222   fSeed(0),
223   fInitCDBCalled(sim.fInitCDBCalled),
224   fInitRunNumberCalled(sim.fInitRunNumberCalled),
225   fSetRunNumberFromDataCalled(sim.fSetRunNumberFromDataCalled),
226   fEmbeddingFlag(sim.fEmbeddingFlag),
227   fRunQA(kTRUE), 
228   fRunHLT(sim.fRunHLT)
229 {
230 // copy constructor
231
232   for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
233     if (!sim.fEventsPerFile[i]) continue;
234     fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
235   }
236
237   fBkgrdFileNames = new TObjArray;
238   for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
239     if (!sim.fBkgrdFileNames->At(i)) continue;
240     fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
241   }
242
243   for (Int_t i = 0; i < sim.fSpecCDBUri.GetEntriesFast(); i++) {
244     if (sim.fSpecCDBUri[i]) fSpecCDBUri.Add(sim.fSpecCDBUri[i]->Clone());
245   }
246   fgInstance = this;
247
248 // for QA
249    for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) 
250         fQACycles[iDet] = sim.fQACycles[iDet];
251 }
252
253 //_____________________________________________________________________________
254 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
255 {
256 // assignment operator
257
258   this->~AliSimulation();
259   new(this) AliSimulation(sim);
260   return *this;
261 }
262
263 //_____________________________________________________________________________
264 AliSimulation::~AliSimulation()
265 {
266 // clean up
267
268   fEventsPerFile.Delete();
269 //  if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
270 //  delete fAlignObjArray; fAlignObjArray=0;
271
272   if (fBkgrdFileNames) {
273     fBkgrdFileNames->Delete();
274     delete fBkgrdFileNames;
275   }
276
277   fSpecCDBUri.Delete();
278   if (fgInstance==this) fgInstance = 0;
279
280   AliCodeTimer::Instance()->Print();
281 }
282
283
284 //_____________________________________________________________________________
285 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
286 {
287 // set the number of events for one run
288
289   fNEvents = nEvents;
290 }
291
292 //_____________________________________________________________________________
293 void AliSimulation::InitCDB()
294 {
295 // activate a default CDB storage
296 // First check if we have any CDB storage set, because it is used 
297 // to retrieve the calibration and alignment constants
298
299   if (fInitCDBCalled) return;
300   fInitCDBCalled = kTRUE;
301
302   AliCDBManager* man = AliCDBManager::Instance();
303   if (man->IsDefaultStorageSet())
304   {
305     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
306     AliWarning("Default CDB storage has been already set !");
307     AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
308     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
309     fCDBUri = man->GetDefaultStorage()->GetURI();
310   }
311   else {
312     if (fCDBUri.Length() > 0) 
313     {
314         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
315         AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
316         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
317     } else {
318         fCDBUri="local://$ALICE_ROOT";
319         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
320         AliWarning("Default CDB storage not yet set !!!!");
321         AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
322         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
323                 
324     }
325     man->SetDefaultStorage(fCDBUri);
326   }
327
328   // Now activate the detector specific CDB storage locations
329   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
330     TObject* obj = fSpecCDBUri[i];
331     if (!obj) continue;
332     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
333     AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
334     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
335     man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
336   }
337       
338 }
339
340 //_____________________________________________________________________________
341 void AliSimulation::InitRunNumber(){
342 // check run number. If not set, set it to 0 !!!!
343   
344   if (fInitRunNumberCalled) return;
345   fInitRunNumberCalled = kTRUE;
346   
347   AliCDBManager* man = AliCDBManager::Instance();
348   if (man->GetRun() >= 0)
349   {
350     if(fRun >= 0) {
351         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
352         AliWarning(Form("Run number is already set in AliCDBManager: %d !", man->GetRun()));
353         AliWarning(Form("Ignoring the run number declared in AliSimulation: %d", fRun));
354         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
355     }
356     fRun = man->GetRun();
357   }
358   else {
359     if(fRun >= 0) {
360         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
361         AliDebug(2, Form("Setting CDB run number to: %d",fRun));
362         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
363     } else {
364         fRun=0;
365         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
366         AliWarning("Run number not yet set !!!!");
367         AliWarning(Form("Setting it now to: %d", fRun));
368         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
369         
370     }
371     man->SetRun(fRun);
372   }
373
374   man->Print();
375
376 }
377
378 //_____________________________________________________________________________
379 void AliSimulation::SetCDBLock() {
380   // Set CDB lock: from now on it is forbidden to reset the run number
381   // or the default storage or to activate any further storage!
382   
383   AliCDBManager::Instance()->SetLock(1);
384 }
385
386 //_____________________________________________________________________________
387 void AliSimulation::SetDefaultStorage(const char* uri) {
388 // Store the desired default CDB storage location
389 // Activate it later within the Run() method
390
391   fCDBUri = uri;
392
393 }
394
395 //_____________________________________________________________________________
396 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
397 // Store a detector-specific CDB storage location
398 // Activate it later within the Run() method
399
400   AliCDBPath aPath(calibType);
401   if(!aPath.IsValid()){
402         AliError(Form("Not a valid path: %s", calibType));
403         return;
404   }
405
406   TObject* obj = fSpecCDBUri.FindObject(calibType);
407   if (obj) fSpecCDBUri.Remove(obj);
408   fSpecCDBUri.Add(new TNamed(calibType, uri));
409
410 }
411
412 //_____________________________________________________________________________
413 void AliSimulation::SetRunNumber(Int_t run)
414 {
415 // sets run number
416 // Activate it later within the Run() method
417
418         fRun = run;
419 }
420
421 //_____________________________________________________________________________
422 void AliSimulation::SetSeed(Int_t seed)
423 {
424 // sets seed number
425 // Activate it later within the Run() method
426
427         fSeed = seed;
428 }
429
430 //_____________________________________________________________________________
431 Bool_t AliSimulation::SetRunNumberFromData()
432 {
433   // Set the CDB manager run number
434   // The run number is retrieved from gAlice
435
436     if (fSetRunNumberFromDataCalled) return kTRUE;
437     fSetRunNumberFromDataCalled = kTRUE;    
438   
439     AliCDBManager* man = AliCDBManager::Instance();
440     Int_t runData = -1, runCDB = -1;
441   
442     AliRunLoader* runLoader = LoadRun("READ");
443     if (!runLoader) return kFALSE;
444     else {
445         runData = runLoader->GetAliRun()->GetHeader()->GetRun();
446         delete runLoader;
447     }
448   
449     runCDB = man->GetRun();
450     if(runCDB >= 0) {
451         if (runCDB != runData) {
452                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
453                 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
454                 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
455                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");    
456         }
457         
458     }
459       
460     man->SetRun(runData);
461     fRun = runData;
462     
463     if(man->GetRun() < 0) {
464         AliError("Run number not properly initalized!");
465         return kFALSE;
466     }
467   
468     man->Print();
469     
470     return kTRUE;
471 }
472
473 //_____________________________________________________________________________
474 void AliSimulation::SetConfigFile(const char* fileName)
475 {
476 // set the name of the config file
477
478   fConfigFileName = fileName;
479 }
480
481 //_____________________________________________________________________________
482 void AliSimulation::SetGAliceFile(const char* fileName)
483 {
484 // set the name of the galice file
485 // the path is converted to an absolute one if it is relative
486
487   fGAliceFileName = fileName;
488   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
489     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
490                                                 fGAliceFileName);
491     fGAliceFileName = absFileName;
492     delete[] absFileName;
493   }
494
495   AliDebug(2, Form("galice file name set to %s", fileName));
496 }
497
498 //_____________________________________________________________________________
499 void AliSimulation::SetEventsPerFile(const char* detector, const char* type, 
500                                      Int_t nEvents)
501 {
502 // set the number of events per file for the given detector and data type
503 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
504
505   TNamed* obj = new TNamed(detector, type);
506   obj->SetUniqueID(nEvents);
507   fEventsPerFile.Add(obj);
508 }
509
510 //_____________________________________________________________________________
511 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
512 {
513   // Read the alignment objects from CDB.
514   // Each detector is supposed to have the
515   // alignment objects in DET/Align/Data CDB path.
516   // All the detector objects are then collected,
517   // sorted by geometry level (starting from ALIC) and
518   // then applied to the TGeo geometry.
519   // Finally an overlaps check is performed.
520
521   if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
522     AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
523     return kFALSE;
524   }  
525   
526   // initialize CDB storage, run number, set CDB lock
527   InitCDB();
528 //  if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
529   SetCDBLock();
530     
531   Bool_t delRunLoader = kFALSE;
532   if (!runLoader) {
533     runLoader = LoadRun("READ");
534     if (!runLoader) return kFALSE;
535     delRunLoader = kTRUE;
536   }
537   
538   // Export ideal geometry 
539   if(!gAlice->IsRootGeometry()) AliGeomManager::GetGeometry()->Export("geometry.root");
540
541   // Load alignment data from CDB and apply to geometry through AliGeomManager
542   if(fLoadAlignFromCDB){
543     
544     TString detStr = fLoadAlObjsListOfDets;
545     TString loadAlObjsListOfDets = "";
546     
547     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
548     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
549       AliModule* det = (AliModule*) detArray->At(iDet);
550       if (!det || !det->IsActive()) continue;
551       if (IsSelected(det->GetName(), detStr)) {
552         //add det to list of dets to be aligned from CDB
553         loadAlObjsListOfDets += det->GetName();
554         loadAlObjsListOfDets += " ";
555       }
556     } // end loop over detectors
557     loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
558     AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
559   }else{
560     // Check if the array with alignment objects was
561     // provided by the user. If yes, apply the objects
562     // to the present TGeo geometry
563     if (fAlignObjArray) {
564       if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
565         AliError("The misalignment of one or more volumes failed!"
566                  "Compare the list of simulated detectors and the list of detector alignment data!");
567         if (delRunLoader) delete runLoader;
568         return kFALSE;
569       }
570     }
571   }
572
573   // Update the internal geometry of modules (ITS needs it)
574   TString detStr = fLoadAlObjsListOfDets;
575   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
576   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
577
578     AliModule* det = (AliModule*) detArray->At(iDet);
579     if (!det || !det->IsActive()) continue;
580     if (IsSelected(det->GetName(), detStr)) {
581       det->UpdateInternalGeometry();
582     }
583   } // end loop over detectors
584
585
586   if (delRunLoader) delete runLoader;
587
588   return kTRUE;
589 }
590
591 //_____________________________________________________________________________
592 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
593 {
594 // add a file with background events for merging
595
596   TObjString* fileNameStr = new TObjString(fileName);
597   fileNameStr->SetUniqueID(nSignalPerBkgrd);
598   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
599   fBkgrdFileNames->Add(fileNameStr);
600 }
601
602 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
603 {
604 // add a file with background events for embeddin
605   MergeWith(fileName, nSignalPerBkgrd);
606   fEmbeddingFlag = kTRUE;
607 }
608
609 //_____________________________________________________________________________
610 Bool_t AliSimulation::Run(Int_t nEvents)
611 {
612 // run the generation, simulation and digitization
613
614  
615   AliCodeTimerAuto("")
616   
617   // Load run number and seed from environmental vars
618   ProcessEnvironmentVars();
619
620   gRandom->SetSeed(fSeed);
621    
622   if (nEvents > 0) fNEvents = nEvents;
623
624   // generation and simulation -> hits
625   if (fRunGeneration) {
626     if (!RunSimulation()) if (fStopOnError) return kFALSE;
627   }
628            
629   // initialize CDB storage from external environment
630   // (either CDB manager or AliSimulation setters),
631   // if not already done in RunSimulation()
632   InitCDB();
633   
634   // Set run number in CDBManager from data 
635   // From this point on the run number must be always loaded from data!
636   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
637   
638   // Set CDB lock: from now on it is forbidden to reset the run number
639   // or the default storage or to activate any further storage!
640   SetCDBLock();
641
642   // If RunSimulation was not called, load the geometry and misalign it
643   if (!AliGeomManager::GetGeometry()) {
644     // Initialize the geometry manager
645     AliGeomManager::LoadGeometry("geometry.root");
646     if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
647     // Misalign geometry
648     if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
649   }
650
651   // hits -> summable digits
652   if (!fMakeSDigits.IsNull()) {
653     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
654  
655   }
656   
657
658   
659   // summable digits -> digits  
660   if (!fMakeDigits.IsNull()) {
661     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
662       if (fStopOnError) return kFALSE;
663     }
664    }
665
666   
667   
668   // hits -> digits
669   if (!fMakeDigitsFromHits.IsNull()) {
670     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
671       AliWarning(Form("Merging and direct creation of digits from hits " 
672                  "was selected for some detectors. "
673                  "No merging will be done for the following detectors: %s",
674                  fMakeDigitsFromHits.Data()));
675     }
676     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
677       if (fStopOnError) return kFALSE;
678     }
679   }
680
681   
682   
683   // digits -> trigger
684   if (!RunTrigger(fMakeTrigger)) {
685     if (fStopOnError) return kFALSE;
686   }
687
688   
689   
690   // digits -> raw data
691   if (!fWriteRawData.IsNull()) {
692     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
693                       fDeleteIntermediateFiles)) {
694       if (fStopOnError) return kFALSE;
695     }
696   }
697
698   
699   
700   // run HLT simulation
701   if (!fRunHLT.IsNull()) {
702     if (!RunHLT()) {
703       if (fStopOnError) return kFALSE;
704     }
705   }
706   
707  // //QA
708 //      if (fRunQA) {
709 //              Bool_t rv = RunQA() ; 
710 //              if (!rv)
711 //                      if (fStopOnError) 
712 //                              return kFALSE ;         
713 //      }
714
715   // Cleanup of CDB manager: cache and active storages!
716   AliCDBManager::Instance()->ClearCache();
717
718   return kTRUE;
719 }
720
721 //_____________________________________________________________________________
722 Bool_t AliSimulation::RunTrigger(const char* descriptors)
723 {
724   // run the trigger
725
726   AliCodeTimerAuto("")
727
728   // initialize CDB storage from external environment
729   // (either CDB manager or AliSimulation setters),
730   // if not already done in RunSimulation()
731   InitCDB();
732   
733   // Set run number in CDBManager from data 
734   // From this point on the run number must be always loaded from data!
735   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
736   
737   // Set CDB lock: from now on it is forbidden to reset the run number
738   // or the default storage or to activate any further storage!
739   SetCDBLock();
740    
741    AliRunLoader* runLoader = LoadRun("READ");
742    if (!runLoader) return kFALSE;
743    TString des = descriptors;
744
745    if (des.IsNull()) {
746      if (gAlice->GetTriggerDescriptor() != "") {
747        des = gAlice->GetTriggerDescriptor();
748      }
749      else {
750        AliWarning("No trigger descriptor is specified. Skipping the trigger simulation...");
751        return kTRUE;
752      }
753    }
754
755    runLoader->MakeTree( "GG" );
756    AliCentralTrigger* aCTP = runLoader->GetTrigger();
757   // Load Descriptors
758    aCTP->LoadDescriptor( des );
759
760   // digits -> trigger
761    if( !aCTP->RunTrigger( runLoader ) ) {
762       if (fStopOnError) {
763     //  delete aCTP;
764          return kFALSE;
765       }
766    }
767
768    delete runLoader;
769
770    return kTRUE;
771 }
772
773 //_____________________________________________________________________________
774 Bool_t AliSimulation::WriteTriggerRawData()
775 {
776   // Writes the CTP (trigger) DDL raw data
777   // Details of the format are given in the
778   // trigger TDR - pages 134 and 135.
779   AliCTPRawData writer;
780   writer.RawData();
781
782   return kTRUE;
783 }
784
785 //_____________________________________________________________________________
786 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
787 {
788 // run the generation and simulation
789
790   AliCodeTimerAuto("")
791
792   // initialize CDB storage and run number from external environment
793   // (either CDB manager or AliSimulation setters)
794   InitCDB();
795   InitRunNumber();
796   SetCDBLock();
797   
798   if (!gAlice) {
799     AliError("no gAlice object. Restart aliroot and try again.");
800     return kFALSE;
801   }
802   if (gAlice->Modules()->GetEntries() > 0) {
803     AliError("gAlice was already run. Restart aliroot and try again.");
804     return kFALSE;
805   }
806
807   AliInfo(Form("initializing gAlice with config file %s",
808           fConfigFileName.Data()));
809   StdoutToAliInfo(StderrToAliError(
810     gAlice->Init(fConfigFileName.Data());
811   ););
812   
813   // Get the trigger descriptor string
814   // Either from AliSimulation or from
815   // gAlice
816   if (fMakeTrigger.IsNull()) {
817     if (gAlice->GetTriggerDescriptor() != "")
818       fMakeTrigger = gAlice->GetTriggerDescriptor();
819   }
820   else
821     gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
822
823   // Set run number in CDBManager
824   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
825
826   AliRunLoader* runLoader = gAlice->GetRunLoader();
827   if (!runLoader) {
828              AliError(Form("gAlice has no run loader object. "
829                              "Check your config file: %s", fConfigFileName.Data()));
830              return kFALSE;
831   }
832   SetGAliceFile(runLoader->GetFileName());
833       
834   // Misalign geometry
835 #if ROOT_VERSION_CODE < 331527
836   AliGeomManager::SetGeometry(gGeoManager);
837   MisalignGeometry(runLoader);
838 #endif
839
840 //   AliRunLoader* runLoader = gAlice->GetRunLoader();
841 //   if (!runLoader) {
842 //     AliError(Form("gAlice has no run loader object. "
843 //                   "Check your config file: %s", fConfigFileName.Data()));
844 //     return kFALSE;
845 //   }
846 //   SetGAliceFile(runLoader->GetFileName());
847
848   if (!gAlice->Generator()) {
849     AliError(Form("gAlice has no generator object. "
850                   "Check your config file: %s", fConfigFileName.Data()));
851     return kFALSE;
852   }
853   if (nEvents <= 0) nEvents = fNEvents;
854
855   // get vertex from background file in case of merging
856   if (fUseBkgrdVertex &&
857       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
858     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
859     const char* fileName = ((TObjString*)
860                             (fBkgrdFileNames->At(0)))->GetName();
861     AliInfo(Form("The vertex will be taken from the background "
862                  "file %s with nSignalPerBackground = %d", 
863                  fileName, signalPerBkgrd));
864     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
865     gAlice->Generator()->SetVertexGenerator(vtxGen);
866   }
867
868   if (!fRunSimulation) {
869     gAlice->Generator()->SetTrackingFlag(0);
870   }
871
872   // set the number of events per file for given detectors and data types
873   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
874     if (!fEventsPerFile[i]) continue;
875     const char* detName = fEventsPerFile[i]->GetName();
876     const char* typeName = fEventsPerFile[i]->GetTitle();
877     TString loaderName(detName);
878     loaderName += "Loader";
879     AliLoader* loader = runLoader->GetLoader(loaderName);
880     if (!loader) {
881       AliError(Form("RunSimulation", "no loader for %s found\n"
882                     "Number of events per file not set for %s %s", 
883                     detName, typeName, detName));
884       continue;
885     }
886     AliDataLoader* dataLoader = 
887       loader->GetDataLoader(typeName);
888     if (!dataLoader) {
889       AliError(Form("no data loader for %s found\n"
890                     "Number of events per file not set for %s %s", 
891                     typeName, detName, typeName));
892       continue;
893     }
894     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
895     AliDebug(1, Form("number of events per file set to %d for %s %s",
896                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
897   }
898
899   AliInfo("running gAlice");
900   StdoutToAliInfo(StderrToAliError(
901     gAlice->Run(nEvents);
902   ););
903
904   delete runLoader;
905
906   return kTRUE;
907 }
908
909 //_____________________________________________________________________________
910 Bool_t AliSimulation::RunSDigitization(const char* detectors)
911 {
912 // run the digitization and produce summable digits
913
914   AliCodeTimerAuto("")
915
916   // initialize CDB storage, run number, set CDB lock
917   InitCDB();
918   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
919   SetCDBLock();
920   
921   AliRunLoader* runLoader = LoadRun();
922   if (!runLoader) return kFALSE;
923
924   TString detStr = detectors;
925   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
926   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
927     AliModule* det = (AliModule*) detArray->At(iDet);
928     if (!det || !det->IsActive()) continue;
929     if (IsSelected(det->GetName(), detStr)) {
930       AliInfo(Form("creating summable digits for %s", det->GetName()));
931       AliCodeTimerAuto(Form("creating summable digits for %s", det->GetName()));
932           
933       det->Hits2SDigits();
934     }
935   }
936
937   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
938     AliError(Form("the following detectors were not found: %s",
939                   detStr.Data()));
940     if (fStopOnError) return kFALSE;
941   }
942
943   delete runLoader;
944
945   return kTRUE;
946 }
947
948
949 //_____________________________________________________________________________
950 Bool_t AliSimulation::RunDigitization(const char* detectors, 
951                                       const char* excludeDetectors)
952 {
953 // run the digitization and produce digits from sdigits
954
955   AliCodeTimerAuto("")
956
957   // initialize CDB storage, run number, set CDB lock
958   InitCDB();
959   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
960   SetCDBLock();
961   
962   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
963   if (gAlice) delete gAlice;
964   gAlice = NULL;
965
966   Int_t nStreams = 1;
967   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
968   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
969   AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
970   // manager->SetEmbeddingFlag(fEmbeddingFlag);
971   manager->SetInputStream(0, fGAliceFileName.Data());
972   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
973     const char* fileName = ((TObjString*)
974                             (fBkgrdFileNames->At(iStream-1)))->GetName();
975     manager->SetInputStream(iStream, fileName);
976   }
977
978   TString detStr = detectors;
979   TString detExcl = excludeDetectors;
980   manager->GetInputStream(0)->ImportgAlice();
981   AliRunLoader* runLoader = 
982     AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
983   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
984   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
985     AliModule* det = (AliModule*) detArray->At(iDet);
986     if (!det || !det->IsActive()) continue;
987     if (IsSelected(det->GetName(), detStr) && 
988         !IsSelected(det->GetName(), detExcl)) {
989       AliDigitizer* digitizer = det->CreateDigitizer(manager);
990       
991       if (!digitizer) {
992         AliError(Form("no digitizer for %s", det->GetName()));
993         if (fStopOnError) return kFALSE;
994       } else {
995         digitizer->SetRegionOfInterest(fRegionOfInterest);
996       }
997     }
998   }
999
1000   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1001     AliError(Form("the following detectors were not found: %s", 
1002                   detStr.Data()));
1003     if (fStopOnError) return kFALSE;
1004   }
1005
1006   if (!manager->GetListOfTasks()->IsEmpty()) {
1007     AliInfo("executing digitization");
1008     manager->Exec("");
1009   }
1010
1011   delete manager;
1012
1013   return kTRUE;
1014 }
1015
1016 //_____________________________________________________________________________
1017 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1018 {
1019 // run the digitization and produce digits from hits
1020
1021   AliCodeTimerAuto("")
1022
1023   // initialize CDB storage, run number, set CDB lock
1024   InitCDB();
1025   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1026   SetCDBLock();
1027   
1028   AliRunLoader* runLoader = LoadRun("READ");
1029   if (!runLoader) return kFALSE;
1030
1031   TString detStr = detectors;
1032   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1033   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1034     AliModule* det = (AliModule*) detArray->At(iDet);
1035     if (!det || !det->IsActive()) continue;
1036     if (IsSelected(det->GetName(), detStr)) {
1037       AliInfo(Form("creating digits from hits for %s", det->GetName()));
1038       det->Hits2Digits();
1039     }
1040   }
1041
1042   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1043     AliError(Form("the following detectors were not found: %s", 
1044                   detStr.Data()));
1045     if (fStopOnError) return kFALSE;
1046   }
1047
1048   delete runLoader;
1049   //PH Temporary fix to avoid interference with the PHOS loder/getter
1050   //PH The problem has to be solved in more general way 09/06/05
1051
1052   return kTRUE;
1053 }
1054
1055 //_____________________________________________________________________________
1056 Bool_t AliSimulation::WriteRawData(const char* detectors, 
1057                                    const char* fileName,
1058                                    Bool_t deleteIntermediateFiles)
1059 {
1060 // convert the digits to raw data
1061 // First DDL raw data files for the given detectors are created.
1062 // If a file name is given, the DDL files are then converted to a DATE file.
1063 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
1064 // afterwards.
1065 // If the file name has the extension ".root", the DATE file is converted
1066 // to a root file.
1067 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1068
1069   AliCodeTimerAuto("")
1070
1071   if (!WriteRawFiles(detectors)) {
1072     if (fStopOnError) return kFALSE;
1073   }
1074
1075   TString dateFileName(fileName);
1076   if (!dateFileName.IsNull()) {
1077     Bool_t rootOutput = dateFileName.EndsWith(".root");
1078     if (rootOutput) dateFileName += ".date";
1079     if (!ConvertRawFilesToDate(dateFileName)) {
1080       if (fStopOnError) return kFALSE;
1081     }
1082     if (deleteIntermediateFiles) {
1083       AliRunLoader* runLoader = LoadRun("READ");
1084       if (runLoader) for (Int_t iEvent = 0; 
1085                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1086         char command[256];
1087         sprintf(command, "rm -r raw%d", iEvent);
1088         gSystem->Exec(command);
1089       }
1090     }
1091
1092     if (rootOutput) {
1093       if (!ConvertDateToRoot(dateFileName, fileName)) {
1094         if (fStopOnError) return kFALSE;
1095       }
1096       if (deleteIntermediateFiles) {
1097         gSystem->Unlink(dateFileName);
1098       }
1099     }
1100   }
1101
1102   return kTRUE;
1103 }
1104
1105 //_____________________________________________________________________________
1106 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1107 {
1108 // convert the digits to raw data DDL files
1109
1110   AliCodeTimerAuto("")
1111   
1112   AliRunLoader* runLoader = LoadRun("READ");
1113   if (!runLoader) return kFALSE;
1114
1115   // write raw data to DDL files
1116   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1117     AliInfo(Form("processing event %d", iEvent));
1118     runLoader->GetEvent(iEvent);
1119     TString baseDir = gSystem->WorkingDirectory();
1120     char dirName[256];
1121     sprintf(dirName, "raw%d", iEvent);
1122     gSystem->MakeDirectory(dirName);
1123     if (!gSystem->ChangeDirectory(dirName)) {
1124       AliError(Form("couldn't change to directory %s", dirName));
1125       if (fStopOnError) return kFALSE; else continue;
1126     }
1127
1128     TString detStr = detectors;
1129     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1130     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1131       AliModule* det = (AliModule*) detArray->At(iDet);
1132       if (!det || !det->IsActive()) continue;
1133       if (IsSelected(det->GetName(), detStr)) {
1134         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1135         det->Digits2Raw();
1136       }
1137     }
1138
1139     if (!WriteTriggerRawData())
1140       if (fStopOnError) return kFALSE;
1141
1142     gSystem->ChangeDirectory(baseDir);
1143     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1144       AliError(Form("the following detectors were not found: %s", 
1145                     detStr.Data()));
1146       if (fStopOnError) return kFALSE;
1147     }
1148   }
1149
1150   delete runLoader;
1151   
1152   return kTRUE;
1153 }
1154
1155 //_____________________________________________________________________________
1156 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1157 {
1158 // convert raw data DDL files to a DATE file with the program "dateStream"
1159
1160   AliCodeTimerAuto("")
1161   
1162   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1163   if (!path) {
1164     AliError("the program dateStream was not found");
1165     if (fStopOnError) return kFALSE;
1166   } else {
1167     delete[] path;
1168   }
1169
1170   AliRunLoader* runLoader = LoadRun("READ");
1171   if (!runLoader) return kFALSE;
1172
1173   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1174   char command[256];
1175   // Note the option -s. It is used in order to avoid
1176   // the generation of SOR/EOR events.
1177   sprintf(command, "dateStream -s -D -o %s -# %d -C", 
1178           dateFileName, runLoader->GetNumberOfEvents());
1179   FILE* pipe = gSystem->OpenPipe(command, "w");
1180
1181   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1182     fprintf(pipe, "GDC\n");
1183     Float_t ldc = 0;
1184     Int_t prevLDC = -1;
1185
1186     // loop over detectors and DDLs
1187     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1188       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1189
1190         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1191         Int_t ldcID = Int_t(ldc + 0.0001);
1192         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1193
1194         char rawFileName[256];
1195         sprintf(rawFileName, "raw%d/%s", 
1196                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1197
1198         // check existence and size of raw data file
1199         FILE* file = fopen(rawFileName, "rb");
1200         if (!file) continue;
1201         fseek(file, 0, SEEK_END);
1202         unsigned long size = ftell(file);
1203         fclose(file);
1204         if (!size) continue;
1205
1206         if (ldcID != prevLDC) {
1207           fprintf(pipe, " LDC Id %d\n", ldcID);
1208           prevLDC = ldcID;
1209         }
1210         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1211       }
1212     }
1213   }
1214
1215   Int_t result = gSystem->ClosePipe(pipe);
1216
1217   delete runLoader;
1218   return (result == 0);
1219 }
1220
1221 //_____________________________________________________________________________
1222 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1223                                         const char* rootFileName)
1224 {
1225 // convert a DATE file to a root file with the program "alimdc"
1226
1227   // ALIMDC setup
1228   const Int_t kDBSize = 2000000000;
1229   const Int_t kTagDBSize = 1000000000;
1230   const Bool_t kFilter = kFALSE;
1231   const Int_t kCompression = 1;
1232
1233   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1234   if (!path) {
1235     AliError("the program alimdc was not found");
1236     if (fStopOnError) return kFALSE;
1237   } else {
1238     delete[] path;
1239   }
1240
1241   AliInfo(Form("converting DATE file %s to root file %s", 
1242                dateFileName, rootFileName));
1243
1244   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1245   const char* tagDBFS    = "/tmp/mdc1/tags";
1246
1247   // User defined file system locations
1248   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1249     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1250   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1251     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1252   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1253     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1254
1255   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1256   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1257   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1258
1259   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1260   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1261   gSystem->Exec(Form("mkdir %s",tagDBFS));
1262
1263   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1264                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1265   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1266
1267   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1268   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1269   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1270
1271   return (result == 0);
1272 }
1273
1274
1275 //_____________________________________________________________________________
1276 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1277 {
1278 // delete existing run loaders, open a new one and load gAlice
1279
1280   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1281   AliRunLoader* runLoader = 
1282     AliRunLoader::Open(fGAliceFileName.Data(), 
1283                        AliConfig::GetDefaultEventFolderName(), mode);
1284   if (!runLoader) {
1285     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1286     return NULL;
1287   }
1288   runLoader->LoadgAlice();
1289   runLoader->LoadHeader();
1290   gAlice = runLoader->GetAliRun();
1291   if (!gAlice) {
1292     AliError(Form("no gAlice object found in file %s", 
1293                   fGAliceFileName.Data()));
1294     return NULL;
1295   }
1296   return runLoader;
1297 }
1298
1299 //_____________________________________________________________________________
1300 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1301 {
1302 // get or calculate the number of signal events per background event
1303
1304   if (!fBkgrdFileNames) return 1;
1305   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1306   if (nBkgrdFiles == 0) return 1;
1307
1308   // get the number of signal events
1309   if (nEvents <= 0) {
1310     AliRunLoader* runLoader = 
1311         AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1312     if (!runLoader) return 1;
1313     
1314     nEvents = runLoader->GetNumberOfEvents();
1315     delete runLoader;
1316   }
1317
1318   Int_t result = 0;
1319   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1320     // get the number of background events
1321     const char* fileName = ((TObjString*)
1322                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1323     AliRunLoader* runLoader =
1324       AliRunLoader::Open(fileName, "BKGRD");
1325     if (!runLoader) continue;
1326     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1327     delete runLoader;
1328   
1329     // get or calculate the number of signal per background events
1330     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1331     if (nSignalPerBkgrd <= 0) {
1332       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1333     } else if (result && (result != nSignalPerBkgrd)) {
1334       AliInfo(Form("the number of signal events per background event "
1335                    "will be changed from %d to %d for stream %d", 
1336                    nSignalPerBkgrd, result, iBkgrdFile+1));
1337       nSignalPerBkgrd = result;
1338     }
1339
1340     if (!result) result = nSignalPerBkgrd;
1341     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1342       AliWarning(Form("not enough background events (%d) for %d signal events "
1343                       "using %d signal per background events for stream %d",
1344                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1345     }
1346   }
1347
1348   return result;
1349 }
1350
1351 //_____________________________________________________________________________
1352 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1353 {
1354 // check whether detName is contained in detectors
1355 // if yes, it is removed from detectors
1356
1357   // check if all detectors are selected
1358   if ((detectors.CompareTo("ALL") == 0) ||
1359       detectors.BeginsWith("ALL ") ||
1360       detectors.EndsWith(" ALL") ||
1361       detectors.Contains(" ALL ")) {
1362     detectors = "ALL";
1363     return kTRUE;
1364   }
1365
1366   // search for the given detector
1367   Bool_t result = kFALSE;
1368   if ((detectors.CompareTo(detName) == 0) ||
1369       detectors.BeginsWith(detName+" ") ||
1370       detectors.EndsWith(" "+detName) ||
1371       detectors.Contains(" "+detName+" ")) {
1372     detectors.ReplaceAll(detName, "");
1373     result = kTRUE;
1374   }
1375
1376   // clean up the detectors string
1377   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1378   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1379   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1380
1381   return result;
1382 }
1383
1384 //_____________________________________________________________________________
1385 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName) 
1386 {
1387 //
1388 // Steering routine  to convert raw data in directory rawDirectory/ to fake SDigits. 
1389 // These can be used for embedding of MC tracks into RAW data using the standard 
1390 // merging procedure.
1391 //
1392 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1393 //
1394     if (!gAlice) {
1395         AliError("no gAlice object. Restart aliroot and try again.");
1396         return kFALSE;
1397     }
1398     if (gAlice->Modules()->GetEntries() > 0) {
1399         AliError("gAlice was already run. Restart aliroot and try again.");
1400         return kFALSE;
1401     }
1402     
1403     AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1404     StdoutToAliInfo(StderrToAliError(gAlice->Init(fConfigFileName.Data());););
1405 //
1406 //  Initialize CDB     
1407     InitCDB();
1408     //AliCDBManager* man = AliCDBManager::Instance();
1409     //man->SetRun(0); // Should this come from rawdata header ?
1410     
1411     Int_t iDet;
1412     //
1413     // Get the runloader
1414     AliRunLoader* runLoader = gAlice->GetRunLoader();
1415     //
1416     // Open esd file if available
1417     TFile* esdFile = TFile::Open(esdFileName);
1418     Bool_t esdOK = (esdFile != 0);
1419     AliESD* esd = new AliESD;
1420     TTree* treeESD = 0;
1421     if (esdOK) {
1422         treeESD = (TTree*) esdFile->Get("esdTree");
1423         if (!treeESD) {
1424             AliWarning("No ESD tree found");
1425             esdOK = kFALSE;
1426         } else {
1427             treeESD->SetBranchAddress("ESD", &esd);
1428         }
1429     }
1430     //
1431     // Create the RawReader
1432     TString fileName(rawDirectory);
1433     AliRawReader* rawReader = 0x0;
1434     if (fileName.EndsWith("/")) {
1435       rawReader = new AliRawReaderFile(fileName);
1436     } else if (fileName.EndsWith(".root")) {
1437       rawReader = new AliRawReaderRoot(fileName);
1438     } else if (!fileName.IsNull()) {
1439       rawReader = new AliRawReaderDate(fileName);
1440       rawReader->SelectEvents(7);
1441     }
1442 //     if (!fEquipIdMap.IsNull() && fRawReader)
1443 //       fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1444     //
1445     // Get list of detectors
1446     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1447     //
1448     // Get Header
1449     AliHeader* header = runLoader->GetHeader();
1450     //
1451     TString detStr = fMakeSDigits;
1452     // Event loop
1453     Int_t nev = 0;
1454     while(kTRUE) {
1455         if (!(rawReader->NextEvent())) break;
1456         //
1457         // Detector loop
1458         for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1459             AliModule* det = (AliModule*) detArray->At(iDet);
1460             if (!det || !det->IsActive()) continue;
1461             if (IsSelected(det->GetName(), detStr)) {
1462               AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1463               det->Raw2SDigits(rawReader);
1464               rawReader->Reset();
1465             }
1466         } // detectors
1467
1468
1469         //
1470         //  If ESD information available obtain reconstructed vertex and store in header.
1471         if (esdOK) {
1472             treeESD->GetEvent(nev);
1473             const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1474             Double_t position[3];
1475             esdVertex->GetXYZ(position);
1476             AliGenEventHeader* mcHeader = new  AliGenEventHeader("ESD");
1477             TArrayF mcV;
1478             mcV.Set(3);
1479             for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1480             mcHeader->SetPrimaryVertex(mcV);
1481             header->Reset(0,nev);
1482             header->SetGenEventHeader(mcHeader);
1483             printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1484         }
1485         nev++;
1486 //
1487 //      Finish the event
1488         runLoader->TreeE()->Fill();
1489         runLoader->SetNextEvent();
1490     } // events
1491  
1492     delete rawReader;
1493 //
1494 //  Finish the run 
1495     runLoader->CdGAFile();
1496     runLoader->WriteHeader("OVERWRITE");
1497     runLoader->WriteRunLoader();
1498
1499     return kTRUE;
1500 }
1501
1502 //_____________________________________________________________________________
1503 Int_t AliSimulation::GetDetIndex(const char* detector)
1504 {
1505   // return the detector index corresponding to detector
1506   Int_t index = -1 ; 
1507   for (index = 0; index < fgkNDetectors ; index++) {
1508     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1509           break ; 
1510   }     
1511   return index ; 
1512 }
1513
1514 //_____________________________________________________________________________
1515 Bool_t AliSimulation::RunHLT()
1516 {
1517   // Run the HLT simulation
1518   // HLT simulation is implemented in HLT/sim/AliHLTSimulation
1519   // Disabled if fRunHLT is empty, default vaule is "default".
1520   // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
1521   // The default simulation depends on the HLT component libraries and their
1522   // corresponding agents which define components and chains to run. See
1523   // http://web.ift.uib.no/~kjeks/doc/alice-hlt/
1524   // http://web.ift.uib.no/~kjeks/doc/alice-hlt/classAliHLTModuleAgent.html
1525   //
1526   // The libraries to be loaded can be specified as an option.
1527   // <pre>
1528   // AliSimulation sim;
1529   // sim.SetRunHLT("libAliHLTSample.so");
1530   // </pre>
1531   // will only load <tt>libAliHLTSample.so</tt>
1532
1533   // Other available options:
1534   // \li loglevel=<i>level</i> <br>
1535   //     logging level for this processing
1536   // \li alilog=off
1537   //     disable redirection of log messages to AliLog class
1538   // \li config=<i>macro</i>
1539   //     configuration macro
1540   // \li localrec=<i>configuration</i>
1541   //     comma separated list of configurations to be run during simulation
1542
1543   int iResult=0;
1544   AliRunLoader* pRunLoader = LoadRun("READ");
1545   if (!pRunLoader) return kFALSE;
1546
1547   // initialize CDB storage, run number, set CDB lock
1548   InitCDB();
1549   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1550   SetCDBLock();
1551   
1552   // load the library dynamically
1553   gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1554
1555   // check for the library version
1556   AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
1557   if (!fctVersion) {
1558     AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
1559     return kFALSE;
1560   }
1561   if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
1562     AliError(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
1563     return kFALSE;
1564   }
1565
1566   // print compile info
1567   typedef void (*CompileInfo)( char*& date, char*& time);
1568   CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
1569   if (fctInfo) {
1570     char* date="";
1571     char* time="";
1572     (*fctInfo)(date, time);
1573     if (!date) date="unknown";
1574     if (!time) time="unknown";
1575     AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
1576   } else {
1577     AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
1578   }
1579
1580   // create instance of the HLT simulation
1581   AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
1582   AliHLTSimulation* pHLT=NULL;
1583   if (fctCreate==NULL || (pHLT=(fctCreate()))==NULL) {
1584     AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
1585     return kFALSE;    
1586   }
1587
1588   // init the HLT simulation
1589   if (fRunHLT.CompareTo("default")==0) fRunHLT="";
1590   AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
1591   if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, fRunHLT.Data())))<0) {
1592     AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
1593   } else {
1594     // run the HLT simulation
1595     AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
1596     if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
1597       AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
1598     }
1599   }
1600
1601   // delete the instance
1602   AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
1603   if (fctDelete==NULL || fctDelete(pHLT)<0) {
1604     AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
1605   }
1606   pHLT=NULL;
1607
1608   return iResult>=0?kTRUE:kFALSE;
1609 }
1610
1611 //_____________________________________________________________________________
1612 Bool_t AliSimulation::RunQA()
1613 {
1614         // run the QA on summable hits, digits or digits
1615
1616         AliQADataMakerSteer qas ; 
1617     qas.SetRunLoader(gAlice->GetRunLoader()) ;
1618
1619         Bool_t rv =  qas.Run("ALL", AliQA::kHITS) ; 
1620 //      qas.Reset() ; 
1621         rv *= qas.Run(fMakeSDigits.Data(), AliQA::kSDIGITS) ;   
1622 //      qas.Reset() ; 
1623         rv *= qas.Run(fMakeDigits.Data(), AliQA::kDIGITS) ;     
1624 //      qas.Reset() ; 
1625         rv *= qas.Run(fMakeDigitsFromHits.Data(), AliQA::kDIGITS) ; 
1626
1627         return rv ; 
1628 }
1629
1630 //_____________________________________________________________________________
1631 void AliSimulation::ProcessEnvironmentVars()
1632 {
1633 // Extract run number and random generator seed from env variables
1634
1635     AliInfo("Processing environment variables");
1636     
1637     // Random Number seed
1638     
1639     // first check that seed is not already set
1640     if (fSeed == 0) {
1641         if (gSystem->Getenv("CONFIG_SEED")) {
1642                 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
1643         }
1644     } else {
1645         if (gSystem->Getenv("CONFIG_SEED")) {
1646                 AliInfo(Form("Seed for random number generation already set (%d)"
1647                              ": CONFIG_SEED variable ignored!", fSeed));
1648         }
1649     }
1650    
1651     AliInfo(Form("Seed for random number generation = %d ", fSeed)); 
1652
1653     // Run Number
1654     
1655     // first check that run number is not already set
1656     if(fRun < 0) {    
1657         if (gSystem->Getenv("DC_RUN")) {
1658                 fRun = atoi(gSystem->Getenv("DC_RUN"));
1659         }
1660     } else {
1661         if (gSystem->Getenv("DC_RUN")) {
1662                 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
1663         }
1664     }
1665     
1666     AliInfo(Form("Run number = %d", fRun)); 
1667 }