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