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