]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
New feature in AliSimulation which allows to write a separate raw-data file with...
[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 #include "AliHLTSimulation.h"
141 #include "AliQADataMakerSteer.h"
142
143 ClassImp(AliSimulation)
144
145 AliSimulation *AliSimulation::fgInstance = 0;
146 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
147
148 //_____________________________________________________________________________
149 AliSimulation::AliSimulation(const char* configFileName,
150                              const char* name, const char* title) :
151   TNamed(name, title),
152
153   fRunGeneration(kTRUE),
154   fRunSimulation(kTRUE),
155   fLoadAlignFromCDB(kTRUE),
156   fLoadAlObjsListOfDets("ALL"),
157   fMakeSDigits("ALL"),
158   fMakeDigits("ALL"),
159   fMakeTrigger(""),
160   fMakeDigitsFromHits(""),
161   fWriteRawData(""),
162   fRawDataFileName(""),
163   fDeleteIntermediateFiles(kFALSE),
164   fWriteSelRawData(kFALSE),
165   fStopOnError(kFALSE),
166
167   fNEvents(1),
168   fConfigFileName(configFileName),
169   fGAliceFileName("galice.root"),
170   fEventsPerFile(),
171   fBkgrdFileNames(NULL),
172   fAlignObjArray(NULL),
173   fUseBkgrdVertex(kTRUE),
174   fRegionOfInterest(kFALSE),
175   fCDBUri(""),
176   fSpecCDBUri(),
177   fRun(-1),
178   fSeed(0),
179   fInitCDBCalled(kFALSE),
180   fInitRunNumberCalled(kFALSE),
181   fSetRunNumberFromDataCalled(kFALSE),
182   fEmbeddingFlag(kFALSE),
183   fRunQA(kTRUE), 
184   fRunHLT("default")
185 {
186 // create simulation object with default parameters
187   fgInstance = this;
188   SetGAliceFile("galice.root");
189   
190 // for QA
191    for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) 
192         fQACycles[iDet] = 999999;
193 }
194
195 //_____________________________________________________________________________
196 AliSimulation::AliSimulation(const AliSimulation& sim) :
197   TNamed(sim),
198
199   fRunGeneration(sim.fRunGeneration),
200   fRunSimulation(sim.fRunSimulation),
201   fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
202   fLoadAlObjsListOfDets(sim.fLoadAlObjsListOfDets),
203   fMakeSDigits(sim.fMakeSDigits),
204   fMakeDigits(sim.fMakeDigits),
205   fMakeTrigger(sim.fMakeTrigger),
206   fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
207   fWriteRawData(sim.fWriteRawData),
208   fRawDataFileName(""),
209   fDeleteIntermediateFiles(kFALSE),
210   fWriteSelRawData(kFALSE),
211   fStopOnError(sim.fStopOnError),
212
213   fNEvents(sim.fNEvents),
214   fConfigFileName(sim.fConfigFileName),
215   fGAliceFileName(sim.fGAliceFileName),
216   fEventsPerFile(),
217   fBkgrdFileNames(NULL),
218   fAlignObjArray(NULL),
219   fUseBkgrdVertex(sim.fUseBkgrdVertex),
220   fRegionOfInterest(sim.fRegionOfInterest),
221   fCDBUri(sim.fCDBUri),
222   fSpecCDBUri(),
223   fRun(-1),
224   fSeed(0),
225   fInitCDBCalled(sim.fInitCDBCalled),
226   fInitRunNumberCalled(sim.fInitRunNumberCalled),
227   fSetRunNumberFromDataCalled(sim.fSetRunNumberFromDataCalled),
228   fEmbeddingFlag(sim.fEmbeddingFlag),
229   fRunQA(kTRUE), 
230   fRunHLT(sim.fRunHLT)
231 {
232 // copy constructor
233
234   for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
235     if (!sim.fEventsPerFile[i]) continue;
236     fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
237   }
238
239   fBkgrdFileNames = new TObjArray;
240   for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
241     if (!sim.fBkgrdFileNames->At(i)) continue;
242     fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
243   }
244
245   for (Int_t i = 0; i < sim.fSpecCDBUri.GetEntriesFast(); i++) {
246     if (sim.fSpecCDBUri[i]) fSpecCDBUri.Add(sim.fSpecCDBUri[i]->Clone());
247   }
248   fgInstance = this;
249
250 // for QA
251    for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) 
252         fQACycles[iDet] = sim.fQACycles[iDet];
253 }
254
255 //_____________________________________________________________________________
256 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
257 {
258 // assignment operator
259
260   this->~AliSimulation();
261   new(this) AliSimulation(sim);
262   return *this;
263 }
264
265 //_____________________________________________________________________________
266 AliSimulation::~AliSimulation()
267 {
268 // clean up
269
270   fEventsPerFile.Delete();
271 //  if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
272 //  delete fAlignObjArray; fAlignObjArray=0;
273
274   if (fBkgrdFileNames) {
275     fBkgrdFileNames->Delete();
276     delete fBkgrdFileNames;
277   }
278
279   fSpecCDBUri.Delete();
280   if (fgInstance==this) fgInstance = 0;
281
282   AliCodeTimer::Instance()->Print();
283 }
284
285
286 //_____________________________________________________________________________
287 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
288 {
289 // set the number of events for one run
290
291   fNEvents = nEvents;
292 }
293
294 //_____________________________________________________________________________
295 void AliSimulation::InitCDB()
296 {
297 // activate a default CDB storage
298 // First check if we have any CDB storage set, because it is used 
299 // to retrieve the calibration and alignment constants
300
301   if (fInitCDBCalled) return;
302   fInitCDBCalled = kTRUE;
303
304   AliCDBManager* man = AliCDBManager::Instance();
305   if (man->IsDefaultStorageSet())
306   {
307     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
308     AliWarning("Default CDB storage has been already set !");
309     AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
310     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
311     fCDBUri = man->GetDefaultStorage()->GetURI();
312   }
313   else {
314     if (fCDBUri.Length() > 0) 
315     {
316         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
317         AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
318         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
319     } else {
320         fCDBUri="local://$ALICE_ROOT";
321         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
322         AliWarning("Default CDB storage not yet set !!!!");
323         AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
324         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
325                 
326     }
327     man->SetDefaultStorage(fCDBUri);
328   }
329
330   // Now activate the detector specific CDB storage locations
331   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
332     TObject* obj = fSpecCDBUri[i];
333     if (!obj) continue;
334     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
335     AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
336     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
337     man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
338   }
339       
340 }
341
342 //_____________________________________________________________________________
343 void AliSimulation::InitRunNumber(){
344 // check run number. If not set, set it to 0 !!!!
345   
346   if (fInitRunNumberCalled) return;
347   fInitRunNumberCalled = kTRUE;
348   
349   AliCDBManager* man = AliCDBManager::Instance();
350   if (man->GetRun() >= 0)
351   {
352         AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
353                         "Use external variable DC_RUN or AliSimulation::SetRun()!"));
354   }
355     
356   if(fRun >= 0) {
357         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
358         AliDebug(2, Form("Setting CDB run number to: %d",fRun));
359         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
360   } else {
361         fRun=0;
362         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
363         AliWarning("Run number not yet set !!!!");
364         AliWarning(Form("Setting it now to: %d", fRun));
365         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
366         
367   }
368   man->SetRun(fRun);
369
370   man->Print();
371
372 }
373
374 //_____________________________________________________________________________
375 void AliSimulation::SetCDBLock() {
376   // Set CDB lock: from now on it is forbidden to reset the run number
377   // or the default storage or to activate any further storage!
378   
379   AliCDBManager::Instance()->SetLock(1);
380 }
381
382 //_____________________________________________________________________________
383 void AliSimulation::SetDefaultStorage(const char* uri) {
384 // Store the desired default CDB storage location
385 // Activate it later within the Run() method
386
387   fCDBUri = uri;
388
389 }
390
391 //_____________________________________________________________________________
392 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
393 // Store a detector-specific CDB storage location
394 // Activate it later within the Run() method
395
396   AliCDBPath aPath(calibType);
397   if(!aPath.IsValid()){
398         AliError(Form("Not a valid path: %s", calibType));
399         return;
400   }
401
402   TObject* obj = fSpecCDBUri.FindObject(calibType);
403   if (obj) fSpecCDBUri.Remove(obj);
404   fSpecCDBUri.Add(new TNamed(calibType, uri));
405
406 }
407
408 //_____________________________________________________________________________
409 void AliSimulation::SetRunNumber(Int_t run)
410 {
411 // sets run number
412 // Activate it later within the Run() method
413
414         fRun = run;
415 }
416
417 //_____________________________________________________________________________
418 void AliSimulation::SetSeed(Int_t seed)
419 {
420 // sets seed number
421 // Activate it later within the Run() method
422
423         fSeed = seed;
424 }
425
426 //_____________________________________________________________________________
427 Bool_t AliSimulation::SetRunNumberFromData()
428 {
429   // Set the CDB manager run number
430   // The run number is retrieved from gAlice
431
432     if (fSetRunNumberFromDataCalled) return kTRUE;
433     fSetRunNumberFromDataCalled = kTRUE;    
434   
435     AliCDBManager* man = AliCDBManager::Instance();
436     Int_t runData = -1, runCDB = -1;
437   
438     AliRunLoader* runLoader = LoadRun("READ");
439     if (!runLoader) return kFALSE;
440     else {
441         runData = runLoader->GetAliRun()->GetHeader()->GetRun();
442         delete runLoader;
443     }
444   
445     runCDB = man->GetRun();
446     if(runCDB >= 0) {
447         if (runCDB != runData) {
448                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
449                 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
450                 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
451                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");    
452         }
453         
454     }
455       
456     man->SetRun(runData);
457     fRun = runData;
458     
459     if(man->GetRun() < 0) {
460         AliError("Run number not properly initalized!");
461         return kFALSE;
462     }
463   
464     man->Print();
465     
466     return kTRUE;
467 }
468
469 //_____________________________________________________________________________
470 void AliSimulation::SetConfigFile(const char* fileName)
471 {
472 // set the name of the config file
473
474   fConfigFileName = fileName;
475 }
476
477 //_____________________________________________________________________________
478 void AliSimulation::SetGAliceFile(const char* fileName)
479 {
480 // set the name of the galice file
481 // the path is converted to an absolute one if it is relative
482
483   fGAliceFileName = fileName;
484   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
485     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
486                                                 fGAliceFileName);
487     fGAliceFileName = absFileName;
488     delete[] absFileName;
489   }
490
491   AliDebug(2, Form("galice file name set to %s", fileName));
492 }
493
494 //_____________________________________________________________________________
495 void AliSimulation::SetEventsPerFile(const char* detector, const char* type, 
496                                      Int_t nEvents)
497 {
498 // set the number of events per file for the given detector and data type
499 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
500
501   TNamed* obj = new TNamed(detector, type);
502   obj->SetUniqueID(nEvents);
503   fEventsPerFile.Add(obj);
504 }
505
506 //_____________________________________________________________________________
507 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
508 {
509   // Read the alignment objects from CDB.
510   // Each detector is supposed to have the
511   // alignment objects in DET/Align/Data CDB path.
512   // All the detector objects are then collected,
513   // sorted by geometry level (starting from ALIC) and
514   // then applied to the TGeo geometry.
515   // Finally an overlaps check is performed.
516
517   if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
518     AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
519     return kFALSE;
520   }  
521   
522   // initialize CDB storage, run number, set CDB lock
523   InitCDB();
524 //  if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
525   SetCDBLock();
526     
527   Bool_t delRunLoader = kFALSE;
528   if (!runLoader) {
529     runLoader = LoadRun("READ");
530     if (!runLoader) return kFALSE;
531     delRunLoader = kTRUE;
532   }
533   
534   // Export ideal geometry 
535   if(!gAlice->IsRootGeometry()) AliGeomManager::GetGeometry()->Export("geometry.root");
536
537   // Load alignment data from CDB and apply to geometry through AliGeomManager
538   if(fLoadAlignFromCDB){
539     
540     TString detStr = fLoadAlObjsListOfDets;
541     TString loadAlObjsListOfDets = "";
542     
543     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
544     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
545       AliModule* det = (AliModule*) detArray->At(iDet);
546       if (!det || !det->IsActive()) continue;
547       if (IsSelected(det->GetName(), detStr)) {
548         //add det to list of dets to be aligned from CDB
549         loadAlObjsListOfDets += det->GetName();
550         loadAlObjsListOfDets += " ";
551       }
552     } // end loop over detectors
553     loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
554     AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
555   }else{
556     // Check if the array with alignment objects was
557     // provided by the user. If yes, apply the objects
558     // to the present TGeo geometry
559     if (fAlignObjArray) {
560       if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
561         AliError("The misalignment of one or more volumes failed!"
562                  "Compare the list of simulated detectors and the list of detector alignment data!");
563         if (delRunLoader) delete runLoader;
564         return kFALSE;
565       }
566     }
567   }
568
569   // Update the internal geometry of modules (ITS needs it)
570   TString detStr = fLoadAlObjsListOfDets;
571   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
572   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
573
574     AliModule* det = (AliModule*) detArray->At(iDet);
575     if (!det || !det->IsActive()) continue;
576     if (IsSelected(det->GetName(), detStr)) {
577       det->UpdateInternalGeometry();
578     }
579   } // end loop over detectors
580
581
582   if (delRunLoader) delete runLoader;
583
584   return kTRUE;
585 }
586
587 //_____________________________________________________________________________
588 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
589 {
590 // add a file with background events for merging
591
592   TObjString* fileNameStr = new TObjString(fileName);
593   fileNameStr->SetUniqueID(nSignalPerBkgrd);
594   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
595   fBkgrdFileNames->Add(fileNameStr);
596 }
597
598 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
599 {
600 // add a file with background events for embeddin
601   MergeWith(fileName, nSignalPerBkgrd);
602   fEmbeddingFlag = kTRUE;
603 }
604
605 //_____________________________________________________________________________
606 Bool_t AliSimulation::Run(Int_t nEvents)
607 {
608 // run the generation, simulation and digitization
609
610  
611   AliCodeTimerAuto("")
612   
613   // Load run number and seed from environmental vars
614   ProcessEnvironmentVars();
615
616   gRandom->SetSeed(fSeed);
617    
618   if (nEvents > 0) fNEvents = nEvents;
619
620   // generation and simulation -> hits
621   if (fRunGeneration) {
622     if (!RunSimulation()) if (fStopOnError) return kFALSE;
623   }
624            
625   // initialize CDB storage from external environment
626   // (either CDB manager or AliSimulation setters),
627   // if not already done in RunSimulation()
628   InitCDB();
629   
630   // Set run number in CDBManager from data 
631   // From this point on the run number must be always loaded from data!
632   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
633   
634   // Set CDB lock: from now on it is forbidden to reset the run number
635   // or the default storage or to activate any further storage!
636   SetCDBLock();
637
638   // If RunSimulation was not called, load the geometry and misalign it
639   if (!AliGeomManager::GetGeometry()) {
640     // Initialize the geometry manager
641     AliGeomManager::LoadGeometry("geometry.root");
642     if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
643     // Misalign geometry
644     if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
645   }
646
647   // hits -> summable digits
648   if (!fMakeSDigits.IsNull()) {
649     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
650  
651   }
652   
653
654   
655   // summable digits -> digits  
656   if (!fMakeDigits.IsNull()) {
657     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
658       if (fStopOnError) return kFALSE;
659     }
660    }
661
662   
663   
664   // hits -> digits
665   if (!fMakeDigitsFromHits.IsNull()) {
666     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
667       AliWarning(Form("Merging and direct creation of digits from hits " 
668                  "was selected for some detectors. "
669                  "No merging will be done for the following detectors: %s",
670                  fMakeDigitsFromHits.Data()));
671     }
672     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
673       if (fStopOnError) return kFALSE;
674     }
675   }
676
677   
678   
679   // digits -> trigger
680   if (!RunTrigger(fMakeTrigger,fMakeDigits)) {
681     if (fStopOnError) return kFALSE;
682   }
683
684   
685   
686   // digits -> raw data
687   if (!fWriteRawData.IsNull()) {
688     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
689                       fDeleteIntermediateFiles,fWriteSelRawData)) {
690       if (fStopOnError) return kFALSE;
691     }
692   }
693
694   
695   
696   // run HLT simulation
697   if (!fRunHLT.IsNull()) {
698     if (!RunHLT()) {
699       if (fStopOnError) return kFALSE;
700     }
701   }
702   
703   //QA
704         if (fRunQA) {
705                 Bool_t rv = RunQA() ; 
706                 if (!rv)
707                         if (fStopOnError) 
708                                 return kFALSE ;         
709         }
710
711   // Cleanup of CDB manager: cache and active storages!
712   AliCDBManager::Instance()->ClearCache();
713
714   return kTRUE;
715 }
716
717 //_____________________________________________________________________________
718 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
719 {
720   // run the trigger
721
722   AliCodeTimerAuto("")
723
724   // initialize CDB storage from external environment
725   // (either CDB manager or AliSimulation setters),
726   // if not already done in RunSimulation()
727   InitCDB();
728   
729   // Set run number in CDBManager from data 
730   // From this point on the run number must be always loaded from data!
731   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
732   
733   // Set CDB lock: from now on it is forbidden to reset the run number
734   // or the default storage or to activate any further storage!
735   SetCDBLock();
736    
737    AliRunLoader* runLoader = LoadRun("READ");
738    if (!runLoader) return kFALSE;
739    TString trconfiguration = config;
740
741    if (trconfiguration.IsNull()) {
742      if (gAlice->GetTriggerDescriptor() != "") {
743        trconfiguration = gAlice->GetTriggerDescriptor();
744      }
745      else
746        AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
747    }
748
749    runLoader->MakeTree( "GG" );
750    AliCentralTrigger* aCTP = runLoader->GetTrigger();
751    // Load Configuration
752    if (!aCTP->LoadConfiguration( trconfiguration ))
753      return kFALSE;
754
755    // digits -> trigger
756    if( !aCTP->RunTrigger( runLoader , detectors ) ) {
757       if (fStopOnError) {
758         //  delete aCTP;
759         return kFALSE;
760       }
761    }
762
763    delete runLoader;
764
765    return kTRUE;
766 }
767
768 //_____________________________________________________________________________
769 Bool_t AliSimulation::WriteTriggerRawData()
770 {
771   // Writes the CTP (trigger) DDL raw data
772   // Details of the format are given in the
773   // trigger TDR - pages 134 and 135.
774   AliCTPRawData writer;
775   writer.RawData();
776
777   return kTRUE;
778 }
779
780 //_____________________________________________________________________________
781 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
782 {
783 // run the generation and simulation
784
785   AliCodeTimerAuto("")
786
787   // initialize CDB storage and run number from external environment
788   // (either CDB manager or AliSimulation setters)
789   InitCDB();
790   InitRunNumber();
791   SetCDBLock();
792   
793   if (!gAlice) {
794     AliError("no gAlice object. Restart aliroot and try again.");
795     return kFALSE;
796   }
797   if (gAlice->Modules()->GetEntries() > 0) {
798     AliError("gAlice was already run. Restart aliroot and try again.");
799     return kFALSE;
800   }
801
802   AliInfo(Form("initializing gAlice with config file %s",
803           fConfigFileName.Data()));
804   StdoutToAliInfo(StderrToAliError(
805     gAlice->Init(fConfigFileName.Data());
806   ););
807   
808   // Get the trigger descriptor string
809   // Either from AliSimulation or from
810   // gAlice
811   if (fMakeTrigger.IsNull()) {
812     if (gAlice->GetTriggerDescriptor() != "")
813       fMakeTrigger = gAlice->GetTriggerDescriptor();
814   }
815   else
816     gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
817
818   // Set run number in CDBManager
819   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
820
821   AliRunLoader* runLoader = gAlice->GetRunLoader();
822   if (!runLoader) {
823              AliError(Form("gAlice has no run loader object. "
824                              "Check your config file: %s", fConfigFileName.Data()));
825              return kFALSE;
826   }
827   SetGAliceFile(runLoader->GetFileName());
828       
829   // Misalign geometry
830 #if ROOT_VERSION_CODE < 331527
831   AliGeomManager::SetGeometry(gGeoManager);
832   MisalignGeometry(runLoader);
833 #endif
834
835 //   AliRunLoader* runLoader = gAlice->GetRunLoader();
836 //   if (!runLoader) {
837 //     AliError(Form("gAlice has no run loader object. "
838 //                   "Check your config file: %s", fConfigFileName.Data()));
839 //     return kFALSE;
840 //   }
841 //   SetGAliceFile(runLoader->GetFileName());
842
843   if (!gAlice->Generator()) {
844     AliError(Form("gAlice has no generator object. "
845                   "Check your config file: %s", fConfigFileName.Data()));
846     return kFALSE;
847   }
848   if (nEvents <= 0) nEvents = fNEvents;
849
850   // get vertex from background file in case of merging
851   if (fUseBkgrdVertex &&
852       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
853     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
854     const char* fileName = ((TObjString*)
855                             (fBkgrdFileNames->At(0)))->GetName();
856     AliInfo(Form("The vertex will be taken from the background "
857                  "file %s with nSignalPerBackground = %d", 
858                  fileName, signalPerBkgrd));
859     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
860     gAlice->Generator()->SetVertexGenerator(vtxGen);
861   }
862
863   if (!fRunSimulation) {
864     gAlice->Generator()->SetTrackingFlag(0);
865   }
866
867   // set the number of events per file for given detectors and data types
868   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
869     if (!fEventsPerFile[i]) continue;
870     const char* detName = fEventsPerFile[i]->GetName();
871     const char* typeName = fEventsPerFile[i]->GetTitle();
872     TString loaderName(detName);
873     loaderName += "Loader";
874     AliLoader* loader = runLoader->GetLoader(loaderName);
875     if (!loader) {
876       AliError(Form("RunSimulation", "no loader for %s found\n"
877                     "Number of events per file not set for %s %s", 
878                     detName, typeName, detName));
879       continue;
880     }
881     AliDataLoader* dataLoader = 
882       loader->GetDataLoader(typeName);
883     if (!dataLoader) {
884       AliError(Form("no data loader for %s found\n"
885                     "Number of events per file not set for %s %s", 
886                     typeName, detName, typeName));
887       continue;
888     }
889     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
890     AliDebug(1, Form("number of events per file set to %d for %s %s",
891                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
892   }
893
894   AliInfo("running gAlice");
895   StdoutToAliInfo(StderrToAliError(
896     gAlice->Run(nEvents);
897   ););
898
899   delete runLoader;
900
901   return kTRUE;
902 }
903
904 //_____________________________________________________________________________
905 Bool_t AliSimulation::RunSDigitization(const char* detectors)
906 {
907 // run the digitization and produce summable digits
908
909   AliCodeTimerAuto("")
910
911   // initialize CDB storage, run number, set CDB lock
912   InitCDB();
913   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
914   SetCDBLock();
915   
916   AliRunLoader* runLoader = LoadRun();
917   if (!runLoader) return kFALSE;
918
919   TString detStr = detectors;
920   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
921   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
922     AliModule* det = (AliModule*) detArray->At(iDet);
923     if (!det || !det->IsActive()) continue;
924     if (IsSelected(det->GetName(), detStr)) {
925       AliInfo(Form("creating summable digits for %s", det->GetName()));
926       AliCodeTimerAuto(Form("creating summable digits for %s", det->GetName()));
927           
928       det->Hits2SDigits();
929     }
930   }
931
932   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
933     AliError(Form("the following detectors were not found: %s",
934                   detStr.Data()));
935     if (fStopOnError) return kFALSE;
936   }
937
938   delete runLoader;
939
940   return kTRUE;
941 }
942
943
944 //_____________________________________________________________________________
945 Bool_t AliSimulation::RunDigitization(const char* detectors, 
946                                       const char* excludeDetectors)
947 {
948 // run the digitization and produce digits from sdigits
949
950   AliCodeTimerAuto("")
951
952   // initialize CDB storage, run number, set CDB lock
953   InitCDB();
954   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
955   SetCDBLock();
956   
957   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
958   if (gAlice) delete gAlice;
959   gAlice = NULL;
960
961   Int_t nStreams = 1;
962   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
963   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
964   AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
965   // manager->SetEmbeddingFlag(fEmbeddingFlag);
966   manager->SetInputStream(0, fGAliceFileName.Data());
967   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
968     const char* fileName = ((TObjString*)
969                             (fBkgrdFileNames->At(iStream-1)))->GetName();
970     manager->SetInputStream(iStream, fileName);
971   }
972
973   TString detStr = detectors;
974   TString detExcl = excludeDetectors;
975   manager->GetInputStream(0)->ImportgAlice();
976   AliRunLoader* runLoader = 
977     AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
978   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
979   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
980     AliModule* det = (AliModule*) detArray->At(iDet);
981     if (!det || !det->IsActive()) continue;
982     if (IsSelected(det->GetName(), detStr) && 
983         !IsSelected(det->GetName(), detExcl)) {
984       AliDigitizer* digitizer = det->CreateDigitizer(manager);
985       
986       if (!digitizer) {
987         AliError(Form("no digitizer for %s", det->GetName()));
988         if (fStopOnError) return kFALSE;
989       } else {
990         digitizer->SetRegionOfInterest(fRegionOfInterest);
991       }
992     }
993   }
994
995   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
996     AliError(Form("the following detectors were not found: %s", 
997                   detStr.Data()));
998     if (fStopOnError) return kFALSE;
999   }
1000
1001   if (!manager->GetListOfTasks()->IsEmpty()) {
1002     AliInfo("executing digitization");
1003     manager->Exec("");
1004   }
1005
1006   delete manager;
1007
1008   return kTRUE;
1009 }
1010
1011 //_____________________________________________________________________________
1012 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1013 {
1014 // run the digitization and produce digits from hits
1015
1016   AliCodeTimerAuto("")
1017
1018   // initialize CDB storage, run number, set CDB lock
1019   InitCDB();
1020   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1021   SetCDBLock();
1022   
1023   AliRunLoader* runLoader = LoadRun("READ");
1024   if (!runLoader) return kFALSE;
1025
1026   TString detStr = detectors;
1027   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1028   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1029     AliModule* det = (AliModule*) detArray->At(iDet);
1030     if (!det || !det->IsActive()) continue;
1031     if (IsSelected(det->GetName(), detStr)) {
1032       AliInfo(Form("creating digits from hits for %s", det->GetName()));
1033       det->Hits2Digits();
1034     }
1035   }
1036
1037   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1038     AliError(Form("the following detectors were not found: %s", 
1039                   detStr.Data()));
1040     if (fStopOnError) return kFALSE;
1041   }
1042
1043   delete runLoader;
1044   //PH Temporary fix to avoid interference with the PHOS loder/getter
1045   //PH The problem has to be solved in more general way 09/06/05
1046
1047   return kTRUE;
1048 }
1049
1050 //_____________________________________________________________________________
1051 Bool_t AliSimulation::WriteRawData(const char* detectors, 
1052                                    const char* fileName,
1053                                    Bool_t deleteIntermediateFiles,
1054                                    Bool_t selrawdata)
1055 {
1056 // convert the digits to raw data
1057 // First DDL raw data files for the given detectors are created.
1058 // If a file name is given, the DDL files are then converted to a DATE file.
1059 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
1060 // afterwards.
1061 // If the file name has the extension ".root", the DATE file is converted
1062 // to a root file.
1063 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1064 // 'selrawdata' flag can be used to enable writing of detectors raw data
1065 // accoring to the trigger cluster.
1066
1067   AliCodeTimerAuto("")
1068
1069   if (!WriteRawFiles(detectors)) {
1070     if (fStopOnError) return kFALSE;
1071   }
1072
1073   TString dateFileName(fileName);
1074   if (!dateFileName.IsNull()) {
1075     Bool_t rootOutput = dateFileName.EndsWith(".root");
1076     if (rootOutput) dateFileName += ".date";
1077     TString selDateFileName;
1078     if (selrawdata) {
1079       selDateFileName = "selected.";
1080       selDateFileName+= dateFileName;
1081     }
1082     if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1083       if (fStopOnError) return kFALSE;
1084     }
1085     if (deleteIntermediateFiles) {
1086       AliRunLoader* runLoader = LoadRun("READ");
1087       if (runLoader) for (Int_t iEvent = 0; 
1088                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1089         char command[256];
1090         sprintf(command, "rm -r raw%d", iEvent);
1091         gSystem->Exec(command);
1092       }
1093     }
1094
1095     if (rootOutput) {
1096       if (!ConvertDateToRoot(dateFileName, fileName)) {
1097         if (fStopOnError) return kFALSE;
1098       }
1099       if (deleteIntermediateFiles) {
1100         gSystem->Unlink(dateFileName);
1101       }
1102       if (selrawdata) {
1103         TString selFileName = "selected.";
1104         selFileName        += fileName;
1105         if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1106           if (fStopOnError) return kFALSE;
1107         }
1108         if (deleteIntermediateFiles) {
1109           gSystem->Unlink(selDateFileName);
1110         }
1111       }
1112     }
1113   }
1114
1115   return kTRUE;
1116 }
1117
1118 //_____________________________________________________________________________
1119 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1120 {
1121 // convert the digits to raw data DDL files
1122
1123   AliCodeTimerAuto("")
1124   
1125   AliRunLoader* runLoader = LoadRun("READ");
1126   if (!runLoader) return kFALSE;
1127
1128   // write raw data to DDL files
1129   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1130     AliInfo(Form("processing event %d", iEvent));
1131     runLoader->GetEvent(iEvent);
1132     TString baseDir = gSystem->WorkingDirectory();
1133     char dirName[256];
1134     sprintf(dirName, "raw%d", iEvent);
1135     gSystem->MakeDirectory(dirName);
1136     if (!gSystem->ChangeDirectory(dirName)) {
1137       AliError(Form("couldn't change to directory %s", dirName));
1138       if (fStopOnError) return kFALSE; else continue;
1139     }
1140
1141     TString detStr = detectors;
1142     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1143     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1144       AliModule* det = (AliModule*) detArray->At(iDet);
1145       if (!det || !det->IsActive()) continue;
1146       if (IsSelected(det->GetName(), detStr)) {
1147         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1148         det->Digits2Raw();
1149       }
1150     }
1151
1152     if (!WriteTriggerRawData())
1153       if (fStopOnError) return kFALSE;
1154
1155     gSystem->ChangeDirectory(baseDir);
1156     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1157       AliError(Form("the following detectors were not found: %s", 
1158                     detStr.Data()));
1159       if (fStopOnError) return kFALSE;
1160     }
1161   }
1162
1163   delete runLoader;
1164   
1165   return kTRUE;
1166 }
1167
1168 //_____________________________________________________________________________
1169 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1170                                             const char* selDateFileName)
1171 {
1172 // convert raw data DDL files to a DATE file with the program "dateStream"
1173 // The second argument is not empty when the user decides to write
1174 // the detectors raw data according to the trigger cluster.
1175
1176   AliCodeTimerAuto("")
1177   
1178   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1179   if (!path) {
1180     AliError("the program dateStream was not found");
1181     if (fStopOnError) return kFALSE;
1182   } else {
1183     delete[] path;
1184   }
1185
1186   AliRunLoader* runLoader = LoadRun("READ");
1187   if (!runLoader) return kFALSE;
1188
1189   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1190   Bool_t selrawdata = kFALSE;
1191   if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1192
1193   char command[256];
1194   // Note the option -s. It is used in order to avoid
1195   // the generation of SOR/EOR events.
1196   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1197           dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1198   FILE* pipe = gSystem->OpenPipe(command, "w");
1199
1200   Int_t selEvents = 0;
1201   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1202     fprintf(pipe, "GDC\n");
1203     Float_t ldc = 0;
1204     Int_t prevLDC = -1;
1205
1206     if (selrawdata) {
1207       // Check if the event was triggered by CTP
1208       runLoader->GetEvent(iEvent);
1209       if (!runLoader->LoadTrigger()) {
1210         AliCentralTrigger *aCTP = runLoader->GetTrigger();
1211         if (aCTP->GetClassMask()) selEvents++;
1212       }
1213       else {
1214         AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1215         selrawdata = kFALSE;
1216       }
1217     }
1218
1219     // loop over detectors and DDLs
1220     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1221       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1222
1223         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1224         Int_t ldcID = Int_t(ldc + 0.0001);
1225         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1226
1227         char rawFileName[256];
1228         sprintf(rawFileName, "raw%d/%s", 
1229                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1230
1231         // check existence and size of raw data file
1232         FILE* file = fopen(rawFileName, "rb");
1233         if (!file) continue;
1234         fseek(file, 0, SEEK_END);
1235         unsigned long size = ftell(file);
1236         fclose(file);
1237         if (!size) continue;
1238
1239         if (ldcID != prevLDC) {
1240           fprintf(pipe, " LDC Id %d\n", ldcID);
1241           prevLDC = ldcID;
1242         }
1243         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1244       }
1245     }
1246   }
1247
1248   Int_t result = gSystem->ClosePipe(pipe);
1249
1250   if (!selrawdata && selEvents > 0) {
1251     delete runLoader;
1252     return (result == 0);
1253   }
1254
1255   AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1256   
1257   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1258           selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1259   FILE* pipe2 = gSystem->OpenPipe(command, "w");
1260
1261   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1262
1263     // Get the trigger decision and cluster
1264     TString detClust;
1265     runLoader->GetEvent(iEvent);
1266     if (!runLoader->LoadTrigger()) {
1267       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1268       if (aCTP->GetClassMask() == 0) continue;
1269       detClust = aCTP->GetTriggeredDetectors();
1270       AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1271     }
1272
1273     fprintf(pipe2, "GDC\n");
1274     Float_t ldc = 0;
1275     Int_t prevLDC = -1;
1276
1277     // loop over detectors and DDLs
1278     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1279       // Write only raw data from detectors that
1280       // are contained in the trigger cluster(s)
1281       if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1282
1283       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1284
1285         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1286         Int_t ldcID = Int_t(ldc + 0.0001);
1287         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1288
1289         char rawFileName[256];
1290         sprintf(rawFileName, "raw%d/%s", 
1291                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1292
1293         // check existence and size of raw data file
1294         FILE* file = fopen(rawFileName, "rb");
1295         if (!file) continue;
1296         fseek(file, 0, SEEK_END);
1297         unsigned long size = ftell(file);
1298         fclose(file);
1299         if (!size) continue;
1300
1301         if (ldcID != prevLDC) {
1302           fprintf(pipe2, " LDC Id %d\n", ldcID);
1303           prevLDC = ldcID;
1304         }
1305         fprintf(pipe2, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1306       }
1307     }
1308   }
1309
1310   Int_t result2 = gSystem->ClosePipe(pipe2);
1311
1312   delete runLoader;
1313   return ((result == 0) && (result2 == 0));
1314 }
1315
1316 //_____________________________________________________________________________
1317 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1318                                         const char* rootFileName)
1319 {
1320 // convert a DATE file to a root file with the program "alimdc"
1321
1322   // ALIMDC setup
1323   const Int_t kDBSize = 2000000000;
1324   const Int_t kTagDBSize = 1000000000;
1325   const Bool_t kFilter = kFALSE;
1326   const Int_t kCompression = 1;
1327
1328   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1329   if (!path) {
1330     AliError("the program alimdc was not found");
1331     if (fStopOnError) return kFALSE;
1332   } else {
1333     delete[] path;
1334   }
1335
1336   AliInfo(Form("converting DATE file %s to root file %s", 
1337                dateFileName, rootFileName));
1338
1339   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1340   const char* tagDBFS    = "/tmp/mdc1/tags";
1341
1342   // User defined file system locations
1343   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1344     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1345   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1346     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1347   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1348     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1349
1350   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1351   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1352   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1353
1354   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1355   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1356   gSystem->Exec(Form("mkdir %s",tagDBFS));
1357
1358   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1359                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1360   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1361
1362   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1363   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1364   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1365
1366   return (result == 0);
1367 }
1368
1369
1370 //_____________________________________________________________________________
1371 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1372 {
1373 // delete existing run loaders, open a new one and load gAlice
1374
1375   while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1376   AliRunLoader* runLoader = 
1377     AliRunLoader::Open(fGAliceFileName.Data(), 
1378                        AliConfig::GetDefaultEventFolderName(), mode);
1379   if (!runLoader) {
1380     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1381     return NULL;
1382   }
1383   runLoader->LoadgAlice();
1384   runLoader->LoadHeader();
1385   gAlice = runLoader->GetAliRun();
1386   if (!gAlice) {
1387     AliError(Form("no gAlice object found in file %s", 
1388                   fGAliceFileName.Data()));
1389     return NULL;
1390   }
1391   return runLoader;
1392 }
1393
1394 //_____________________________________________________________________________
1395 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1396 {
1397 // get or calculate the number of signal events per background event
1398
1399   if (!fBkgrdFileNames) return 1;
1400   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1401   if (nBkgrdFiles == 0) return 1;
1402
1403   // get the number of signal events
1404   if (nEvents <= 0) {
1405     AliRunLoader* runLoader = 
1406         AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1407     if (!runLoader) return 1;
1408     
1409     nEvents = runLoader->GetNumberOfEvents();
1410     delete runLoader;
1411   }
1412
1413   Int_t result = 0;
1414   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1415     // get the number of background events
1416     const char* fileName = ((TObjString*)
1417                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1418     AliRunLoader* runLoader =
1419       AliRunLoader::Open(fileName, "BKGRD");
1420     if (!runLoader) continue;
1421     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1422     delete runLoader;
1423   
1424     // get or calculate the number of signal per background events
1425     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1426     if (nSignalPerBkgrd <= 0) {
1427       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1428     } else if (result && (result != nSignalPerBkgrd)) {
1429       AliInfo(Form("the number of signal events per background event "
1430                    "will be changed from %d to %d for stream %d", 
1431                    nSignalPerBkgrd, result, iBkgrdFile+1));
1432       nSignalPerBkgrd = result;
1433     }
1434
1435     if (!result) result = nSignalPerBkgrd;
1436     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1437       AliWarning(Form("not enough background events (%d) for %d signal events "
1438                       "using %d signal per background events for stream %d",
1439                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1440     }
1441   }
1442
1443   return result;
1444 }
1445
1446 //_____________________________________________________________________________
1447 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1448 {
1449 // check whether detName is contained in detectors
1450 // if yes, it is removed from detectors
1451
1452   // check if all detectors are selected
1453   if ((detectors.CompareTo("ALL") == 0) ||
1454       detectors.BeginsWith("ALL ") ||
1455       detectors.EndsWith(" ALL") ||
1456       detectors.Contains(" ALL ")) {
1457     detectors = "ALL";
1458     return kTRUE;
1459   }
1460
1461   // search for the given detector
1462   Bool_t result = kFALSE;
1463   if ((detectors.CompareTo(detName) == 0) ||
1464       detectors.BeginsWith(detName+" ") ||
1465       detectors.EndsWith(" "+detName) ||
1466       detectors.Contains(" "+detName+" ")) {
1467     detectors.ReplaceAll(detName, "");
1468     result = kTRUE;
1469   }
1470
1471   // clean up the detectors string
1472   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1473   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1474   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1475
1476   return result;
1477 }
1478
1479 //_____________________________________________________________________________
1480 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName) 
1481 {
1482 //
1483 // Steering routine  to convert raw data in directory rawDirectory/ to fake SDigits. 
1484 // These can be used for embedding of MC tracks into RAW data using the standard 
1485 // merging procedure.
1486 //
1487 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1488 //
1489     if (!gAlice) {
1490         AliError("no gAlice object. Restart aliroot and try again.");
1491         return kFALSE;
1492     }
1493     if (gAlice->Modules()->GetEntries() > 0) {
1494         AliError("gAlice was already run. Restart aliroot and try again.");
1495         return kFALSE;
1496     }
1497     
1498     AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1499     StdoutToAliInfo(StderrToAliError(gAlice->Init(fConfigFileName.Data());););
1500 //
1501 //  Initialize CDB     
1502     InitCDB();
1503     //AliCDBManager* man = AliCDBManager::Instance();
1504     //man->SetRun(0); // Should this come from rawdata header ?
1505     
1506     Int_t iDet;
1507     //
1508     // Get the runloader
1509     AliRunLoader* runLoader = gAlice->GetRunLoader();
1510     //
1511     // Open esd file if available
1512     TFile* esdFile = TFile::Open(esdFileName);
1513     Bool_t esdOK = (esdFile != 0);
1514     AliESD* esd = new AliESD;
1515     TTree* treeESD = 0;
1516     if (esdOK) {
1517         treeESD = (TTree*) esdFile->Get("esdTree");
1518         if (!treeESD) {
1519             AliWarning("No ESD tree found");
1520             esdOK = kFALSE;
1521         } else {
1522             treeESD->SetBranchAddress("ESD", &esd);
1523         }
1524     }
1525     //
1526     // Create the RawReader
1527     TString fileName(rawDirectory);
1528     AliRawReader* rawReader = 0x0;
1529     if (fileName.EndsWith("/")) {
1530       rawReader = new AliRawReaderFile(fileName);
1531     } else if (fileName.EndsWith(".root")) {
1532       rawReader = new AliRawReaderRoot(fileName);
1533     } else if (!fileName.IsNull()) {
1534       rawReader = new AliRawReaderDate(fileName);
1535       rawReader->SelectEvents(7);
1536     }
1537 //     if (!fEquipIdMap.IsNull() && fRawReader)
1538 //       fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1539     //
1540     // Get list of detectors
1541     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1542     //
1543     // Get Header
1544     AliHeader* header = runLoader->GetHeader();
1545     //
1546     TString detStr = fMakeSDigits;
1547     // Event loop
1548     Int_t nev = 0;
1549     while(kTRUE) {
1550         if (!(rawReader->NextEvent())) break;
1551         //
1552         // Detector loop
1553         for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1554             AliModule* det = (AliModule*) detArray->At(iDet);
1555             if (!det || !det->IsActive()) continue;
1556             if (IsSelected(det->GetName(), detStr)) {
1557               AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1558               det->Raw2SDigits(rawReader);
1559               rawReader->Reset();
1560             }
1561         } // detectors
1562
1563
1564         //
1565         //  If ESD information available obtain reconstructed vertex and store in header.
1566         if (esdOK) {
1567             treeESD->GetEvent(nev);
1568             const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1569             Double_t position[3];
1570             esdVertex->GetXYZ(position);
1571             AliGenEventHeader* mcHeader = new  AliGenEventHeader("ESD");
1572             TArrayF mcV;
1573             mcV.Set(3);
1574             for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1575             mcHeader->SetPrimaryVertex(mcV);
1576             header->Reset(0,nev);
1577             header->SetGenEventHeader(mcHeader);
1578             printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1579         }
1580         nev++;
1581 //
1582 //      Finish the event
1583         runLoader->TreeE()->Fill();
1584         runLoader->SetNextEvent();
1585     } // events
1586  
1587     delete rawReader;
1588 //
1589 //  Finish the run 
1590     runLoader->CdGAFile();
1591     runLoader->WriteHeader("OVERWRITE");
1592     runLoader->WriteRunLoader();
1593
1594     return kTRUE;
1595 }
1596
1597 //_____________________________________________________________________________
1598 Int_t AliSimulation::GetDetIndex(const char* detector)
1599 {
1600   // return the detector index corresponding to detector
1601   Int_t index = -1 ; 
1602   for (index = 0; index < fgkNDetectors ; index++) {
1603     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1604           break ; 
1605   }     
1606   return index ; 
1607 }
1608
1609 //_____________________________________________________________________________
1610 Bool_t AliSimulation::RunHLT()
1611 {
1612   // Run the HLT simulation
1613   // HLT simulation is implemented in HLT/sim/AliHLTSimulation
1614   // Disabled if fRunHLT is empty, default vaule is "default".
1615   // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
1616   // The default simulation depends on the HLT component libraries and their
1617   // corresponding agents which define components and chains to run. See
1618   // http://web.ift.uib.no/~kjeks/doc/alice-hlt/
1619   // http://web.ift.uib.no/~kjeks/doc/alice-hlt/classAliHLTModuleAgent.html
1620   //
1621   // The libraries to be loaded can be specified as an option.
1622   // <pre>
1623   // AliSimulation sim;
1624   // sim.SetRunHLT("libAliHLTSample.so");
1625   // </pre>
1626   // will only load <tt>libAliHLTSample.so</tt>
1627
1628   // Other available options:
1629   // \li loglevel=<i>level</i> <br>
1630   //     logging level for this processing
1631   // \li alilog=off
1632   //     disable redirection of log messages to AliLog class
1633   // \li config=<i>macro</i>
1634   //     configuration macro
1635   // \li localrec=<i>configuration</i>
1636   //     comma separated list of configurations to be run during simulation
1637
1638   int iResult=0;
1639   AliRunLoader* pRunLoader = LoadRun("READ");
1640   if (!pRunLoader) return kFALSE;
1641
1642   // initialize CDB storage, run number, set CDB lock
1643   InitCDB();
1644   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1645   SetCDBLock();
1646   
1647   // load the library dynamically
1648   gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1649
1650   // check for the library version
1651   AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
1652   if (!fctVersion) {
1653     AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
1654     return kFALSE;
1655   }
1656   if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
1657     AliError(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
1658     return kFALSE;
1659   }
1660
1661   // print compile info
1662   typedef void (*CompileInfo)( char*& date, char*& time);
1663   CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
1664   if (fctInfo) {
1665     char* date="";
1666     char* time="";
1667     (*fctInfo)(date, time);
1668     if (!date) date="unknown";
1669     if (!time) time="unknown";
1670     AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
1671   } else {
1672     AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
1673   }
1674
1675   // create instance of the HLT simulation
1676   AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
1677   AliHLTSimulation* pHLT=NULL;
1678   if (fctCreate==NULL || (pHLT=(fctCreate()))==NULL) {
1679     AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
1680     return kFALSE;    
1681   }
1682
1683   // init the HLT simulation
1684   if (fRunHLT.CompareTo("default")==0) fRunHLT="";
1685   AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
1686   if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, fRunHLT.Data())))<0) {
1687     AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
1688   } else {
1689     // run the HLT simulation
1690     AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
1691     if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
1692       AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
1693     }
1694   }
1695
1696   // delete the instance
1697   AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
1698   if (fctDelete==NULL || fctDelete(pHLT)<0) {
1699     AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
1700   }
1701   pHLT=NULL;
1702
1703   return iResult>=0?kTRUE:kFALSE;
1704 }
1705
1706 //_____________________________________________________________________________
1707 Bool_t AliSimulation::RunQA()
1708 {
1709         // run the QA on summable hits, digits or digits
1710
1711         AliQADataMakerSteer qas ; 
1712     qas.SetRunLoader(gAlice->GetRunLoader()) ;
1713
1714         Bool_t rv =  qas.Run("ALL", AliQA::kHITS) ; 
1715 //      qas.Reset() ; 
1716         rv *= qas.Run(fMakeSDigits.Data(), AliQA::kSDIGITS) ;   
1717 //      qas.Reset() ; 
1718         rv *= qas.Run(fMakeDigits.Data(), AliQA::kDIGITS) ;     
1719 //      qas.Reset() ; 
1720         rv *= qas.Run(fMakeDigitsFromHits.Data(), AliQA::kDIGITS) ; 
1721
1722         return rv ; 
1723 }
1724
1725 //_____________________________________________________________________________
1726 void AliSimulation::ProcessEnvironmentVars()
1727 {
1728 // Extract run number and random generator seed from env variables
1729
1730     AliInfo("Processing environment variables");
1731     
1732     // Random Number seed
1733     
1734     // first check that seed is not already set
1735     if (fSeed == 0) {
1736         if (gSystem->Getenv("CONFIG_SEED")) {
1737                 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
1738         }
1739     } else {
1740         if (gSystem->Getenv("CONFIG_SEED")) {
1741                 AliInfo(Form("Seed for random number generation already set (%d)"
1742                              ": CONFIG_SEED variable ignored!", fSeed));
1743         }
1744     }
1745    
1746     AliInfo(Form("Seed for random number generation = %d ", fSeed)); 
1747
1748     // Run Number
1749     
1750     // first check that run number is not already set
1751     if(fRun < 0) {    
1752         if (gSystem->Getenv("DC_RUN")) {
1753                 fRun = atoi(gSystem->Getenv("DC_RUN"));
1754         }
1755     } else {
1756         if (gSystem->Getenv("DC_RUN")) {
1757                 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
1758         }
1759     }
1760     
1761     AliInfo(Form("Run number = %d", fRun)); 
1762 }