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