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