]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
Several fixes related to the generation of DATE and ROOT raw data files. 1. Root...
[u/mrichter/AliRoot.git] / STEER / AliSimulation.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // class for running generation, simulation and digitization                 //
21 //                                                                           //
22 // Hits, sdigits and digits are created for all detectors by typing:         //
23 //                                                                           //
24 //   AliSimulation sim;                                                      //
25 //   sim.Run();                                                              //
26 //                                                                           //
27 // The Run method returns kTRUE in case of successful execution.             //
28 // The number of events can be given as argument to the Run method or it     //
29 // can be set by                                                             //
30 //                                                                           //
31 //   sim.SetNumberOfEvents(n);                                               //
32 //                                                                           //
33 // The name of the configuration file can be passed as argument to the       //
34 // AliSimulation constructor or can be specified by                          //
35 //                                                                           //
36 //   sim.SetConfigFile("...");                                               //
37 //                                                                           //
38 // The generation of particles and the simulation of detector hits can be    //
39 // switched on or off by                                                     //
40 //                                                                           //
41 //   sim.SetRunGeneration(kTRUE);   // generation of primary particles       //
42 //   sim.SetRunSimulation(kFALSE);  // but no tracking                       //
43 //                                                                           //
44 // For which detectors sdigits and digits will be created, can be steered    //
45 // by                                                                        //
46 //                                                                           //
47 //   sim.SetMakeSDigits("ALL");     // make sdigits for all detectors        //
48 //   sim.SetMakeDigits("ITS TPC");  // make digits only for ITS and TPC      //
49 //                                                                           //
50 // The argument is a (case sensitive) string with the names of the           //
51 // detectors separated by a space. An empty string ("") can be used to       //
52 // disable the creation of sdigits or digits. The special string "ALL"       //
53 // selects all available detectors. This is the default.                     //
54 //                                                                           //
55 // The creation of digits from hits instead of from sdigits can be selected  //
56 // by                                                                        //
57 //                                                                           //
58 //   sim.SetMakeDigitsFromHits("TRD");                                       //
59 //                                                                           //
60 // The argument is again a string with the selected detectors. Be aware that //
61 // this feature is not available for all detectors and that merging is not   //
62 // possible, when digits are created directly from hits.                     //
63 //                                                                           //
64 // Background events can be merged by calling                                //
65 //                                                                           //
66 //   sim.MergeWith("background/galice.root", 2);                             //
67 //                                                                           //
68 // The first argument is the file name of the background galice file. The    //
69 // second argument is the number of signal events per background event.      //
70 // By default this number is calculated from the number of available         //
71 // background events. MergeWith can be called several times to merge more    //
72 // than two event streams. It is assumed that the sdigits were already       //
73 // produced for the background events.                                       //
74 //                                                                           //
75 // The output of raw data can be switched on by calling                      //
76 //                                                                           //
77 //   sim.SetWriteRawData("MUON");   // write raw data for MUON               //
78 //                                                                           //
79 // The default output format of the raw data are DDL files. They are         //
80 // converted to a DATE file, if a file name is given as second argument.     //
81 // For this conversion the program "dateStream" is required. If the file     //
82 // name has the extension ".root", the DATE file is converted to a root      //
83 // file. The program "alimdc" is used for this purpose. For the conversion   //
84 // to DATE and root format the two conversion programs have to be installed. //
85 // Only the raw data in the final format is kept if the third argument is    //
86 // kTRUE.                                                                    //
87 //                                                                           //
88 // The methods RunSimulation, RunSDigitization, RunDigitization,             //
89 // RunHitsDigitization and WriteRawData can be used to run only parts of     //
90 // the full simulation chain. The creation of raw data DDL files and their   //
91 // conversion to the DATE or root format can be run directly by calling      //
92 // the methods WriteRawFiles, ConvertRawFilesToDate and ConvertDateToRoot.   //
93 //                                                                           //
94 // The default number of events per file, which is usually set in the        //
95 // config file, can be changed for individual detectors and data types       //
96 // by calling                                                                //
97 //                                                                           //
98 //   sim.SetEventsPerFile("PHOS", "Reconstructed Points", 3);                //
99 //                                                                           //
100 // The first argument is the detector, the second one the data type and the  //
101 // last one the number of events per file. Valid data types are "Hits",      //
102 // "Summable Digits", "Digits", "Reconstructed Points" and "Tracks".         //
103 // The number of events per file has to be set before the simulation of      //
104 // hits. Otherwise it has no effect.                                         //
105 //                                                                           //
106 ///////////////////////////////////////////////////////////////////////////////
107
108 #include <TGeoManager.h>
109 #include <TObjString.h>
110 #include <TStopwatch.h>
111 #include <TSystem.h>
112 #include <TFile.h>
113
114 #include "AliCDBStorage.h"
115 #include "AliCDBEntry.h"
116 #include "AliCDBManager.h"
117 #include "AliAlignObj.h"
118 #include "AliCentralTrigger.h"
119 #include "AliDAQ.h"
120 #include "AliDigitizer.h"
121 #include "AliGenerator.h"
122 #include "AliLog.h"
123 #include "AliModule.h"
124 #include "AliRun.h"
125 #include "AliRunDigitizer.h"
126 #include "AliRunLoader.h"
127 #include "AliSimulation.h"
128 #include "AliVertexGenFile.h"
129 #include "AliCentralTrigger.h"
130 #include "AliCTPRawData.h"
131 #include "AliRawReaderFile.h"
132 #include "AliESD.h"
133 #include "AliHeader.h"
134 #include "AliGenEventHeader.h"
135
136 ClassImp(AliSimulation)
137
138
139 //_____________________________________________________________________________
140 AliSimulation::AliSimulation(const char* configFileName, const char* cdbUri,
141                              const char* name, const char* title) :
142   TNamed(name, title),
143
144   fRunGeneration(kTRUE),
145   fRunSimulation(kTRUE),
146   fLoadAlignFromCDB(kTRUE),
147   fLoadAlignData("ALL"),
148   fMakeSDigits("ALL"),
149   fMakeDigits("ALL"),
150   fMakeTrigger(""),
151   fMakeDigitsFromHits(""),
152   fWriteRawData(""),
153   fRawDataFileName(""),
154   fDeleteIntermediateFiles(kFALSE),
155   fStopOnError(kFALSE),
156
157   fNEvents(1),
158   fConfigFileName(configFileName),
159   fGAliceFileName("galice.root"),
160   fEventsPerFile(),
161   fBkgrdFileNames(NULL),
162   fAlignObjArray(NULL),
163   fUseBkgrdVertex(kTRUE),
164   fRegionOfInterest(kFALSE),
165   fCDBUri(cdbUri),
166   fSpecCDBUri(),
167   fEmbeddingFlag(kFALSE)
168 {
169 // create simulation object with default parameters
170
171   SetGAliceFile("galice.root");
172 }
173
174 //_____________________________________________________________________________
175 AliSimulation::AliSimulation(const AliSimulation& sim) :
176   TNamed(sim),
177
178   fRunGeneration(sim.fRunGeneration),
179   fRunSimulation(sim.fRunSimulation),
180   fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
181   fLoadAlignData(sim.fLoadAlignData),
182   fMakeSDigits(sim.fMakeSDigits),
183   fMakeDigits(sim.fMakeDigits),
184   fMakeTrigger(sim.fMakeTrigger),
185   fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
186   fWriteRawData(sim.fWriteRawData),
187   fRawDataFileName(""),
188   fDeleteIntermediateFiles(kFALSE),
189   fStopOnError(sim.fStopOnError),
190
191   fNEvents(sim.fNEvents),
192   fConfigFileName(sim.fConfigFileName),
193   fGAliceFileName(sim.fGAliceFileName),
194   fEventsPerFile(),
195   fBkgrdFileNames(NULL),
196   fAlignObjArray(NULL),
197   fUseBkgrdVertex(sim.fUseBkgrdVertex),
198   fRegionOfInterest(sim.fRegionOfInterest),
199   fCDBUri(sim.fCDBUri),
200   fSpecCDBUri(),
201   fEmbeddingFlag(sim.fEmbeddingFlag)
202 {
203 // copy constructor
204
205   for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
206     if (!sim.fEventsPerFile[i]) continue;
207     fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
208   }
209
210   fBkgrdFileNames = new TObjArray;
211   for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
212     if (!sim.fBkgrdFileNames->At(i)) continue;
213     fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
214   }
215
216   for (Int_t i = 0; i < sim.fSpecCDBUri.GetEntriesFast(); i++) {
217     if (sim.fSpecCDBUri[i]) fSpecCDBUri.Add(sim.fSpecCDBUri[i]->Clone());
218   }
219
220 }
221
222 //_____________________________________________________________________________
223 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
224 {
225 // assignment operator
226
227   this->~AliSimulation();
228   new(this) AliSimulation(sim);
229   return *this;
230 }
231
232 //_____________________________________________________________________________
233 AliSimulation::~AliSimulation()
234 {
235 // clean up
236
237   fEventsPerFile.Delete();
238 //  if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
239 //  delete fAlignObjArray; fAlignObjArray=0;
240
241   if (fBkgrdFileNames) {
242     fBkgrdFileNames->Delete();
243     delete fBkgrdFileNames;
244   }
245
246   fSpecCDBUri.Delete();
247 }
248
249
250 //_____________________________________________________________________________
251 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
252 {
253 // set the number of events for one run
254
255   fNEvents = nEvents;
256 }
257
258 //_____________________________________________________________________________
259 void AliSimulation::InitCDBStorage()
260 {
261 // activate a default CDB storage
262 // First check if we have any CDB storage set, because it is used 
263 // to retrieve the calibration and alignment constants
264
265   AliCDBManager* man = AliCDBManager::Instance();
266   if (man->IsDefaultStorageSet())
267   {
268     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
269     AliWarning("Default CDB storage has been already set !");
270     AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
271     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
272     fCDBUri = "";
273   }
274   else {
275     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
276     AliWarning(Form("Default CDB storage is set to: %s",fCDBUri.Data()));
277     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
278     man->SetDefaultStorage(fCDBUri);
279   }
280
281   // Now activate the detector specific CDB storage locations
282   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
283     TObject* obj = fSpecCDBUri[i];
284     if (!obj) continue;
285     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
286     AliWarning(Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
287     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
288     man->SetSpecificStorage(obj->GetName(),obj->GetTitle());
289   }
290   man->Print();
291 }
292
293 //_____________________________________________________________________________
294 void AliSimulation::SetDefaultStorage(const char* uri) {
295 // Store the desired default CDB storage location
296 // Activate it later within the Run() method
297
298   fCDBUri = uri;
299
300 }
301
302 //_____________________________________________________________________________
303 void AliSimulation::SetSpecificStorage(const char* detName, const char* uri) {
304 // Store a detector-specific CDB storage location
305 // Activate it later within the Run() method
306
307   TObject* obj = fSpecCDBUri.FindObject(detName);
308   if (obj) fSpecCDBUri.Remove(obj);
309   fSpecCDBUri.Add(new TNamed(detName, uri));
310
311 }
312
313 //_____________________________________________________________________________
314 void AliSimulation::SetConfigFile(const char* fileName)
315 {
316 // set the name of the config file
317
318   fConfigFileName = fileName;
319 }
320
321 //_____________________________________________________________________________
322 void AliSimulation::SetGAliceFile(const char* fileName)
323 {
324 // set the name of the galice file
325 // the path is converted to an absolute one if it is relative
326
327   fGAliceFileName = fileName;
328   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
329     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
330                                                 fGAliceFileName);
331     fGAliceFileName = absFileName;
332     delete[] absFileName;
333   }
334
335   AliDebug(2, Form("galice file name set to %s", fileName));
336 }
337
338 //_____________________________________________________________________________
339 void AliSimulation::SetEventsPerFile(const char* detector, const char* type, 
340                                      Int_t nEvents)
341 {
342 // set the number of events per file for the given detector and data type
343 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
344
345   TNamed* obj = new TNamed(detector, type);
346   obj->SetUniqueID(nEvents);
347   fEventsPerFile.Add(obj);
348 }
349
350 //_____________________________________________________________________________
351 Bool_t AliSimulation::ApplyAlignObjsToGeom(TObjArray* alObjArray)
352 {
353   // Read collection of alignment objects (AliAlignObj derived) saved
354   // in the TClonesArray ClArrayName and apply them to the geometry
355   // manager singleton.
356   //
357   alObjArray->Sort();
358   Int_t nvols = alObjArray->GetEntriesFast();
359
360   Bool_t flag = kTRUE;
361
362   for(Int_t j=0; j<nvols; j++)
363     {
364       AliAlignObj* alobj = (AliAlignObj*) alObjArray->UncheckedAt(j);
365       if (alobj->ApplyToGeometry() == kFALSE) flag = kFALSE;
366     }
367
368   if (AliDebugLevelClass() >= 1) {
369     gGeoManager->GetTopNode()->CheckOverlaps(20);
370     TObjArray* ovexlist = gGeoManager->GetListOfOverlaps();
371     if(ovexlist->GetEntriesFast()){  
372       AliError("The application of alignment objects to the geometry caused huge overlaps/extrusions!");
373    }
374   }
375
376   return flag;
377
378 }
379
380 //_____________________________________________________________________________
381 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* fileName, const char* clArrayName)
382 {
383   // read collection of alignment objects (AliAlignObj derived) saved
384   // in the TClonesArray ClArrayName in the file fileName and apply
385   // them to the TGeo geometry passed as argument
386   //
387
388   TFile* inFile = TFile::Open(fileName,"READ");
389   if (!inFile || !inFile->IsOpen()) {
390     AliErrorClass(Form("Could not open file %s !",fileName));
391     return kFALSE;
392   }
393
394   TClonesArray* alObjArray = ((TClonesArray*) inFile->Get(clArrayName));
395   inFile->Close();
396   if (!alObjArray) {
397     AliErrorClass(Form("Could not get array (%s) from file (%s) !",clArrayName,fileName));
398     return kFALSE;
399   }
400
401   return ApplyAlignObjsToGeom(alObjArray);
402
403 }
404
405 //_____________________________________________________________________________
406 Bool_t AliSimulation::ApplyAlignObjsToGeom(AliCDBParam* param, AliCDBId& Id)
407 {
408   // read collection of alignment objects (AliAlignObj derived) saved
409   // in the TClonesArray ClArrayName in the AliCDBEntry identified by
410   // param (to get the AliCDBStorage) and Id; apply the alignment objects
411   // to the TGeo geometry passed as argument
412   //
413
414   AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage(param);
415   AliCDBEntry* entry = storage->Get(Id);
416   TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
417
418   return ApplyAlignObjsToGeom(AlObjArray);
419
420 }
421
422 //_____________________________________________________________________________
423 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* uri, const char* path, Int_t runnum, Int_t version, Int_t sversion)
424 {
425   // read collection of alignment objects (AliAlignObj derived) saved
426   // in the TClonesArray ClArrayName in the AliCDBEntry identified by
427   // param (to get the AliCDBStorage) and Id; apply the alignment objects
428   // to the TGeo geometry passed as argument
429   //
430
431   AliCDBParam* param = AliCDBManager::Instance()->CreateParameter(uri);
432   AliCDBId id(path, runnum, runnum, version, sversion);
433
434   return ApplyAlignObjsToGeom(param, id);
435
436 }
437
438 //_____________________________________________________________________________
439 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* detName, Int_t runnum, Int_t version, Int_t sversion)
440 {
441   // read collection of alignment objects (AliAlignObj derived) saved
442   // in the TClonesArray ClArrayName in the AliCDBEntry identified by
443   // param (to get the AliCDBStorage) and Id; apply the alignment objects
444   // to the TGeo geometry passed as argument
445   //
446
447   AliCDBPath path(detName,"Align","Data");
448   AliCDBEntry* entry = AliCDBManager::Instance()->Get(path.GetPath(),runnum,version,sversion);
449
450   if(!entry) return kFALSE;
451   TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
452
453   return ApplyAlignObjsToGeom(AlObjArray);
454 }
455
456 //_____________________________________________________________________________
457 Bool_t AliSimulation::SetAlignObjArraySingleDet(const char* detName)
458 {
459   // Fills array of single detector's alignable objects from CDB
460   
461   AliDebug(2, Form("Loading alignment data for detector: %s",detName));
462   
463   AliCDBEntry *entry;
464         
465   AliCDBPath path(detName,"Align","Data");
466         
467   entry=AliCDBManager::Instance()->Get(path.GetPath());
468   if(!entry){ 
469         AliDebug(2,Form("Couldn't load alignment data for detector %s",detName));
470         return kFALSE;
471   }
472   entry->SetOwner(1);
473   TClonesArray *alignArray = (TClonesArray*) entry->GetObject();        
474   alignArray->SetOwner(0);
475   AliDebug(2,Form("Found %d alignment objects for %s",
476                         alignArray->GetEntries(),detName));
477
478   AliAlignObj *alignObj=0;
479   TIter iter(alignArray);
480         
481   // loop over align objects in detector
482   while( ( alignObj=(AliAlignObj *) iter.Next() ) ){
483         fAlignObjArray->Add(alignObj);
484   }
485   // delete entry --- Don't delete, it is cached!
486         
487   AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
488   return kTRUE;
489
490 }
491
492 //_____________________________________________________________________________
493 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
494 {
495   // Read the alignment objects from CDB.
496   // Each detector is supposed to have the
497   // alignment objects in DET/Align/Data CDB path.
498   // All the detector objects are then collected,
499   // sorted by geometry level (starting from ALIC) and
500   // then applied to the TGeo geometry.
501   // Finally an overlaps check is performed.
502
503   Bool_t delRunLoader = kFALSE;
504   if (!runLoader) {
505     runLoader = LoadRun("READ");
506     if (!runLoader) return kFALSE;
507     delRunLoader = kTRUE;
508   }
509
510   // Load alignment data from CDB and fill fAlignObjArray 
511   if(fLoadAlignFromCDB){
512         if(!fAlignObjArray) fAlignObjArray = new TObjArray();
513         
514         //fAlignObjArray->RemoveAll(); 
515         fAlignObjArray->Clear();        
516         fAlignObjArray->SetOwner(0);
517  
518         TString detStr = fLoadAlignData;
519         TString dataNotLoaded="";
520         TString dataLoaded="";
521   
522         TObjArray* detArray = runLoader->GetAliRun()->Detectors();
523         for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
524                 AliModule* det = (AliModule*) detArray->At(iDet);
525                 if (!det) continue;
526                 if (IsSelected(det->GetName(), detStr)) {
527                         if(!SetAlignObjArraySingleDet(det->GetName())){
528                                 dataNotLoaded += det->GetName();
529                                 dataNotLoaded += " ";
530                         } else {
531                                 dataLoaded += det->GetName();
532                                 dataLoaded += " ";
533                         }
534                 }
535         } // end loop over detectors
536   
537         if ((detStr.CompareTo("ALL") == 0)) detStr = "";
538         dataNotLoaded += detStr;
539         AliInfo(Form("Alignment data loaded for: %s",
540                           dataLoaded.Data()));
541         AliInfo(Form("Didn't/couldn't load alignment data for: %s",
542                           dataNotLoaded.Data()));
543   } // fLoadAlignFromCDB flag
544  
545   // Check if the array with alignment objects was
546   // provided by the user. If yes, apply the objects
547   // to the present TGeo geometry
548   if (fAlignObjArray) {
549     if (gGeoManager && gGeoManager->IsClosed()) {
550       if (ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
551         AliError("The misalignment of one or more volumes failed!"
552                  "Compare the list of simulated detectors and the list of detector alignment data!");
553         if (delRunLoader) delete runLoader;
554         return kFALSE;
555       }
556     }
557     else {
558       AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
559       if (delRunLoader) delete runLoader;
560       return kFALSE;
561     }
562   }
563
564
565   // Update the internal geometry of modules (ITS needs it)
566   TString detStr = fLoadAlignData;
567   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
568   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
569
570     AliModule* det = (AliModule*) detArray->At(iDet);
571     if (!det || !det->IsActive()) continue;
572     if (IsSelected(det->GetName(), detStr)) {
573       det->UpdateInternalGeometry();
574     }
575   } // end loop over detectors
576
577
578   if (delRunLoader) delete runLoader;
579
580   return kTRUE;
581 }
582
583
584 //_____________________________________________________________________________
585 Bool_t AliSimulation::SetRunNumber()
586 {
587   // Set the CDB manager run number
588   // The run number is retrieved from gAlice
589
590   if(AliCDBManager::Instance()->GetRun() < 0) {
591     AliRunLoader* runLoader = LoadRun("READ");
592     if (!runLoader) return kFALSE;
593     else {
594       AliCDBManager::Instance()->SetRun(runLoader->GetAliRun()->GetRunNumber());
595       AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
596       delete runLoader;
597     }
598   }
599   return kTRUE;
600 }
601
602 //_____________________________________________________________________________
603 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
604 {
605 // add a file with background events for merging
606
607   TObjString* fileNameStr = new TObjString(fileName);
608   fileNameStr->SetUniqueID(nSignalPerBkgrd);
609   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
610   fBkgrdFileNames->Add(fileNameStr);
611 }
612
613 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
614 {
615 // add a file with background events for embeddin
616   MergeWith(fileName, nSignalPerBkgrd);
617   fEmbeddingFlag = kTRUE;
618 }
619
620 //_____________________________________________________________________________
621 Bool_t AliSimulation::Run(Int_t nEvents)
622 {
623 // run the generation, simulation and digitization
624
625   InitCDBStorage();
626
627   if (nEvents > 0) fNEvents = nEvents;
628
629   // generation and simulation -> hits
630   if (fRunGeneration) {
631     if (!RunSimulation()) if (fStopOnError) return kFALSE;
632   }
633
634   // Set run number in CDBManager (if it is not already set in RunSimulation)
635   if (!SetRunNumber()) if (fStopOnError) return kFALSE;
636
637   // Load and misalign the geometry
638   if (!gGeoManager) {
639     TGeoManager::Import("geometry.root");
640     if (!gGeoManager) if (fStopOnError) return kFALSE;
641     if (!MisalignGeometry()) if (fStopOnError) return kFALSE;
642   }
643   
644   // hits -> summable digits
645   if (!fMakeSDigits.IsNull()) {
646     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
647   }
648
649   // summable digits -> digits
650   if (!fMakeDigits.IsNull()) {
651     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
652       if (fStopOnError) return kFALSE;
653     }
654   }
655
656   // hits -> digits
657   if (!fMakeDigitsFromHits.IsNull()) {
658     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
659       AliWarning(Form("Merging and direct creation of digits from hits " 
660                  "was selected for some detectors. "
661                  "No merging will be done for the following detectors: %s",
662                  fMakeDigitsFromHits.Data()));
663     }
664     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
665       if (fStopOnError) return kFALSE;
666     }
667   }
668
669   // digits -> trigger
670   if (!RunTrigger(fMakeTrigger)) {
671     if (fStopOnError) return kFALSE;
672   }
673
674   // digits -> raw data
675   if (!fWriteRawData.IsNull()) {
676     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
677                       fDeleteIntermediateFiles)) {
678       if (fStopOnError) return kFALSE;
679     }
680   }
681
682   return kTRUE;
683 }
684
685 //_____________________________________________________________________________
686 Bool_t AliSimulation::RunTrigger(const char* descriptors)
687 {
688   // run the trigger
689
690    TStopwatch stopwatch;
691    stopwatch.Start();
692
693    AliRunLoader* runLoader = LoadRun("READ");
694    if (!runLoader) return kFALSE;
695    TString des = descriptors;
696
697    if (des.IsNull()) {
698      if (gAlice->GetTriggerDescriptor() != "") {
699        des = gAlice->GetTriggerDescriptor();
700      }
701      else {
702        AliWarning("No trigger descriptor is specified. Skipping the trigger simulation...");
703        return kTRUE;
704      }
705    }
706
707    runLoader->MakeTree( "CT" );
708    AliCentralTrigger* aCTP = runLoader->GetTrigger();
709   // Load Descriptors
710    aCTP->LoadDescriptor( des );
711
712   // digits -> trigger
713    if( !aCTP->RunTrigger( runLoader ) ) {
714       if (fStopOnError) {
715     //  delete aCTP;
716          return kFALSE;
717       }
718    }
719
720    AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
721            stopwatch.RealTime(),stopwatch.CpuTime()));
722
723    delete runLoader;
724
725    return kTRUE;
726 }
727
728 //_____________________________________________________________________________
729 Bool_t AliSimulation::WriteTriggerRawData()
730 {
731   // Writes the CTP (trigger) DDL raw data
732   // Details of the format are given in the
733   // trigger TDR - pages 134 and 135.
734   AliCTPRawData writer;
735   writer.RawData();
736
737   return kTRUE;
738 }
739
740 //_____________________________________________________________________________
741 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
742 {
743 // run the generation and simulation
744
745   TStopwatch stopwatch;
746   stopwatch.Start();
747
748   if (!gAlice) {
749     AliError("no gAlice object. Restart aliroot and try again.");
750     return kFALSE;
751   }
752   if (gAlice->Modules()->GetEntries() > 0) {
753     AliError("gAlice was already run. Restart aliroot and try again.");
754     return kFALSE;
755   }
756
757   AliInfo(Form("initializing gAlice with config file %s",
758           fConfigFileName.Data()));
759   StdoutToAliInfo(StderrToAliError(
760     gAlice->Init(fConfigFileName.Data());
761   ););
762
763   // Get the trigger descriptor string
764   // Either from AliSimulation or from
765   // gAlice
766   if (fMakeTrigger.IsNull()) {
767     if (gAlice->GetTriggerDescriptor() != "")
768       fMakeTrigger = gAlice->GetTriggerDescriptor();
769   }
770   else
771     gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
772
773   // Set run number in CDBManager
774   AliCDBManager::Instance()->SetRun(gAlice->GetRunNumber());
775   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
776
777   AliRunLoader* runLoader = gAlice->GetRunLoader();
778   if (!runLoader) {
779              AliError(Form("gAlice has no run loader object. "
780                              "Check your config file: %s", fConfigFileName.Data()));
781              return kFALSE;
782   }
783   SetGAliceFile(runLoader->GetFileName());
784  
785   // Export ideal geometry 
786   if (gGeoManager) gGeoManager->Export("geometry.root");
787
788   // Misalign geometry
789 //   if (!MisalignGeometry(runLoader)) {
790 //     delete runLoader;
791 //     return kFALSE;
792 //   }
793   MisalignGeometry(runLoader);
794
795   // Export (mis)aligned geometry 
796   if (gGeoManager) gGeoManager->Export("misaligned_geometry.root");
797
798 //   AliRunLoader* runLoader = gAlice->GetRunLoader();
799 //   if (!runLoader) {
800 //     AliError(Form("gAlice has no run loader object. "
801 //                   "Check your config file: %s", fConfigFileName.Data()));
802 //     return kFALSE;
803 //   }
804 //   SetGAliceFile(runLoader->GetFileName());
805
806   if (!gAlice->Generator()) {
807     AliError(Form("gAlice has no generator object. "
808                   "Check your config file: %s", fConfigFileName.Data()));
809     return kFALSE;
810   }
811   if (nEvents <= 0) nEvents = fNEvents;
812
813   // get vertex from background file in case of merging
814   if (fUseBkgrdVertex &&
815       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
816     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
817     const char* fileName = ((TObjString*)
818                             (fBkgrdFileNames->At(0)))->GetName();
819     AliInfo(Form("The vertex will be taken from the background "
820                  "file %s with nSignalPerBackground = %d", 
821                  fileName, signalPerBkgrd));
822     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
823     gAlice->Generator()->SetVertexGenerator(vtxGen);
824   }
825
826   if (!fRunSimulation) {
827     gAlice->Generator()->SetTrackingFlag(0);
828   }
829
830   // set the number of events per file for given detectors and data types
831   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
832     if (!fEventsPerFile[i]) continue;
833     const char* detName = fEventsPerFile[i]->GetName();
834     const char* typeName = fEventsPerFile[i]->GetTitle();
835     TString loaderName(detName);
836     loaderName += "Loader";
837     AliLoader* loader = runLoader->GetLoader(loaderName);
838     if (!loader) {
839       AliError(Form("RunSimulation", "no loader for %s found\n"
840                     "Number of events per file not set for %s %s", 
841                     detName, typeName, detName));
842       continue;
843     }
844     AliDataLoader* dataLoader = 
845       loader->GetDataLoader(typeName);
846     if (!dataLoader) {
847       AliError(Form("no data loader for %s found\n"
848                     "Number of events per file not set for %s %s", 
849                     typeName, detName, typeName));
850       continue;
851     }
852     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
853     AliDebug(1, Form("number of events per file set to %d for %s %s",
854                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
855   }
856
857   AliInfo("running gAlice");
858   StdoutToAliInfo(StderrToAliError(
859     gAlice->Run(nEvents);
860   ););
861
862   delete runLoader;
863
864   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
865                stopwatch.RealTime(),stopwatch.CpuTime()));
866
867   return kTRUE;
868 }
869
870 //_____________________________________________________________________________
871 Bool_t AliSimulation::RunSDigitization(const char* detectors)
872 {
873 // run the digitization and produce summable digits
874
875   TStopwatch stopwatch;
876   stopwatch.Start();
877
878   AliRunLoader* runLoader = LoadRun();
879   if (!runLoader) return kFALSE;
880
881   TString detStr = detectors;
882   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
883   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
884     AliModule* det = (AliModule*) detArray->At(iDet);
885     if (!det || !det->IsActive()) continue;
886     if (IsSelected(det->GetName(), detStr)) {
887       AliInfo(Form("creating summable digits for %s", det->GetName()));
888       TStopwatch stopwatchDet;
889       stopwatchDet.Start();
890       det->Hits2SDigits();
891       AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
892            det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
893     }
894   }
895
896   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
897     AliError(Form("the following detectors were not found: %s",
898                   detStr.Data()));
899     if (fStopOnError) return kFALSE;
900   }
901
902   delete runLoader;
903
904   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
905            stopwatch.RealTime(),stopwatch.CpuTime()));
906
907   return kTRUE;
908 }
909
910
911 //_____________________________________________________________________________
912 Bool_t AliSimulation::RunDigitization(const char* detectors, 
913                                       const char* excludeDetectors)
914 {
915 // run the digitization and produce digits from sdigits
916
917   TStopwatch stopwatch;
918   stopwatch.Start();
919
920   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
921   if (gAlice) delete gAlice;
922   gAlice = NULL;
923
924   Int_t nStreams = 1;
925   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
926   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
927   AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
928   // manager->SetEmbeddingFlag(fEmbeddingFlag);
929   manager->SetInputStream(0, fGAliceFileName.Data());
930   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
931     const char* fileName = ((TObjString*)
932                             (fBkgrdFileNames->At(iStream-1)))->GetName();
933     manager->SetInputStream(iStream, fileName);
934   }
935
936   TString detStr = detectors;
937   TString detExcl = excludeDetectors;
938   manager->GetInputStream(0)->ImportgAlice();
939   AliRunLoader* runLoader = 
940     AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
941   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
942   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
943     AliModule* det = (AliModule*) detArray->At(iDet);
944     if (!det || !det->IsActive()) continue;
945     if (IsSelected(det->GetName(), detStr) && 
946         !IsSelected(det->GetName(), detExcl)) {
947       AliDigitizer* digitizer = det->CreateDigitizer(manager);
948       
949       if (!digitizer) {
950         AliError(Form("no digitizer for %s", det->GetName()));
951         if (fStopOnError) return kFALSE;
952       } else {
953         digitizer->SetRegionOfInterest(fRegionOfInterest);
954       }
955     }
956   }
957
958   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
959     AliError(Form("the following detectors were not found: %s", 
960                   detStr.Data()));
961     if (fStopOnError) return kFALSE;
962   }
963
964   if (!manager->GetListOfTasks()->IsEmpty()) {
965     AliInfo("executing digitization");
966     manager->Exec("");
967   }
968
969   delete manager;
970
971   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
972                stopwatch.RealTime(),stopwatch.CpuTime()));
973   
974   return kTRUE;
975 }
976
977 //_____________________________________________________________________________
978 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
979 {
980 // run the digitization and produce digits from hits
981
982   TStopwatch stopwatch;
983   stopwatch.Start();
984
985   AliRunLoader* runLoader = LoadRun("READ");
986   if (!runLoader) return kFALSE;
987
988   TString detStr = detectors;
989   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
990   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
991     AliModule* det = (AliModule*) detArray->At(iDet);
992     if (!det || !det->IsActive()) continue;
993     if (IsSelected(det->GetName(), detStr)) {
994       AliInfo(Form("creating digits from hits for %s", det->GetName()));
995       det->Hits2Digits();
996     }
997   }
998
999   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1000     AliError(Form("the following detectors were not found: %s", 
1001                   detStr.Data()));
1002     if (fStopOnError) return kFALSE;
1003   }
1004
1005   delete runLoader;
1006   //PH Temporary fix to avoid interference with the PHOS loder/getter
1007   //PH The problem has to be solved in more general way 09/06/05
1008
1009   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1010                stopwatch.RealTime(),stopwatch.CpuTime()));
1011
1012   return kTRUE;
1013 }
1014
1015 //_____________________________________________________________________________
1016 Bool_t AliSimulation::WriteRawData(const char* detectors, 
1017                                    const char* fileName,
1018                                    Bool_t deleteIntermediateFiles)
1019 {
1020 // convert the digits to raw data
1021 // First DDL raw data files for the given detectors are created.
1022 // If a file name is given, the DDL files are then converted to a DATE file.
1023 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
1024 // afterwards.
1025 // If the file name has the extension ".root", the DATE file is converted
1026 // to a root file.
1027 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1028
1029   TStopwatch stopwatch;
1030   stopwatch.Start();
1031
1032   if (!WriteRawFiles(detectors)) {
1033     if (fStopOnError) return kFALSE;
1034   }
1035
1036   TString dateFileName(fileName);
1037   if (!dateFileName.IsNull()) {
1038     Bool_t rootOutput = dateFileName.EndsWith(".root");
1039     if (rootOutput) dateFileName += ".date";
1040     if (!ConvertRawFilesToDate(dateFileName)) {
1041       if (fStopOnError) return kFALSE;
1042     }
1043     if (deleteIntermediateFiles) {
1044       AliRunLoader* runLoader = LoadRun("READ");
1045       if (runLoader) for (Int_t iEvent = 0; 
1046                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1047         char command[256];
1048         sprintf(command, "rm -r raw%d", iEvent);
1049         gSystem->Exec(command);
1050       }
1051     }
1052
1053     if (rootOutput) {
1054       if (!ConvertDateToRoot(dateFileName, fileName)) {
1055         if (fStopOnError) return kFALSE;
1056       }
1057       if (deleteIntermediateFiles) {
1058         gSystem->Unlink(dateFileName);
1059       }
1060     }
1061   }
1062
1063   AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1064                stopwatch.RealTime(),stopwatch.CpuTime()));
1065
1066   return kTRUE;
1067 }
1068
1069 //_____________________________________________________________________________
1070 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1071 {
1072 // convert the digits to raw data DDL files
1073
1074   AliRunLoader* runLoader = LoadRun("READ");
1075   if (!runLoader) return kFALSE;
1076
1077   // write raw data to DDL files
1078   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1079     AliInfo(Form("processing event %d", iEvent));
1080     runLoader->GetEvent(iEvent);
1081     TString baseDir = gSystem->WorkingDirectory();
1082     char dirName[256];
1083     sprintf(dirName, "raw%d", iEvent);
1084     gSystem->MakeDirectory(dirName);
1085     if (!gSystem->ChangeDirectory(dirName)) {
1086       AliError(Form("couldn't change to directory %s", dirName));
1087       if (fStopOnError) return kFALSE; else continue;
1088     }
1089
1090     TString detStr = detectors;
1091     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1092     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1093       AliModule* det = (AliModule*) detArray->At(iDet);
1094       if (!det || !det->IsActive()) continue;
1095       if (IsSelected(det->GetName(), detStr)) {
1096         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1097         det->Digits2Raw();
1098       }
1099     }
1100
1101     if (!WriteTriggerRawData())
1102       if (fStopOnError) return kFALSE;
1103
1104     gSystem->ChangeDirectory(baseDir);
1105     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1106       AliError(Form("the following detectors were not found: %s", 
1107                     detStr.Data()));
1108       if (fStopOnError) return kFALSE;
1109     }
1110   }
1111
1112   delete runLoader;
1113   return kTRUE;
1114 }
1115
1116 //_____________________________________________________________________________
1117 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1118 {
1119 // convert raw data DDL files to a DATE file with the program "dateStream"
1120
1121   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1122   if (!path) {
1123     AliError("the program dateStream was not found");
1124     if (fStopOnError) return kFALSE;
1125   } else {
1126     delete[] path;
1127   }
1128
1129   AliRunLoader* runLoader = LoadRun("READ");
1130   if (!runLoader) return kFALSE;
1131
1132   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1133   char command[256];
1134   // Note the option -s. It is used in order to avoid
1135   // the generation of SOR/EOR events.
1136   sprintf(command, "dateStream -s -D -o %s -# %d -C", 
1137           dateFileName, runLoader->GetNumberOfEvents());
1138   FILE* pipe = gSystem->OpenPipe(command, "w");
1139
1140   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1141     fprintf(pipe, "GDC\n");
1142     Float_t ldc = 0;
1143     Int_t prevLDC = -1;
1144
1145     // loop over detectors and DDLs
1146     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1147       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1148
1149         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1150         Int_t ldcID = Int_t(ldc + 0.0001);
1151         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1152
1153         char rawFileName[256];
1154         sprintf(rawFileName, "raw%d/%s", 
1155                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1156
1157         // check existence and size of raw data file
1158         FILE* file = fopen(rawFileName, "rb");
1159         if (!file) continue;
1160         fseek(file, 0, SEEK_END);
1161         unsigned long size = ftell(file);
1162         fclose(file);
1163         if (!size) continue;
1164
1165         if (ldcID != prevLDC) {
1166           fprintf(pipe, " LDC Id %d\n", ldcID);
1167           prevLDC = ldcID;
1168         }
1169         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1170       }
1171     }
1172   }
1173
1174   Int_t result = gSystem->ClosePipe(pipe);
1175
1176   delete runLoader;
1177   return (result == 0);
1178 }
1179
1180 //_____________________________________________________________________________
1181 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1182                                         const char* rootFileName)
1183 {
1184 // convert a DATE file to a root file with the program "alimdc"
1185
1186   // ALIMDC setup
1187   const Int_t kDBSize = 2000000000;
1188   const Int_t kTagDBSize = 1000000000;
1189   const Bool_t kFilter = kFALSE;
1190   const Int_t kCompression = 0;
1191
1192   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1193   if (!path) {
1194     AliError("the program alimdc was not found");
1195     if (fStopOnError) return kFALSE;
1196   } else {
1197     delete[] path;
1198   }
1199
1200   AliInfo(Form("converting DATE file %s to root file %s", 
1201                dateFileName, rootFileName));
1202
1203   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1204   const char* tagDBFS    = "/tmp/mdc1/tags";
1205   const char* runDBFS    = "/tmp/mdc1/meta";
1206
1207   // User defined file system locations
1208   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1209     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1210   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1211     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1212   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1213     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1214   if (gSystem->Getenv("ALIMDC_RUNDB")) 
1215     runDBFS = gSystem->Getenv("ALIMDC_RUNDB");
1216
1217   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1218   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1219   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1220   gSystem->Exec(Form("rm -rf %s",runDBFS));
1221
1222   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1223   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1224   gSystem->Exec(Form("mkdir %s",tagDBFS));
1225   gSystem->Exec(Form("mkdir %s",runDBFS));
1226
1227   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1228                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1229   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1230
1231   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1232   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1233   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1234   gSystem->Exec(Form("rm -rf %s",runDBFS));
1235
1236   return (result == 0);
1237 }
1238
1239
1240 //_____________________________________________________________________________
1241 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1242 {
1243 // delete existing run loaders, open a new one and load gAlice
1244
1245   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1246   AliRunLoader* runLoader = 
1247     AliRunLoader::Open(fGAliceFileName.Data(), 
1248                        AliConfig::GetDefaultEventFolderName(), mode);
1249   if (!runLoader) {
1250     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1251     return NULL;
1252   }
1253   runLoader->LoadgAlice();
1254   gAlice = runLoader->GetAliRun();
1255   if (!gAlice) {
1256     AliError(Form("no gAlice object found in file %s", 
1257                   fGAliceFileName.Data()));
1258     return NULL;
1259   }
1260   return runLoader;
1261 }
1262
1263 //_____________________________________________________________________________
1264 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1265 {
1266 // get or calculate the number of signal events per background event
1267
1268   if (!fBkgrdFileNames) return 1;
1269   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1270   if (nBkgrdFiles == 0) return 1;
1271
1272   // get the number of signal events
1273   if (nEvents <= 0) {
1274     AliRunLoader* runLoader = 
1275         AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1276     if (!runLoader) return 1;
1277     
1278     nEvents = runLoader->GetNumberOfEvents();
1279     delete runLoader;
1280   }
1281
1282   Int_t result = 0;
1283   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1284     // get the number of background events
1285     const char* fileName = ((TObjString*)
1286                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1287     AliRunLoader* runLoader =
1288       AliRunLoader::Open(fileName, "BKGRD");
1289     if (!runLoader) continue;
1290     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1291     delete runLoader;
1292   
1293     // get or calculate the number of signal per background events
1294     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1295     if (nSignalPerBkgrd <= 0) {
1296       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1297     } else if (result && (result != nSignalPerBkgrd)) {
1298       AliInfo(Form("the number of signal events per background event "
1299                    "will be changed from %d to %d for stream %d", 
1300                    nSignalPerBkgrd, result, iBkgrdFile+1));
1301       nSignalPerBkgrd = result;
1302     }
1303
1304     if (!result) result = nSignalPerBkgrd;
1305     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1306       AliWarning(Form("not enough background events (%d) for %d signal events "
1307                       "using %d signal per background events for stream %d",
1308                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1309     }
1310   }
1311
1312   return result;
1313 }
1314
1315 //_____________________________________________________________________________
1316 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1317 {
1318 // check whether detName is contained in detectors
1319 // if yes, it is removed from detectors
1320
1321   // check if all detectors are selected
1322   if ((detectors.CompareTo("ALL") == 0) ||
1323       detectors.BeginsWith("ALL ") ||
1324       detectors.EndsWith(" ALL") ||
1325       detectors.Contains(" ALL ")) {
1326     detectors = "ALL";
1327     return kTRUE;
1328   }
1329
1330   // search for the given detector
1331   Bool_t result = kFALSE;
1332   if ((detectors.CompareTo(detName) == 0) ||
1333       detectors.BeginsWith(detName+" ") ||
1334       detectors.EndsWith(" "+detName) ||
1335       detectors.Contains(" "+detName+" ")) {
1336     detectors.ReplaceAll(detName, "");
1337     result = kTRUE;
1338   }
1339
1340   // clean up the detectors string
1341   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1342   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1343   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1344
1345   return result;
1346 }
1347
1348 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName) 
1349 {
1350 //
1351 // Steering routine  to convert raw data in directory rawDirectory/ to fake SDigits. 
1352 // These can be used for embedding of MC tracks into RAW data using the standard 
1353 // merging procedure.
1354 //
1355 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1356 //
1357     if (!gAlice) {
1358         AliError("no gAlice object. Restart aliroot and try again.");
1359         return kFALSE;
1360     }
1361     if (gAlice->Modules()->GetEntries() > 0) {
1362         AliError("gAlice was already run. Restart aliroot and try again.");
1363         return kFALSE;
1364     }
1365     
1366     AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1367     StdoutToAliInfo(StderrToAliError(gAlice->Init(fConfigFileName.Data());););
1368 //
1369 //  Initialize CDB     
1370     InitCDBStorage();
1371     AliCDBManager* man = AliCDBManager::Instance();
1372     man->SetRun(0); // Should this come from rawdata header ?
1373     
1374     Int_t iDet;
1375     //
1376     // Get the runloader
1377     AliRunLoader* runLoader = gAlice->GetRunLoader();
1378     //
1379     // Open esd file if available
1380     TFile* esdFile = TFile::Open(esdFileName);
1381     Bool_t esdOK = (esdFile != 0);
1382     AliESD* esd = new AliESD;
1383     TTree* treeESD = 0;
1384     if (esdOK) {
1385         treeESD = (TTree*) esdFile->Get("esdTree");
1386         if (!treeESD) {
1387             AliWarning("No ESD tree found");
1388             esdOK = kFALSE;
1389         } else {
1390             treeESD->SetBranchAddress("ESD", &esd);
1391         }
1392     }
1393     //
1394     // Create the RawReader
1395     AliRawReaderFile* rawReader = new AliRawReaderFile(rawDirectory);
1396     //
1397     // Get list of detectors
1398     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1399     //
1400     // Get Header
1401     AliHeader* header = runLoader->GetHeader();
1402     //
1403     // Event loop
1404     Int_t nev = 0;
1405     while(kTRUE) {
1406         if (!(rawReader->NextEvent())) break;
1407         //
1408         // Detector loop
1409         for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1410             AliModule* det = (AliModule*) detArray->At(iDet);
1411             AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1412             det->Raw2SDigits(rawReader);
1413             rawReader->Reset();
1414         } // detectors
1415
1416         //
1417         //  If ESD information available obtain reconstructed vertex and store in header.
1418         if (esdOK) {
1419             treeESD->GetEvent(nev);
1420             const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1421             Double_t position[3];
1422             esdVertex->GetXYZ(position);
1423             AliGenEventHeader* mcHeader = new  AliGenEventHeader("ESD");
1424             TArrayF mcV;
1425             mcV.Set(3);
1426             for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1427             mcHeader->SetPrimaryVertex(mcV);
1428             header->Reset(0,nev);
1429             header->SetGenEventHeader(mcHeader);
1430             printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1431         }
1432         nev++;
1433 //
1434 //      Finish the event
1435         runLoader->TreeE()->Fill();
1436         runLoader->SetNextEvent();
1437     } // events
1438  
1439     delete rawReader;
1440 //
1441 //  Finish the run 
1442     runLoader->CdGAFile();
1443     runLoader->WriteHeader("OVERWRITE");
1444     runLoader->WriteRunLoader();
1445
1446     return kTRUE;
1447 }