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