]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
Do not export the misaligned geometry (Raffaele)
[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 <TStopwatch.h>
112 #include <TSystem.h>
113 #include <TFile.h>
114
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 "AliESD.h"
135 #include "AliHeader.h"
136 #include "AliGenEventHeader.h"
137 #include "AliMC.h"
138
139 ClassImp(AliSimulation)
140
141 AliSimulation *AliSimulation::fgInstance = 0;
142
143 //_____________________________________________________________________________
144 AliSimulation::AliSimulation(const char* configFileName, const char* cdbUri,
145                              const char* name, const char* title) :
146   TNamed(name, title),
147
148   fRunGeneration(kTRUE),
149   fRunSimulation(kTRUE),
150   fLoadAlignFromCDB(kTRUE),
151   fLoadAlObjsListOfDets("ALL"),
152   fMakeSDigits("ALL"),
153   fMakeDigits("ALL"),
154   fMakeTrigger(""),
155   fMakeDigitsFromHits(""),
156   fWriteRawData(""),
157   fRawDataFileName(""),
158   fDeleteIntermediateFiles(kFALSE),
159   fStopOnError(kFALSE),
160
161   fNEvents(1),
162   fConfigFileName(configFileName),
163   fGAliceFileName("galice.root"),
164   fEventsPerFile(),
165   fBkgrdFileNames(NULL),
166   fAlignObjArray(NULL),
167   fUseBkgrdVertex(kTRUE),
168   fRegionOfInterest(kFALSE),
169   fCDBUri(cdbUri),
170   fSpecCDBUri(),
171   fEmbeddingFlag(kFALSE)
172 {
173 // create simulation object with default parameters
174   fgInstance = this;
175   SetGAliceFile("galice.root");
176 }
177
178 //_____________________________________________________________________________
179 AliSimulation::AliSimulation(const AliSimulation& sim) :
180   TNamed(sim),
181
182   fRunGeneration(sim.fRunGeneration),
183   fRunSimulation(sim.fRunSimulation),
184   fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
185   fLoadAlObjsListOfDets(sim.fLoadAlObjsListOfDets),
186   fMakeSDigits(sim.fMakeSDigits),
187   fMakeDigits(sim.fMakeDigits),
188   fMakeTrigger(sim.fMakeTrigger),
189   fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
190   fWriteRawData(sim.fWriteRawData),
191   fRawDataFileName(""),
192   fDeleteIntermediateFiles(kFALSE),
193   fStopOnError(sim.fStopOnError),
194
195   fNEvents(sim.fNEvents),
196   fConfigFileName(sim.fConfigFileName),
197   fGAliceFileName(sim.fGAliceFileName),
198   fEventsPerFile(),
199   fBkgrdFileNames(NULL),
200   fAlignObjArray(NULL),
201   fUseBkgrdVertex(sim.fUseBkgrdVertex),
202   fRegionOfInterest(sim.fRegionOfInterest),
203   fCDBUri(sim.fCDBUri),
204   fSpecCDBUri(),
205   fEmbeddingFlag(sim.fEmbeddingFlag)
206 {
207 // copy constructor
208
209   for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
210     if (!sim.fEventsPerFile[i]) continue;
211     fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
212   }
213
214   fBkgrdFileNames = new TObjArray;
215   for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
216     if (!sim.fBkgrdFileNames->At(i)) continue;
217     fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
218   }
219
220   for (Int_t i = 0; i < sim.fSpecCDBUri.GetEntriesFast(); i++) {
221     if (sim.fSpecCDBUri[i]) fSpecCDBUri.Add(sim.fSpecCDBUri[i]->Clone());
222   }
223   fgInstance = this;
224 }
225
226 //_____________________________________________________________________________
227 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
228 {
229 // assignment operator
230
231   this->~AliSimulation();
232   new(this) AliSimulation(sim);
233   return *this;
234 }
235
236 //_____________________________________________________________________________
237 AliSimulation::~AliSimulation()
238 {
239 // clean up
240
241   fEventsPerFile.Delete();
242 //  if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
243 //  delete fAlignObjArray; fAlignObjArray=0;
244
245   if (fBkgrdFileNames) {
246     fBkgrdFileNames->Delete();
247     delete fBkgrdFileNames;
248   }
249
250   fSpecCDBUri.Delete();
251   if (fgInstance==this) fgInstance = 0;
252 }
253
254
255 //_____________________________________________________________________________
256 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
257 {
258 // set the number of events for one run
259
260   fNEvents = nEvents;
261 }
262
263 //_____________________________________________________________________________
264 void AliSimulation::InitCDBStorage()
265 {
266 // activate a default CDB storage
267 // First check if we have any CDB storage set, because it is used 
268 // to retrieve the calibration and alignment constants
269
270   AliCDBManager* man = AliCDBManager::Instance();
271   if (man->IsDefaultStorageSet())
272   {
273     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
274     AliWarning("Default CDB storage has been already set !");
275     AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
276     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
277     fCDBUri = "";
278   }
279   else {
280     AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
281     AliDebug(2, Form("Default CDB storage is set to: %s",fCDBUri.Data()));
282     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
283     man->SetDefaultStorage(fCDBUri);
284   }
285
286   // Now activate the detector specific CDB storage locations
287   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
288     TObject* obj = fSpecCDBUri[i];
289     if (!obj) continue;
290     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
291     AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
292     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
293     man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
294   }
295   man->Print();
296 }
297
298 //_____________________________________________________________________________
299 void AliSimulation::SetDefaultStorage(const char* uri) {
300 // Store the desired default CDB storage location
301 // Activate it later within the Run() method
302
303   fCDBUri = uri;
304
305 }
306
307 //_____________________________________________________________________________
308 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
309 // Store a detector-specific CDB storage location
310 // Activate it later within the Run() method
311
312   AliCDBPath aPath(calibType);
313   if(!aPath.IsValid()){
314         AliError(Form("Not a valid path: %s", calibType));
315         return;
316   }
317
318   TObject* obj = fSpecCDBUri.FindObject(calibType);
319   if (obj) fSpecCDBUri.Remove(obj);
320   fSpecCDBUri.Add(new TNamed(calibType, uri));
321
322 }
323
324 //_____________________________________________________________________________
325 void AliSimulation::SetConfigFile(const char* fileName)
326 {
327 // set the name of the config file
328
329   fConfigFileName = fileName;
330 }
331
332 //_____________________________________________________________________________
333 void AliSimulation::SetGAliceFile(const char* fileName)
334 {
335 // set the name of the galice file
336 // the path is converted to an absolute one if it is relative
337
338   fGAliceFileName = fileName;
339   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
340     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
341                                                 fGAliceFileName);
342     fGAliceFileName = absFileName;
343     delete[] absFileName;
344   }
345
346   AliDebug(2, Form("galice file name set to %s", fileName));
347 }
348
349 //_____________________________________________________________________________
350 void AliSimulation::SetEventsPerFile(const char* detector, const char* type, 
351                                      Int_t nEvents)
352 {
353 // set the number of events per file for the given detector and data type
354 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
355
356   TNamed* obj = new TNamed(detector, type);
357   obj->SetUniqueID(nEvents);
358   fEventsPerFile.Add(obj);
359 }
360
361 //_____________________________________________________________________________
362 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
363 {
364   // Read the alignment objects from CDB.
365   // Each detector is supposed to have the
366   // alignment objects in DET/Align/Data CDB path.
367   // All the detector objects are then collected,
368   // sorted by geometry level (starting from ALIC) and
369   // then applied to the TGeo geometry.
370   // Finally an overlaps check is performed.
371
372   if (!gGeoManager || !gGeoManager->IsClosed()) {
373     AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
374     return kFALSE;
375   }  
376   Bool_t delRunLoader = kFALSE;
377   if (!runLoader) {
378     runLoader = LoadRun("READ");
379     if (!runLoader) return kFALSE;
380     delRunLoader = kTRUE;
381   }
382
383   // Export ideal geometry 
384   gGeoManager->Export("geometry.root");
385
386   // Load alignment data from CDB and apply to geometry through AliGeomManager
387   if(fLoadAlignFromCDB){
388     
389     TString detStr = fLoadAlObjsListOfDets;
390     TString loadAlObjsListOfDets = "";
391     
392     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
393     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
394       AliModule* det = (AliModule*) detArray->At(iDet);
395       if (!det || !det->IsActive()) continue;
396       if (IsSelected(det->GetName(), detStr)) {
397         //add det to list of dets to be aligned from CDB
398         loadAlObjsListOfDets += det->GetName();
399         loadAlObjsListOfDets += " ";
400       }
401     } // end loop over detectors
402     (AliGeomManager::Instance())->ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
403   }else{
404     // Check if the array with alignment objects was
405     // provided by the user. If yes, apply the objects
406     // to the present TGeo geometry
407     if (fAlignObjArray) {
408       if ((AliGeomManager::Instance())->ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
409         AliError("The misalignment of one or more volumes failed!"
410                  "Compare the list of simulated detectors and the list of detector alignment data!");
411         if (delRunLoader) delete runLoader;
412         return kFALSE;
413       }
414     }
415   }
416
417   // Update the internal geometry of modules (ITS needs it)
418   TString detStr = fLoadAlObjsListOfDets;
419   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
420   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
421
422     AliModule* det = (AliModule*) detArray->At(iDet);
423     if (!det || !det->IsActive()) continue;
424     if (IsSelected(det->GetName(), detStr)) {
425       det->UpdateInternalGeometry();
426     }
427   } // end loop over detectors
428
429
430   if (delRunLoader) delete runLoader;
431
432   // Update the TGeoPhysicalNodes
433   gGeoManager->RefreshPhysicalNodes();
434
435   return kTRUE;
436 }
437
438
439 //_____________________________________________________________________________
440 Bool_t AliSimulation::SetRunNumber()
441 {
442   // Set the CDB manager run number
443   // The run number is retrieved from gAlice
444
445   if(AliCDBManager::Instance()->GetRun() < 0) {
446     AliRunLoader* runLoader = LoadRun("READ");
447     if (!runLoader) return kFALSE;
448     else {
449       AliCDBManager::Instance()->SetRun(runLoader->GetAliRun()->GetRunNumber());
450       AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
451       delete runLoader;
452     }
453   }
454   return kTRUE;
455 }
456
457 //_____________________________________________________________________________
458 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
459 {
460 // add a file with background events for merging
461
462   TObjString* fileNameStr = new TObjString(fileName);
463   fileNameStr->SetUniqueID(nSignalPerBkgrd);
464   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
465   fBkgrdFileNames->Add(fileNameStr);
466 }
467
468 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
469 {
470 // add a file with background events for embeddin
471   MergeWith(fileName, nSignalPerBkgrd);
472   fEmbeddingFlag = kTRUE;
473 }
474
475 //_____________________________________________________________________________
476 Bool_t AliSimulation::Run(Int_t nEvents)
477 {
478 // run the generation, simulation and digitization
479
480   InitCDBStorage();
481
482   if (nEvents > 0) fNEvents = nEvents;
483
484   // generation and simulation -> hits
485   if (fRunGeneration) {
486     if (!RunSimulation()) if (fStopOnError) return kFALSE;
487   }
488
489   // Set run number in CDBManager (if it is not already set in RunSimulation)
490   if (!SetRunNumber()) if (fStopOnError) return kFALSE;
491
492   // If RunSimulation was not called, load the geometry and misalign it
493   if (!gGeoManager) {
494     TGeoManager::Import("geometry.root");
495     if (!gGeoManager) if (fStopOnError) return kFALSE;
496     // Initialize the geometry manager (if not already done)
497     if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
498   }
499
500   // hits -> summable digits
501   if (!fMakeSDigits.IsNull()) {
502     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
503   }
504
505   // summable digits -> digits
506   if (!fMakeDigits.IsNull()) {
507     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
508       if (fStopOnError) return kFALSE;
509     }
510   }
511
512   // hits -> digits
513   if (!fMakeDigitsFromHits.IsNull()) {
514     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
515       AliWarning(Form("Merging and direct creation of digits from hits " 
516                  "was selected for some detectors. "
517                  "No merging will be done for the following detectors: %s",
518                  fMakeDigitsFromHits.Data()));
519     }
520     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
521       if (fStopOnError) return kFALSE;
522     }
523   }
524
525   // digits -> trigger
526   if (!RunTrigger(fMakeTrigger)) {
527     if (fStopOnError) return kFALSE;
528   }
529
530   // digits -> raw data
531   if (!fWriteRawData.IsNull()) {
532     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
533                       fDeleteIntermediateFiles)) {
534       if (fStopOnError) return kFALSE;
535     }
536   }
537
538   return kTRUE;
539 }
540
541 //_____________________________________________________________________________
542 Bool_t AliSimulation::RunTrigger(const char* descriptors)
543 {
544   // run the trigger
545
546    TStopwatch stopwatch;
547    stopwatch.Start();
548
549    AliRunLoader* runLoader = LoadRun("READ");
550    if (!runLoader) return kFALSE;
551    TString des = descriptors;
552
553    if (des.IsNull()) {
554      if (gAlice->GetTriggerDescriptor() != "") {
555        des = gAlice->GetTriggerDescriptor();
556      }
557      else {
558        AliWarning("No trigger descriptor is specified. Skipping the trigger simulation...");
559        return kTRUE;
560      }
561    }
562
563    runLoader->MakeTree( "GG" );
564    AliCentralTrigger* aCTP = runLoader->GetTrigger();
565   // Load Descriptors
566    aCTP->LoadDescriptor( des );
567
568   // digits -> trigger
569    if( !aCTP->RunTrigger( runLoader ) ) {
570       if (fStopOnError) {
571     //  delete aCTP;
572          return kFALSE;
573       }
574    }
575
576    AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
577            stopwatch.RealTime(),stopwatch.CpuTime()));
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   TStopwatch stopwatch;
602   stopwatch.Start();
603
604   if (!gAlice) {
605     AliError("no gAlice object. Restart aliroot and try again.");
606     return kFALSE;
607   }
608   if (gAlice->Modules()->GetEntries() > 0) {
609     AliError("gAlice was already run. Restart aliroot and try again.");
610     return kFALSE;
611   }
612
613   AliInfo(Form("initializing gAlice with config file %s",
614           fConfigFileName.Data()));
615   StdoutToAliInfo(StderrToAliError(
616     gAlice->Init(fConfigFileName.Data());
617   ););
618
619   // Get the trigger descriptor string
620   // Either from AliSimulation or from
621   // gAlice
622   if (fMakeTrigger.IsNull()) {
623     if (gAlice->GetTriggerDescriptor() != "")
624       fMakeTrigger = gAlice->GetTriggerDescriptor();
625   }
626   else
627     gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
628
629   // Set run number in CDBManager
630   AliCDBManager::Instance()->SetRun(gAlice->GetRunNumber());
631   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
632
633   AliRunLoader* runLoader = gAlice->GetRunLoader();
634   if (!runLoader) {
635              AliError(Form("gAlice has no run loader object. "
636                              "Check your config file: %s", fConfigFileName.Data()));
637              return kFALSE;
638   }
639   SetGAliceFile(runLoader->GetFileName());
640  
641   // Misalign geometry
642 #if ROOT_VERSION_CODE < 331527
643   MisalignGeometry(runLoader);
644 #endif
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 }