]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
use Terminate function for producing summary plots (Markus)
[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 <TCint.h>
109 #include <TFile.h>
110 #include <TGeoGlobalMagField.h>
111 #include <TGeoManager.h>
112 #include <TObjString.h>
113 #include <TROOT.h>
114 #include <TSystem.h>
115 #include <TVirtualMC.h>
116 #include <TVirtualMCApplication.h>
117 #include <TDatime.h>
118
119 #include "AliAlignObj.h"
120 #include "AliCDBEntry.h"
121 #include "AliCDBManager.h"
122 #include "AliGRPManager.h"
123 #include "AliCDBStorage.h"
124 #include "AliCTPRawData.h"
125 #include "AliCentralTrigger.h"
126 #include "AliCentralTrigger.h"
127 #include "AliCodeTimer.h"
128 #include "AliDAQ.h"
129 #include "AliDigitizer.h"
130 #include "AliESDEvent.h"
131 #include "AliGRPObject.h"
132 #include "AliGenEventHeader.h"
133 #include "AliGenerator.h"
134 #include "AliGeomManager.h"
135 #include "AliHLTSimulation.h"
136 #include "AliHeader.h"
137 #include "AliLego.h"
138 #include "AliLegoGenerator.h"
139 #include "AliLog.h"
140 #include "AliMC.h"
141 #include "AliMagF.h"
142 #include "AliModule.h"
143 #include "AliPDG.h"
144 #include "AliRawReaderDate.h"
145 #include "AliRawReaderFile.h"
146 #include "AliRawReaderRoot.h"
147 #include "AliRun.h"
148 #include "AliRunDigitizer.h"
149 #include "AliRunLoader.h"
150 #include "AliSimulation.h"
151 #include "AliSysInfo.h"
152 #include "AliVertexGenFile.h"
153
154 ClassImp(AliSimulation)
155
156 AliSimulation *AliSimulation::fgInstance = 0;
157 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
158
159 //_____________________________________________________________________________
160 AliSimulation::AliSimulation(const char* configFileName,
161                              const char* name, const char* title) :
162   TNamed(name, title),
163
164   fRunGeneration(kTRUE),
165   fRunSimulation(kTRUE),
166   fLoadAlignFromCDB(kTRUE),
167   fLoadAlObjsListOfDets("ALL"),
168   fMakeSDigits("ALL"),
169   fMakeDigits("ALL"),
170   fTriggerConfig(""),
171   fMakeDigitsFromHits(""),
172   fWriteRawData(""),
173   fRawDataFileName(""),
174   fDeleteIntermediateFiles(kFALSE),
175   fWriteSelRawData(kFALSE),
176   fStopOnError(kFALSE),
177   fNEvents(1),
178   fConfigFileName(configFileName),
179   fGAliceFileName("galice.root"),
180   fEventsPerFile(),
181   fBkgrdFileNames(NULL),
182   fAlignObjArray(NULL),
183   fUseBkgrdVertex(kTRUE),
184   fRegionOfInterest(kFALSE),
185   fCDBUri(""),
186   fQARefUri(""), 
187   fSpecCDBUri(),
188   fRun(-1),
189   fSeed(0),
190   fInitCDBCalled(kFALSE),
191   fInitRunNumberCalled(kFALSE),
192   fSetRunNumberFromDataCalled(kFALSE),
193   fEmbeddingFlag(kFALSE),
194   fLego(NULL),
195   fKey(0),
196   fUseVertexFromCDB(0),
197   fUseMagFieldFromGRP(0),
198   fGRPWriteLocation(Form("local://%s", gSystem->pwd())),
199   fQADetectors("ALL"),                  
200   fQATasks("ALL"),      
201   fRunQA(kTRUE), 
202   fEventSpecie(AliRecoParam::kDefault),
203   fWriteQAExpertData(kTRUE), 
204   fGeometryFile(),
205   fRunHLT("default"),
206   fpHLT(NULL),
207   fWriteGRPEntry(kTRUE)
208 {
209 // create simulation object with default parameters
210   fgInstance = this;
211   SetGAliceFile("galice.root");
212   
213 // for QA
214         AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ; 
215         qam->SetActiveDetectors(fQADetectors) ; 
216         fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ; 
217         qam->SetTasks(fQATasks) ;       
218 }
219
220 //_____________________________________________________________________________
221 AliSimulation::~AliSimulation()
222 {
223 // clean up
224
225   fEventsPerFile.Delete();
226 //  if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
227 //  delete fAlignObjArray; fAlignObjArray=0;
228
229   if (fBkgrdFileNames) {
230     fBkgrdFileNames->Delete();
231     delete fBkgrdFileNames;
232   }
233
234   fSpecCDBUri.Delete();
235   if (fgInstance==this) fgInstance = 0;
236
237   AliQAManager::QAManager()->ShowQA() ; 
238   AliQAManager::Destroy() ;     
239   AliCodeTimer::Instance()->Print();
240 }
241
242
243 //_____________________________________________________________________________
244 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
245 {
246 // set the number of events for one run
247
248   fNEvents = nEvents;
249 }
250
251 //_____________________________________________________________________________
252 void AliSimulation::InitQA()
253 {
254   // activate a default CDB storage
255   // First check if we have any CDB storage set, because it is used 
256   // to retrieve the calibration and alignment constants
257   
258   if (fInitCDBCalled) return;
259   fInitCDBCalled = kTRUE;
260
261   AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ; 
262   qam->SetActiveDetectors(fQADetectors) ; 
263   fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ; 
264   qam->SetTasks(fQATasks) ;
265         if (fWriteQAExpertData)
266     qam->SetWriteExpert() ; 
267   
268   if (qam->IsDefaultStorageSet()) {
269     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
270     AliWarning("Default QA reference storage has been already set !");
271     AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fQARefUri.Data()));
272     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
273     fQARefUri = qam->GetDefaultStorage()->GetURI();
274   } else {
275       if (fQARefUri.Length() > 0) {
276         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
277         AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
278         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
279       } else {
280         fQARefUri="local://$ALICE_ROOT/QARef";
281         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
282         AliWarning("Default QA reference storage not yet set !!!!");
283         AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
284         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
285       }
286     qam->SetDefaultStorage(fQARefUri);
287   }
288 }
289
290 //_____________________________________________________________________________
291 void AliSimulation::InitCDB()
292 {
293 // activate a default CDB storage
294 // First check if we have any CDB storage set, because it is used 
295 // to retrieve the calibration and alignment constants
296
297   if (fInitCDBCalled) return;
298   fInitCDBCalled = kTRUE;
299
300   AliCDBManager* man = AliCDBManager::Instance();
301   if (man->IsDefaultStorageSet())
302   {
303     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
304     AliWarning("Default CDB storage has been already set !");
305     AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
306     AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
307     fCDBUri = man->GetDefaultStorage()->GetURI();
308   }
309   else {
310     if (fCDBUri.Length() > 0) 
311     {
312         AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
313         AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
314         AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
315     } else {
316         fCDBUri="local://$ALICE_ROOT/OCDB";
317         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
318         AliWarning("Default CDB storage not yet set !!!!");
319         AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
320         AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
321                 
322     }
323     man->SetDefaultStorage(fCDBUri);
324   }
325
326   // Now activate the detector specific CDB storage locations
327   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
328     TObject* obj = fSpecCDBUri[i];
329     if (!obj) continue;
330     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
331     AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
332     AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
333     man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
334   }
335       
336 }
337
338 //_____________________________________________________________________________
339 void AliSimulation::InitRunNumber(){
340 // check run number. If not set, set it to 0 !!!!
341   
342   if (fInitRunNumberCalled) return;
343   fInitRunNumberCalled = kTRUE;
344   
345   AliCDBManager* man = AliCDBManager::Instance();
346   if (man->GetRun() >= 0)
347   {
348         AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
349                         "Use external variable DC_RUN or AliSimulation::SetRun()!"));
350   }
351     
352   if(fRun >= 0) {
353         AliDebug(2, Form("Setting CDB run number to: %d",fRun));
354   } else {
355         fRun=0;
356         AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
357                         fRun));
358   }
359   man->SetRun(fRun);
360
361   man->Print();
362
363 }
364
365 //_____________________________________________________________________________
366 void AliSimulation::SetCDBLock() {
367   // Set CDB lock: from now on it is forbidden to reset the run number
368   // or the default storage or to activate any further storage!
369   
370   ULong_t key = AliCDBManager::Instance()->SetLock(1);
371   if (key) fKey = key;
372 }
373
374 //_____________________________________________________________________________
375 void AliSimulation::SetDefaultStorage(const char* uri) {
376   // Store the desired default CDB storage location
377   // Activate it later within the Run() method
378   
379   fCDBUri = uri;
380   
381 }
382
383 //_____________________________________________________________________________
384 void AliSimulation::SetQARefDefaultStorage(const char* uri) {
385   // Store the desired default CDB storage location
386   // Activate it later within the Run() method
387   
388   fQARefUri = uri;
389   AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
390 }
391
392 //_____________________________________________________________________________
393 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
394 // Store a detector-specific CDB storage location
395 // Activate it later within the Run() method
396
397   AliCDBPath aPath(calibType);
398   if(!aPath.IsValid()){
399         AliError(Form("Not a valid path: %s", calibType));
400         return;
401   }
402
403   TObject* obj = fSpecCDBUri.FindObject(calibType);
404   if (obj) fSpecCDBUri.Remove(obj);
405   fSpecCDBUri.Add(new TNamed(calibType, uri));
406
407 }
408
409 //_____________________________________________________________________________
410 void AliSimulation::SetRunNumber(Int_t run)
411 {
412 // sets run number
413 // Activate it later within the Run() method
414
415         fRun = run;
416 }
417
418 //_____________________________________________________________________________
419 void AliSimulation::SetSeed(Int_t seed)
420 {
421 // sets seed number
422 // Activate it later within the Run() method
423
424         fSeed = seed;
425 }
426
427 //_____________________________________________________________________________
428 Bool_t AliSimulation::SetRunNumberFromData()
429 {
430   // Set the CDB manager run number
431   // The run number is retrieved from gAlice
432
433     if (fSetRunNumberFromDataCalled) return kTRUE;
434     fSetRunNumberFromDataCalled = kTRUE;    
435   
436     AliCDBManager* man = AliCDBManager::Instance();
437     Int_t runData = -1, runCDB = -1;
438   
439     AliRunLoader* runLoader = LoadRun("READ");
440     if (!runLoader) return kFALSE;
441     else {
442         runData = runLoader->GetHeader()->GetRun();
443         delete runLoader;
444     }
445   
446     runCDB = man->GetRun();
447     if(runCDB >= 0) {
448         if (runCDB != runData) {
449                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
450                 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
451                 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
452                 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");    
453         }
454         
455     }
456       
457     man->SetRun(runData);
458     fRun = runData;
459     
460     if(man->GetRun() < 0) {
461         AliError("Run number not properly initalized!");
462         return kFALSE;
463     }
464   
465     man->Print();
466     
467     return kTRUE;
468 }
469
470 //_____________________________________________________________________________
471 void AliSimulation::SetConfigFile(const char* fileName)
472 {
473 // set the name of the config file
474
475   fConfigFileName = fileName;
476 }
477
478 //_____________________________________________________________________________
479 void AliSimulation::SetGAliceFile(const char* fileName)
480 {
481 // set the name of the galice file
482 // the path is converted to an absolute one if it is relative
483
484   fGAliceFileName = fileName;
485   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
486     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
487                                                 fGAliceFileName);
488     fGAliceFileName = absFileName;
489     delete[] absFileName;
490   }
491
492   AliDebug(2, Form("galice file name set to %s", fileName));
493 }
494
495 //_____________________________________________________________________________
496 void AliSimulation::SetEventsPerFile(const char* detector, const char* type, 
497                                      Int_t nEvents)
498 {
499 // set the number of events per file for the given detector and data type
500 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
501
502   TNamed* obj = new TNamed(detector, type);
503   obj->SetUniqueID(nEvents);
504   fEventsPerFile.Add(obj);
505 }
506
507 //_____________________________________________________________________________
508 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
509 {
510   // Read the alignment objects from CDB.
511   // Each detector is supposed to have the
512   // alignment objects in DET/Align/Data CDB path.
513   // All the detector objects are then collected,
514   // sorted by geometry level (starting from ALIC) and
515   // then applied to the TGeo geometry.
516   // Finally an overlaps check is performed.
517
518   if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
519     AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
520     return kFALSE;
521   }  
522   
523   // initialize CDB storage, run number, set CDB lock
524   InitCDB();
525 //  if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
526   SetCDBLock();
527     
528   Bool_t delRunLoader = kFALSE;
529   if (!runLoader) {
530     runLoader = LoadRun("READ");
531     if (!runLoader) return kFALSE;
532     delRunLoader = kTRUE;
533   }
534   
535   // Export ideal geometry 
536   if(!IsGeometryFromFile()) AliGeomManager::GetGeometry()->Export("geometry.root");
537
538   // Load alignment data from CDB and apply to geometry through AliGeomManager
539   if(fLoadAlignFromCDB){
540     
541     TString detStr = fLoadAlObjsListOfDets;
542     TString loadAlObjsListOfDets = "";
543     
544     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
545     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
546       AliModule* det = (AliModule*) detArray->At(iDet);
547       if (!det || !det->IsActive()) continue;
548       if (IsSelected(det->GetName(), detStr)) {
549         //add det to list of dets to be aligned from CDB
550         loadAlObjsListOfDets += det->GetName();
551         loadAlObjsListOfDets += " ";
552       }
553     } // end loop over detectors
554     loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
555     AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
556   }else{
557     // Check if the array with alignment objects was
558     // provided by the user. If yes, apply the objects
559     // to the present TGeo geometry
560     if (fAlignObjArray) {
561       if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
562         AliError("The misalignment of one or more volumes failed!"
563                  "Compare the list of simulated detectors and the list of detector alignment data!");
564         if (delRunLoader) delete runLoader;
565         return kFALSE;
566       }
567     }
568   }
569
570   // Update the internal geometry of modules (ITS needs it)
571   TString detStr = fLoadAlObjsListOfDets;
572   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
573   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
574
575     AliModule* det = (AliModule*) detArray->At(iDet);
576     if (!det || !det->IsActive()) continue;
577     if (IsSelected(det->GetName(), detStr)) {
578       det->UpdateInternalGeometry();
579     }
580   } // end loop over detectors
581
582
583   if (delRunLoader) delete runLoader;
584
585   return kTRUE;
586 }
587
588 //_____________________________________________________________________________
589 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
590 {
591 // add a file with background events for merging
592
593   TObjString* fileNameStr = new TObjString(fileName);
594   fileNameStr->SetUniqueID(nSignalPerBkgrd);
595   if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
596   fBkgrdFileNames->Add(fileNameStr);
597 }
598
599 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
600 {
601 // add a file with background events for embeddin
602   MergeWith(fileName, nSignalPerBkgrd);
603   fEmbeddingFlag = kTRUE;
604 }
605
606 //_____________________________________________________________________________
607 Bool_t AliSimulation::Run(Int_t nEvents)
608 {
609 // run the generation, simulation and digitization
610
611  
612   AliCodeTimerAuto("",0)
613   AliSysInfo::AddStamp("Start_Run");
614   
615   // Load run number and seed from environmental vars
616   ProcessEnvironmentVars();
617   AliSysInfo::AddStamp("ProcessEnvironmentVars");
618
619   gRandom->SetSeed(fSeed);
620    
621   if (nEvents > 0) fNEvents = nEvents;
622
623   // create and setup the HLT instance
624   if (!fRunHLT.IsNull() && !CreateHLT()) {
625     if (fStopOnError) return kFALSE;
626     // disable HLT
627     fRunHLT="";
628   }
629   
630   // generation and simulation -> hits
631   if (fRunGeneration) {
632     if (!RunSimulation()) if (fStopOnError) return kFALSE;
633   }
634   AliSysInfo::AddStamp("RunSimulation");
635            
636   // initialize CDB storage from external environment
637   // (either CDB manager or AliSimulation setters),
638   // if not already done in RunSimulation()
639   InitCDB();
640   AliSysInfo::AddStamp("InitCDB");
641   
642   // Set run number in CDBManager from data 
643   // From this point on the run number must be always loaded from data!
644   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
645   
646   // Set CDB lock: from now on it is forbidden to reset the run number
647   // or the default storage or to activate any further storage!
648   SetCDBLock();
649
650   // If RunSimulation was not called, load the geometry and misalign it
651   if (!AliGeomManager::GetGeometry()) {
652     // Initialize the geometry manager
653     AliGeomManager::LoadGeometry("geometry.root");
654   AliSysInfo::AddStamp("GetGeometry");
655
656     
657 //    // Check that the consistency of symbolic names for the activated subdetectors
658 //    // in the geometry loaded by AliGeomManager
659 //    AliRunLoader* runLoader = LoadRun("READ");
660 //    if (!runLoader) return kFALSE;
661 //
662 //    TString detsToBeChecked = "";
663 //    TObjArray* detArray = runLoader->GetAliRun()->Detectors();
664 //    for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
665 //      AliModule* det = (AliModule*) detArray->At(iDet);
666 //      if (!det || !det->IsActive()) continue;
667 //      detsToBeChecked += det->GetName();
668 //      detsToBeChecked += " ";
669 //    } // end loop over detectors
670 //    if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
671     if(!AliGeomManager::CheckSymNamesLUT("ALL"))
672         AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
673         
674     if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
675     // Misalign geometry
676     if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
677   }
678   AliSysInfo::AddStamp("MissalignGeometry");
679
680
681   // hits -> summable digits
682   AliSysInfo::AddStamp("Start_sdigitization");
683   if (!fMakeSDigits.IsNull()) {
684     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
685  
686   }
687   AliSysInfo::AddStamp("Stop_sdigitization");
688   
689   AliSysInfo::AddStamp("Start_digitization");  
690   // summable digits -> digits  
691   if (!fMakeDigits.IsNull()) {
692     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
693       if (fStopOnError) return kFALSE;
694     }
695    }
696   AliSysInfo::AddStamp("Stop_digitization");
697
698   
699   
700   // hits -> digits
701   if (!fMakeDigitsFromHits.IsNull()) {
702     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
703       AliWarning(Form("Merging and direct creation of digits from hits " 
704                  "was selected for some detectors. "
705                  "No merging will be done for the following detectors: %s",
706                  fMakeDigitsFromHits.Data()));
707     }
708     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
709       if (fStopOnError) return kFALSE;
710     }
711   }
712
713   AliSysInfo::AddStamp("Hits2Digits");
714   
715   
716   // digits -> trigger
717   if (!RunTrigger(fTriggerConfig,fMakeDigits)) {
718     if (fStopOnError) return kFALSE;
719   }
720
721   AliSysInfo::AddStamp("RunTrigger");
722   
723   
724   // digits -> raw data
725   if (!fWriteRawData.IsNull()) {
726     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
727                       fDeleteIntermediateFiles,fWriteSelRawData)) {
728       if (fStopOnError) return kFALSE;
729     }
730   }
731
732   AliSysInfo::AddStamp("WriteRaw");
733   
734   // run HLT simulation on simulated digit data if raw data is not
735   // simulated, otherwise its called as part of WriteRawData
736   if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
737     if (!RunHLT()) {
738       if (fStopOnError) return kFALSE;
739     }
740   }
741
742   AliSysInfo::AddStamp("RunHLT");
743   
744   //QA
745         if (fRunQA) {
746                 Bool_t rv = RunQA() ; 
747                 if (!rv)
748                         if (fStopOnError) 
749                                 return kFALSE ;         
750         }
751
752   AliSysInfo::AddStamp("RunQA");
753
754   // Cleanup of CDB manager: cache and active storages!
755   AliCDBManager::Instance()->ClearCache();
756
757   return kTRUE;
758 }
759
760 //_______________________________________________________________________
761 Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
762                      Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
763                      Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
764 {
765   //
766   // Generates lego plots of:
767   //    - radiation length map phi vs theta
768   //    - radiation length map phi vs eta
769   //    - interaction length map
770   //    - g/cm2 length map
771   //
772   //  ntheta    bins in theta, eta
773   //  themin    minimum angle in theta (degrees)
774   //  themax    maximum angle in theta (degrees)
775   //  nphi      bins in phi
776   //  phimin    minimum angle in phi (degrees)
777   //  phimax    maximum angle in phi (degrees)
778   //  rmin      minimum radius
779   //  rmax      maximum radius
780   //  
781   //
782   //  The number of events generated = ntheta*nphi
783   //  run input parameters in macro setup (default="Config.C")
784   //
785   //  Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
786   //Begin_Html
787   /*
788     <img src="picts/AliRunLego1.gif">
789   */
790   //End_Html
791   //Begin_Html
792   /*
793     <img src="picts/AliRunLego2.gif">
794   */
795   //End_Html
796   //Begin_Html
797   /*
798     <img src="picts/AliRunLego3.gif">
799   */
800   //End_Html
801   //
802
803 // run the generation and simulation
804
805   AliCodeTimerAuto("",0)
806
807   // initialize CDB storage and run number from external environment
808   // (either CDB manager or AliSimulation setters)
809   InitCDB();
810   InitRunNumber();
811   SetCDBLock();
812   
813   if (!gAlice) {
814     AliError("no gAlice object. Restart aliroot and try again.");
815     return kFALSE;
816   }
817   if (gAlice->Modules()->GetEntries() > 0) {
818     AliError("gAlice was already run. Restart aliroot and try again.");
819     return kFALSE;
820   }
821
822   AliInfo(Form("initializing gAlice with config file %s",
823           fConfigFileName.Data()));
824
825   // Number of events 
826     if (nev == -1) nev  = nc1 * nc2;
827     
828   // check if initialisation has been done
829   // If runloader has been initialized, set the number of events per file to nc1 * nc2
830     
831   // Set new generator
832   if (!gener) gener  = new AliLegoGenerator();
833   //
834   // Configure Generator
835
836   gener->SetRadiusRange(rmin, rmax);
837   gener->SetZMax(zmax);
838   gener->SetCoor1Range(nc1, c1min, c1max);
839   gener->SetCoor2Range(nc2, c2min, c2max);
840   
841   
842   //Create Lego object  
843   fLego = new AliLego("lego",gener);
844
845   //__________________________________________________________________________
846
847   gAlice->Announce();
848
849   gROOT->LoadMacro(setup);
850   gInterpreter->ProcessLine(gAlice->GetConfigFunction());
851
852   if(AliCDBManager::Instance()->GetRun() >= 0) { 
853     SetRunNumber(AliCDBManager::Instance()->GetRun());
854   } else {
855     AliWarning("Run number not initialized!!");
856   }
857
858   AliRunLoader::Instance()->CdGAFile();
859   
860   AliPDG::AddParticlesToPdgDataBase();  
861   
862   gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
863
864   gAlice->GetMCApp()->Init();
865   
866   
867   //Must be here because some MCs (G4) adds detectors here and not in Config.C
868   gAlice->InitLoaders();
869   AliRunLoader::Instance()->MakeTree("E");
870   
871   //
872   // Save stuff at the beginning of the file to avoid file corruption
873   AliRunLoader::Instance()->CdGAFile();
874   gAlice->Write();
875
876   //Save current generator
877   AliGenerator *gen=gAlice->GetMCApp()->Generator();
878   gAlice->GetMCApp()->ResetGenerator(gener);
879   //Prepare MC for Lego Run
880   gMC->InitLego();
881   
882   //Run Lego Object
883   
884   
885   AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
886   gMC->ProcessRun(nev);
887   
888   // End of this run, close files
889   FinishRun();
890   // Restore current generator
891   gAlice->GetMCApp()->ResetGenerator(gen);
892   // Delete Lego Object
893   delete fLego;
894
895   return kTRUE;
896 }
897
898 //_____________________________________________________________________________
899 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
900 {
901   // run the trigger
902
903   AliCodeTimerAuto("",0)
904
905   // initialize CDB storage from external environment
906   // (either CDB manager or AliSimulation setters),
907   // if not already done in RunSimulation()
908   InitCDB();
909   
910   // Set run number in CDBManager from data 
911   // From this point on the run number must be always loaded from data!
912   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
913   
914   // Set CDB lock: from now on it is forbidden to reset the run number
915   // or the default storage or to activate any further storage!
916   SetCDBLock();
917    
918    AliRunLoader* runLoader = LoadRun("READ");
919    if (!runLoader) return kFALSE;
920    TString trconfiguration = config;
921
922    if (trconfiguration.IsNull()) {
923      if(!fTriggerConfig.IsNull()) {
924        trconfiguration = fTriggerConfig;
925      }
926      else
927        AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
928    }
929
930    runLoader->MakeTree( "GG" );
931    AliCentralTrigger* aCTP = runLoader->GetTrigger();
932    // Load Configuration
933    if (!aCTP->LoadConfiguration( trconfiguration ))
934      return kFALSE;
935
936    // digits -> trigger
937    if( !aCTP->RunTrigger( runLoader , detectors ) ) {
938       if (fStopOnError) {
939         //  delete aCTP;
940         return kFALSE;
941       }
942    }
943
944    delete runLoader;
945
946    return kTRUE;
947 }
948
949 //_____________________________________________________________________________
950 Bool_t AliSimulation::WriteTriggerRawData()
951 {
952   // Writes the CTP (trigger) DDL raw data
953   // Details of the format are given in the
954   // trigger TDR - pages 134 and 135.
955   AliCTPRawData writer;
956   writer.RawData();
957
958   return kTRUE;
959 }
960
961 //_____________________________________________________________________________
962 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
963 {
964 // run the generation and simulation
965
966   AliCodeTimerAuto("",0)
967
968   // initialize CDB storage and run number from external environment
969   // (either CDB manager or AliSimulation setters)
970   AliSysInfo::AddStamp("RunSimulation_Begin");
971   InitCDB();
972   AliSysInfo::AddStamp("RunSimulation_InitCDB");
973   InitRunNumber();
974   SetCDBLock();
975   AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
976   
977   if (!gAlice) {
978     AliError("no gAlice object. Restart aliroot and try again.");
979     return kFALSE;
980   }
981   if (gAlice->Modules()->GetEntries() > 0) {
982     AliError("gAlice was already run. Restart aliroot and try again.");
983     return kFALSE;
984   }
985
986   AliInfo(Form("initializing gAlice with config file %s",
987           fConfigFileName.Data()));
988
989   //
990   // Initialize ALICE Simulation run
991   //
992   gAlice->Announce();
993
994   //
995   // If requested set the mag. field from the GRP entry.
996   // After this the field is loccked and cannot be changed by Config.C
997   if (fUseMagFieldFromGRP) {
998       AliGRPManager grpM;
999       grpM.ReadGRPEntry();
1000       grpM.SetMagField();
1001       AliInfo("Field is locked now. It cannot be changed in Config.C");
1002   }
1003 //
1004 // Execute Config.C
1005   gROOT->LoadMacro(fConfigFileName.Data());
1006   gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1007   AliSysInfo::AddStamp("RunSimulation_Config");
1008
1009 //
1010 // If requested obtain the vertex position and vertex sigma_z from the CDB
1011 // This overwrites the settings from the Config.C  
1012   if (fUseVertexFromCDB) {
1013       Double_t vtxPos[3] = {0., 0., 0.}; 
1014       Double_t vtxSig[3] = {0., 0., 0.};
1015       AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1016       AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1017       if (vertex) {
1018           vertex->GetXYZ(vtxPos);
1019           vertex->GetSigmaXYZ(vtxSig);
1020           AliInfo("Overwriting Config.C vertex settings !");
1021           AliInfo(Form("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
1022                  vtxPos[0], vtxPos[1], vtxPos[2], vtxSig[2]));
1023
1024           AliGenerator *gen = gAlice->GetMCApp()->Generator();
1025           gen->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]);   // vertex position
1026           gen->SetSigmaZ(vtxSig[2]);
1027       }
1028   }
1029   
1030   if(AliCDBManager::Instance()->GetRun() >= 0) { 
1031     AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
1032     AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
1033   } else {
1034         AliWarning("Run number not initialized!!");
1035   }
1036   
1037    AliRunLoader::Instance()->CdGAFile();
1038    
1039
1040    AliPDG::AddParticlesToPdgDataBase();  
1041
1042    gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1043    AliSysInfo::AddStamp("RunSimulation_GetField");
1044    
1045    gAlice->GetMCApp()->Init();
1046    AliSysInfo::AddStamp("RunSimulation_InitMCApp");
1047
1048    //Must be here because some MCs (G4) adds detectors here and not in Config.C
1049    gAlice->InitLoaders();
1050    AliRunLoader::Instance()->MakeTree("E");
1051    AliRunLoader::Instance()->LoadKinematics("RECREATE");
1052    AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1053    AliRunLoader::Instance()->LoadHits("all","RECREATE");
1054    //
1055    // Save stuff at the beginning of the file to avoid file corruption
1056    AliRunLoader::Instance()->CdGAFile();
1057    gAlice->Write();
1058    gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1059    AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1060   //___________________________________________________________________________________________
1061   
1062   AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
1063
1064   // Set run number in CDBManager
1065   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
1066
1067   AliRunLoader* runLoader = AliRunLoader::Instance();
1068   if (!runLoader) {
1069              AliError(Form("gAlice has no run loader object. "
1070                              "Check your config file: %s", fConfigFileName.Data()));
1071              return kFALSE;
1072   }
1073   SetGAliceFile(runLoader->GetFileName());
1074       
1075   // Misalign geometry
1076 #if ROOT_VERSION_CODE < 331527
1077   AliGeomManager::SetGeometry(gGeoManager);
1078   
1079   // Check that the consistency of symbolic names for the activated subdetectors
1080   // in the geometry loaded by AliGeomManager
1081   TString detsToBeChecked = "";
1082   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1083   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1084     AliModule* det = (AliModule*) detArray->At(iDet);
1085     if (!det || !det->IsActive()) continue;
1086     detsToBeChecked += det->GetName();
1087     detsToBeChecked += " ";
1088   } // end loop over detectors
1089   if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
1090     AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
1091   MisalignGeometry(runLoader);
1092   AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
1093 #endif
1094
1095 //   AliRunLoader* runLoader = AliRunLoader::Instance();
1096 //   if (!runLoader) {
1097 //     AliError(Form("gAlice has no run loader object. "
1098 //                   "Check your config file: %s", fConfigFileName.Data()));
1099 //     return kFALSE;
1100 //   }
1101 //   SetGAliceFile(runLoader->GetFileName());
1102
1103   if (!gAlice->GetMCApp()->Generator()) {
1104     AliError(Form("gAlice has no generator object. "
1105                   "Check your config file: %s", fConfigFileName.Data()));
1106     return kFALSE;
1107   }
1108
1109   // Write GRP entry corresponding to the setting found in Cofig.C
1110   if (fWriteGRPEntry)
1111     WriteGRPEntry();
1112   AliSysInfo::AddStamp("RunSimulation_WriteGRP");
1113
1114   if (nEvents <= 0) nEvents = fNEvents;
1115
1116   // get vertex from background file in case of merging
1117   if (fUseBkgrdVertex &&
1118       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
1119     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
1120     const char* fileName = ((TObjString*)
1121                             (fBkgrdFileNames->At(0)))->GetName();
1122     AliInfo(Form("The vertex will be taken from the background "
1123                  "file %s with nSignalPerBackground = %d", 
1124                  fileName, signalPerBkgrd));
1125     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
1126     gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
1127   }
1128
1129   if (!fRunSimulation) {
1130     gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
1131   }
1132
1133   // set the number of events per file for given detectors and data types
1134   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
1135     if (!fEventsPerFile[i]) continue;
1136     const char* detName = fEventsPerFile[i]->GetName();
1137     const char* typeName = fEventsPerFile[i]->GetTitle();
1138     TString loaderName(detName);
1139     loaderName += "Loader";
1140     AliLoader* loader = runLoader->GetLoader(loaderName);
1141     if (!loader) {
1142       AliError(Form("RunSimulation", "no loader for %s found\n"
1143                     "Number of events per file not set for %s %s", 
1144                     detName, typeName, detName));
1145       continue;
1146     }
1147     AliDataLoader* dataLoader = 
1148       loader->GetDataLoader(typeName);
1149     if (!dataLoader) {
1150       AliError(Form("no data loader for %s found\n"
1151                     "Number of events per file not set for %s %s", 
1152                     typeName, detName, typeName));
1153       continue;
1154     }
1155     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
1156     AliDebug(1, Form("number of events per file set to %d for %s %s",
1157                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
1158   }
1159
1160   AliInfo("running gAlice");
1161   AliSysInfo::AddStamp("Start_ProcessRun");
1162
1163   // Create the Root Tree with one branch per detector
1164   //Hits moved to begin event -> now we are crating separate tree for each event
1165
1166   gMC->ProcessRun(nEvents);
1167
1168   // End of this run, close files
1169   if(nEvents>0) FinishRun();
1170
1171   AliSysInfo::AddStamp("Stop_ProcessRun");
1172   delete runLoader;
1173
1174   return kTRUE;
1175 }
1176
1177 //_____________________________________________________________________________
1178 Bool_t AliSimulation::RunSDigitization(const char* detectors)
1179 {
1180 // run the digitization and produce summable digits
1181   static Int_t eventNr=0;
1182   AliCodeTimerAuto("",0) ;
1183
1184   // initialize CDB storage, run number, set CDB lock
1185   InitCDB();
1186   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1187   SetCDBLock();
1188   
1189   AliRunLoader* runLoader = LoadRun();
1190   if (!runLoader) return kFALSE;
1191
1192   TString detStr = detectors;
1193   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1194   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1195     AliModule* det = (AliModule*) detArray->At(iDet);
1196     if (!det || !det->IsActive()) continue;
1197     if (IsSelected(det->GetName(), detStr)) {
1198       AliInfo(Form("creating summable digits for %s", det->GetName()));
1199       AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
1200       det->Hits2SDigits();
1201       AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
1202       AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
1203     }
1204   }
1205
1206   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1207     AliError(Form("the following detectors were not found: %s",
1208                   detStr.Data()));
1209     if (fStopOnError) return kFALSE;
1210   }
1211   eventNr++;
1212   delete runLoader;
1213
1214   return kTRUE;
1215 }
1216
1217
1218 //_____________________________________________________________________________
1219 Bool_t AliSimulation::RunDigitization(const char* detectors, 
1220                                       const char* excludeDetectors)
1221 {
1222 // run the digitization and produce digits from sdigits
1223
1224   AliCodeTimerAuto("",0)
1225
1226   // initialize CDB storage, run number, set CDB lock
1227   InitCDB();
1228   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1229   SetCDBLock();
1230   
1231   delete AliRunLoader::Instance();
1232   delete gAlice;
1233   gAlice = NULL;
1234
1235   Int_t nStreams = 1;
1236   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1237   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1238   AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
1239   // manager->SetEmbeddingFlag(fEmbeddingFlag);
1240   manager->SetInputStream(0, fGAliceFileName.Data());
1241   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1242     const char* fileName = ((TObjString*)
1243                             (fBkgrdFileNames->At(iStream-1)))->GetName();
1244     manager->SetInputStream(iStream, fileName);
1245   }
1246
1247   TString detStr = detectors;
1248   TString detExcl = excludeDetectors;
1249   manager->GetInputStream(0)->ImportgAlice();
1250   AliRunLoader* runLoader = 
1251     AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
1252   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1253   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1254     AliModule* det = (AliModule*) detArray->At(iDet);
1255     if (!det || !det->IsActive()) continue;
1256     if (IsSelected(det->GetName(), detStr) && 
1257         !IsSelected(det->GetName(), detExcl)) {
1258       AliDigitizer* digitizer = det->CreateDigitizer(manager);
1259       
1260       if (!digitizer) {
1261         AliError(Form("no digitizer for %s", det->GetName()));
1262         if (fStopOnError) return kFALSE;
1263       } else {
1264         digitizer->SetRegionOfInterest(fRegionOfInterest);
1265       }
1266     }
1267   }
1268
1269   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1270     AliError(Form("the following detectors were not found: %s", 
1271                   detStr.Data()));
1272     if (fStopOnError) return kFALSE;
1273   }
1274
1275   if (!manager->GetListOfTasks()->IsEmpty()) {
1276     AliInfo("executing digitization");
1277     manager->Exec("");
1278   }
1279
1280   delete manager;
1281
1282   return kTRUE;
1283 }
1284
1285 //_____________________________________________________________________________
1286 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1287 {
1288 // run the digitization and produce digits from hits
1289
1290   AliCodeTimerAuto("",0)
1291
1292   // initialize CDB storage, run number, set CDB lock
1293   InitCDB();
1294   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1295   SetCDBLock();
1296   
1297   AliRunLoader* runLoader = LoadRun("READ");
1298   if (!runLoader) return kFALSE;
1299
1300   TString detStr = detectors;
1301   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1302   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1303     AliModule* det = (AliModule*) detArray->At(iDet);
1304     if (!det || !det->IsActive()) continue;
1305     if (IsSelected(det->GetName(), detStr)) {
1306       AliInfo(Form("creating digits from hits for %s", det->GetName()));
1307       det->Hits2Digits();
1308     }
1309   }
1310
1311   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1312     AliError(Form("the following detectors were not found: %s", 
1313                   detStr.Data()));
1314     if (fStopOnError) return kFALSE;
1315   }
1316
1317   return kTRUE;
1318 }
1319
1320 //_____________________________________________________________________________
1321 Bool_t AliSimulation::WriteRawData(const char* detectors, 
1322                                    const char* fileName,
1323                                    Bool_t deleteIntermediateFiles,
1324                                    Bool_t selrawdata)
1325 {
1326 // convert the digits to raw data
1327 // First DDL raw data files for the given detectors are created.
1328 // If a file name is given, the DDL files are then converted to a DATE file.
1329 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
1330 // afterwards.
1331 // If the file name has the extension ".root", the DATE file is converted
1332 // to a root file.
1333 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1334 // 'selrawdata' flag can be used to enable writing of detectors raw data
1335 // accoring to the trigger cluster.
1336
1337   AliCodeTimerAuto("",0)
1338   AliSysInfo::AddStamp("WriteRawData_Start");
1339   
1340   TString detStr = detectors;
1341   if (!WriteRawFiles(detStr.Data())) {
1342     if (fStopOnError) return kFALSE;
1343   }
1344   AliSysInfo::AddStamp("WriteRawFiles");
1345
1346   // run HLT simulation on simulated DDL raw files
1347   // and produce HLT ddl raw files to be included in date/root file
1348   // bugfix 2009-06-26: the decision whether to write HLT raw data
1349   // is taken in RunHLT. Here HLT always needs to be run in order to
1350   // create HLT digits, unless its switched off. This is due to the
1351   // special placement of the HLT between the generation of DDL files
1352   // and conversion to DATE/Root file.
1353   detStr.ReplaceAll("HLT", "");
1354   if (!fRunHLT.IsNull()) {
1355     if (!RunHLT()) {
1356       if (fStopOnError) return kFALSE;
1357     }
1358   }
1359   AliSysInfo::AddStamp("WriteRawData_RunHLT");
1360
1361   TString dateFileName(fileName);
1362   if (!dateFileName.IsNull()) {
1363     Bool_t rootOutput = dateFileName.EndsWith(".root");
1364     if (rootOutput) dateFileName += ".date";
1365     TString selDateFileName;
1366     if (selrawdata) {
1367       selDateFileName = "selected.";
1368       selDateFileName+= dateFileName;
1369     }
1370     if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1371       if (fStopOnError) return kFALSE;
1372     }
1373     AliSysInfo::AddStamp("ConvertRawFilesToDate");
1374     if (deleteIntermediateFiles) {
1375       AliRunLoader* runLoader = LoadRun("READ");
1376       if (runLoader) for (Int_t iEvent = 0; 
1377                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1378         char command[256];
1379         sprintf(command, "rm -r raw%d", iEvent);
1380         gSystem->Exec(command);
1381       }
1382       delete runLoader;
1383     }
1384
1385     if (rootOutput) {
1386       if (!ConvertDateToRoot(dateFileName, fileName)) {
1387         if (fStopOnError) return kFALSE;
1388       }
1389       AliSysInfo::AddStamp("ConvertDateToRoot");
1390       if (deleteIntermediateFiles) {
1391         gSystem->Unlink(dateFileName);
1392       }
1393       if (selrawdata) {
1394         TString selFileName = "selected.";
1395         selFileName        += fileName;
1396         if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1397           if (fStopOnError) return kFALSE;
1398         }
1399         if (deleteIntermediateFiles) {
1400           gSystem->Unlink(selDateFileName);
1401         }
1402       }
1403     }
1404   }
1405
1406   return kTRUE;
1407 }
1408
1409 //_____________________________________________________________________________
1410 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1411 {
1412 // convert the digits to raw data DDL files
1413
1414   AliCodeTimerAuto("",0)
1415   
1416   AliRunLoader* runLoader = LoadRun("READ");
1417   if (!runLoader) return kFALSE;
1418
1419   // write raw data to DDL files
1420   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1421     AliInfo(Form("processing event %d", iEvent));
1422     runLoader->GetEvent(iEvent);
1423     TString baseDir = gSystem->WorkingDirectory();
1424     char dirName[256];
1425     sprintf(dirName, "raw%d", iEvent);
1426     gSystem->MakeDirectory(dirName);
1427     if (!gSystem->ChangeDirectory(dirName)) {
1428       AliError(Form("couldn't change to directory %s", dirName));
1429       if (fStopOnError) return kFALSE; else continue;
1430     }
1431
1432     ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1433     runNbFile.close();
1434
1435     TString detStr = detectors;
1436     if (IsSelected("HLT", detStr)) {
1437       // Do nothing. "HLT" will be removed from detStr and HLT raw
1438       // data files are generated in RunHLT.
1439     }
1440
1441     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1442     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1443       AliModule* det = (AliModule*) detArray->At(iDet);
1444       if (!det || !det->IsActive()) continue;
1445       if (IsSelected(det->GetName(), detStr)) {
1446         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1447         det->Digits2Raw();
1448       }
1449     }
1450
1451     if (!WriteTriggerRawData())
1452       if (fStopOnError) return kFALSE;
1453
1454     gSystem->ChangeDirectory(baseDir);
1455     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1456       AliError(Form("the following detectors were not found: %s", 
1457                     detStr.Data()));
1458       if (fStopOnError) return kFALSE;
1459     }
1460   }
1461
1462   delete runLoader;
1463   
1464   return kTRUE;
1465 }
1466
1467 //_____________________________________________________________________________
1468 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1469                                             const char* selDateFileName)
1470 {
1471 // convert raw data DDL files to a DATE file with the program "dateStream"
1472 // The second argument is not empty when the user decides to write
1473 // the detectors raw data according to the trigger cluster.
1474
1475   AliCodeTimerAuto("",0)
1476   
1477   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1478   if (!path) {
1479     AliError("the program dateStream was not found");
1480     if (fStopOnError) return kFALSE;
1481   } else {
1482     delete[] path;
1483   }
1484
1485   AliRunLoader* runLoader = LoadRun("READ");
1486   if (!runLoader) return kFALSE;
1487
1488   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1489   Bool_t selrawdata = kFALSE;
1490   if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1491
1492   char command[256];
1493   // Note the option -s. It is used in order to avoid
1494   // the generation of SOR/EOR events.
1495   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1496           dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1497   FILE* pipe = gSystem->OpenPipe(command, "w");
1498
1499   if (!pipe) {
1500     AliError(Form("Cannot execute command: %s",command));
1501     return kFALSE;
1502   }
1503
1504   Int_t selEvents = 0;
1505   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1506
1507     UInt_t detectorPattern = 0;
1508     runLoader->GetEvent(iEvent);
1509     if (!runLoader->LoadTrigger()) {
1510       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1511       detectorPattern = aCTP->GetClusterMask();
1512       // Check if the event was triggered by CTP
1513       if (selrawdata) {
1514         if (aCTP->GetClassMask()) selEvents++;
1515       }
1516     }
1517     else {
1518       AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1519       if (selrawdata) {
1520         AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1521         selrawdata = kFALSE;
1522       }
1523     }
1524
1525     fprintf(pipe, "GDC DetectorPattern %u\n", detectorPattern);
1526     Float_t ldc = 0;
1527     Int_t prevLDC = -1;
1528
1529     // loop over detectors and DDLs
1530     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1531       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1532
1533         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1534         Int_t ldcID = Int_t(ldc + 0.0001);
1535         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1536
1537         char rawFileName[256];
1538         sprintf(rawFileName, "raw%d/%s", 
1539                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1540
1541         // check existence and size of raw data file
1542         FILE* file = fopen(rawFileName, "rb");
1543         if (!file) continue;
1544         fseek(file, 0, SEEK_END);
1545         unsigned long size = ftell(file);
1546         fclose(file);
1547         if (!size) continue;
1548
1549         if (ldcID != prevLDC) {
1550           fprintf(pipe, " LDC Id %d\n", ldcID);
1551           prevLDC = ldcID;
1552         }
1553         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1554       }
1555     }
1556   }
1557
1558   Int_t result = gSystem->ClosePipe(pipe);
1559
1560   if (!(selrawdata && selEvents > 0)) {
1561     delete runLoader;
1562     return (result == 0);
1563   }
1564
1565   AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1566   
1567   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1568           selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1569   FILE* pipe2 = gSystem->OpenPipe(command, "w");
1570
1571   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1572
1573     // Get the trigger decision and cluster
1574     UInt_t detectorPattern = 0;
1575     TString detClust;
1576     runLoader->GetEvent(iEvent);
1577     if (!runLoader->LoadTrigger()) {
1578       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1579       if (aCTP->GetClassMask() == 0) continue;
1580       detectorPattern = aCTP->GetClusterMask();
1581       detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1582       AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1583     }
1584
1585     fprintf(pipe2, "GDC DetectorPattern %u\n", detectorPattern);
1586     Float_t ldc = 0;
1587     Int_t prevLDC = -1;
1588
1589     // loop over detectors and DDLs
1590     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1591       // Write only raw data from detectors that
1592       // are contained in the trigger cluster(s)
1593       if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1594
1595       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1596
1597         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1598         Int_t ldcID = Int_t(ldc + 0.0001);
1599         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1600
1601         char rawFileName[256];
1602         sprintf(rawFileName, "raw%d/%s", 
1603                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1604
1605         // check existence and size of raw data file
1606         FILE* file = fopen(rawFileName, "rb");
1607         if (!file) continue;
1608         fseek(file, 0, SEEK_END);
1609         unsigned long size = ftell(file);
1610         fclose(file);
1611         if (!size) continue;
1612
1613         if (ldcID != prevLDC) {
1614           fprintf(pipe2, " LDC Id %d\n", ldcID);
1615           prevLDC = ldcID;
1616         }
1617         fprintf(pipe2, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1618       }
1619     }
1620   }
1621
1622   Int_t result2 = gSystem->ClosePipe(pipe2);
1623
1624   delete runLoader;
1625   return ((result == 0) && (result2 == 0));
1626 }
1627
1628 //_____________________________________________________________________________
1629 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1630                                         const char* rootFileName)
1631 {
1632 // convert a DATE file to a root file with the program "alimdc"
1633
1634   // ALIMDC setup
1635   const Int_t kDBSize = 2000000000;
1636   const Int_t kTagDBSize = 1000000000;
1637   const Bool_t kFilter = kFALSE;
1638   const Int_t kCompression = 1;
1639
1640   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1641   if (!path) {
1642     AliError("the program alimdc was not found");
1643     if (fStopOnError) return kFALSE;
1644   } else {
1645     delete[] path;
1646   }
1647
1648   AliInfo(Form("converting DATE file %s to root file %s", 
1649                dateFileName, rootFileName));
1650
1651   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1652   const char* tagDBFS    = "/tmp/mdc1/tags";
1653
1654   // User defined file system locations
1655   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1656     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1657   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1658     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1659   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1660     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1661
1662   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1663   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1664   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1665
1666   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1667   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1668   gSystem->Exec(Form("mkdir %s",tagDBFS));
1669
1670   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1671                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1672   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1673
1674   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1675   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1676   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1677
1678   return (result == 0);
1679 }
1680
1681
1682 //_____________________________________________________________________________
1683 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1684 {
1685 // delete existing run loaders, open a new one and load gAlice
1686
1687   delete AliRunLoader::Instance();
1688   AliRunLoader* runLoader = 
1689     AliRunLoader::Open(fGAliceFileName.Data(), 
1690                        AliConfig::GetDefaultEventFolderName(), mode);
1691   if (!runLoader) {
1692     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1693     return NULL;
1694   }
1695   runLoader->LoadgAlice();
1696   runLoader->LoadHeader();
1697   gAlice = runLoader->GetAliRun();
1698   if (!gAlice) {
1699     AliError(Form("no gAlice object found in file %s", 
1700                   fGAliceFileName.Data()));
1701     return NULL;
1702   }
1703   return runLoader;
1704 }
1705
1706 //_____________________________________________________________________________
1707 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1708 {
1709 // get or calculate the number of signal events per background event
1710
1711   if (!fBkgrdFileNames) return 1;
1712   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1713   if (nBkgrdFiles == 0) return 1;
1714
1715   // get the number of signal events
1716   if (nEvents <= 0) {
1717     AliRunLoader* runLoader = 
1718         AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1719     if (!runLoader) return 1;
1720     
1721     nEvents = runLoader->GetNumberOfEvents();
1722     delete runLoader;
1723   }
1724
1725   Int_t result = 0;
1726   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1727     // get the number of background events
1728     const char* fileName = ((TObjString*)
1729                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1730     AliRunLoader* runLoader =
1731       AliRunLoader::Open(fileName, "BKGRD");
1732     if (!runLoader) continue;
1733     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1734     delete runLoader;
1735   
1736     // get or calculate the number of signal per background events
1737     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1738     if (nSignalPerBkgrd <= 0) {
1739       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1740     } else if (result && (result != nSignalPerBkgrd)) {
1741       AliInfo(Form("the number of signal events per background event "
1742                    "will be changed from %d to %d for stream %d", 
1743                    nSignalPerBkgrd, result, iBkgrdFile+1));
1744       nSignalPerBkgrd = result;
1745     }
1746
1747     if (!result) result = nSignalPerBkgrd;
1748     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1749       AliWarning(Form("not enough background events (%d) for %d signal events "
1750                       "using %d signal per background events for stream %d",
1751                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1752     }
1753   }
1754
1755   return result;
1756 }
1757
1758 //_____________________________________________________________________________
1759 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1760 {
1761 // check whether detName is contained in detectors
1762 // if yes, it is removed from detectors
1763
1764   // check if all detectors are selected
1765   if ((detectors.CompareTo("ALL") == 0) ||
1766       detectors.BeginsWith("ALL ") ||
1767       detectors.EndsWith(" ALL") ||
1768       detectors.Contains(" ALL ")) {
1769     detectors = "ALL";
1770     return kTRUE;
1771   }
1772
1773   // search for the given detector
1774   Bool_t result = kFALSE;
1775   if ((detectors.CompareTo(detName) == 0) ||
1776       detectors.BeginsWith(detName+" ") ||
1777       detectors.EndsWith(" "+detName) ||
1778       detectors.Contains(" "+detName+" ")) {
1779     detectors.ReplaceAll(detName, "");
1780     result = kTRUE;
1781   }
1782
1783   // clean up the detectors string
1784   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1785   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1786   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1787
1788   return result;
1789 }
1790
1791 //_____________________________________________________________________________
1792 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName) 
1793 {
1794 //
1795 // Steering routine  to convert raw data in directory rawDirectory/ to fake SDigits. 
1796 // These can be used for embedding of MC tracks into RAW data using the standard 
1797 // merging procedure.
1798 //
1799 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1800 //
1801   if (!gAlice) {
1802     AliError("no gAlice object. Restart aliroot and try again.");
1803     return kFALSE;
1804   }
1805   if (gAlice->Modules()->GetEntries() > 0) {
1806     AliError("gAlice was already run. Restart aliroot and try again.");
1807     return kFALSE;
1808   }
1809   
1810   AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1811   
1812   gAlice->Announce();
1813   
1814   gROOT->LoadMacro(fConfigFileName.Data());
1815   gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1816   
1817   if(AliCDBManager::Instance()->GetRun() >= 0) { 
1818         SetRunNumber(AliCDBManager::Instance()->GetRun());
1819   } else {
1820         AliWarning("Run number not initialized!!");
1821   }
1822   
1823    AliRunLoader::Instance()->CdGAFile();
1824     
1825    AliPDG::AddParticlesToPdgDataBase();  
1826
1827    gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1828    
1829    gAlice->GetMCApp()->Init();
1830
1831    //Must be here because some MCs (G4) adds detectors here and not in Config.C
1832    gAlice->InitLoaders();
1833    AliRunLoader::Instance()->MakeTree("E");
1834    AliRunLoader::Instance()->LoadKinematics("RECREATE");
1835    AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1836    AliRunLoader::Instance()->LoadHits("all","RECREATE");
1837    //
1838    // Save stuff at the beginning of the file to avoid file corruption
1839    AliRunLoader::Instance()->CdGAFile();
1840    gAlice->Write();
1841 //
1842 //  Initialize CDB     
1843     InitCDB();
1844     //AliCDBManager* man = AliCDBManager::Instance();
1845     //man->SetRun(0); // Should this come from rawdata header ?
1846     
1847     Int_t iDet;
1848     //
1849     // Get the runloader
1850     AliRunLoader* runLoader = AliRunLoader::Instance();
1851     //
1852     // Open esd file if available
1853     TFile* esdFile = TFile::Open(esdFileName);
1854     TTree* treeESD = 0;
1855     AliESDEvent* esd = new AliESDEvent();
1856     esdFile->GetObject("esdTree", treeESD);
1857     if (treeESD) esd->ReadFromTree(treeESD);
1858
1859     //
1860     // Create the RawReader
1861     TString fileName(rawDirectory);
1862     AliRawReader* rawReader = 0x0;
1863     if (fileName.EndsWith("/")) {
1864       rawReader = new AliRawReaderFile(fileName);
1865     } else if (fileName.EndsWith(".root")) {
1866       rawReader = new AliRawReaderRoot(fileName);
1867     } else if (!fileName.IsNull()) {
1868       rawReader = new AliRawReaderDate(fileName);
1869     }
1870 //     if (!fEquipIdMap.IsNull() && fRawReader)
1871 //       fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1872     //
1873     // Get list of detectors
1874     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1875     //
1876     // Get Header
1877     AliHeader* header = runLoader->GetHeader();
1878     //
1879     TString detStr = fMakeSDigits;
1880     // Event loop
1881     Int_t nev = 0;
1882     while(kTRUE) {
1883         if (!(rawReader->NextEvent())) break;
1884         //
1885         // Detector loop
1886         for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1887             AliModule* det = (AliModule*) detArray->At(iDet);
1888             if (!det || !det->IsActive()) continue;
1889             if (IsSelected(det->GetName(), detStr)) {
1890               AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1891               det->Raw2SDigits(rawReader);
1892               rawReader->Reset();
1893             }
1894         } // detectors
1895
1896
1897         //
1898         //  If ESD information available obtain reconstructed vertex and store in header.
1899         if (treeESD) {
1900             treeESD->GetEvent(nev);
1901             const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1902             Double_t position[3];
1903             esdVertex->GetXYZ(position);
1904             AliGenEventHeader* mcHeader = new  AliGenEventHeader("ESD");
1905             TArrayF mcV;
1906             mcV.Set(3);
1907             for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1908             mcHeader->SetPrimaryVertex(mcV);
1909             header->Reset(0,nev);
1910             header->SetGenEventHeader(mcHeader);
1911             printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1912         }
1913         nev++;
1914 //
1915 //      Finish the event
1916         runLoader->TreeE()->Fill();
1917         runLoader->SetNextEvent();
1918     } // events
1919  
1920     delete rawReader;
1921 //
1922 //  Finish the run 
1923     runLoader->CdGAFile();
1924     runLoader->WriteHeader("OVERWRITE");
1925     runLoader->WriteRunLoader();
1926
1927     return kTRUE;
1928 }
1929
1930 //_____________________________________________________________________________
1931 void AliSimulation::FinishRun()
1932 {
1933   //
1934   // Called at the end of the run.
1935   //
1936
1937   if(IsLegoRun()) 
1938    {
1939     AliDebug(1, "Finish Lego");
1940     AliRunLoader::Instance()->CdGAFile();
1941     fLego->FinishRun();
1942    }
1943   
1944   // Clean detector information
1945   TIter next(gAlice->Modules());
1946   AliModule *detector;
1947   while((detector = dynamic_cast<AliModule*>(next()))) {
1948     AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
1949     detector->FinishRun();
1950   }
1951   
1952   AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
1953   AliRunLoader::Instance()->WriteHeader("OVERWRITE");
1954
1955   // Write AliRun info and all detectors parameters
1956   AliRunLoader::Instance()->CdGAFile();
1957   gAlice->Write(0,TObject::kOverwrite);//write AliRun
1958   AliRunLoader::Instance()->Write(0,TObject::kOverwrite);//write RunLoader itself
1959   
1960   if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();  
1961   AliRunLoader::Instance()->Synchronize();
1962 }
1963
1964 //_____________________________________________________________________________
1965 Int_t AliSimulation::GetDetIndex(const char* detector)
1966 {
1967   // return the detector index corresponding to detector
1968   Int_t index = -1 ; 
1969   for (index = 0; index < fgkNDetectors ; index++) {
1970     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1971           break ; 
1972   }     
1973   return index ; 
1974 }
1975
1976 //_____________________________________________________________________________
1977 Bool_t AliSimulation::CreateHLT()
1978 {
1979   // Init the HLT simulation.
1980   // The function  loads the library and creates the instance of AliHLTSimulation.
1981   // the main reason for the decoupled creation is to set the transient OCDB
1982   // objects before the OCDB is locked
1983
1984   // load the library dynamically
1985   gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1986
1987   // check for the library version
1988   AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
1989   if (!fctVersion) {
1990     AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
1991     return kFALSE;
1992   }
1993   if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
1994     AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
1995   }
1996
1997   // print compile info
1998   typedef void (*CompileInfo)( const char*& date, const char*& time);
1999   CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
2000   if (fctInfo) {
2001     const char* date="";
2002     const char* time="";
2003     (*fctInfo)(date, time);
2004     if (!date) date="unknown";
2005     if (!time) time="unknown";
2006     AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
2007   } else {
2008     AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
2009   }
2010
2011   // create instance of the HLT simulation
2012   AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
2013   if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
2014     AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
2015     return kFALSE;    
2016   }
2017
2018   TString specObjects;
2019   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
2020     if (specObjects.Length()>0) specObjects+=" ";
2021     specObjects+=fSpecCDBUri[i]->GetName();
2022   }
2023
2024   AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
2025   if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
2026     AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
2027   }
2028
2029   return kTRUE;
2030 }
2031
2032 //_____________________________________________________________________________
2033 Bool_t AliSimulation::RunHLT()
2034 {
2035   // Run the HLT simulation
2036   // HLT simulation is implemented in HLT/sim/AliHLTSimulation
2037   // Disabled if fRunHLT is empty, default vaule is "default".
2038   // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
2039   // The default simulation depends on the HLT component libraries and their
2040   // corresponding agents which define components and chains to run. See
2041   // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
2042   // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
2043   //
2044   // The libraries to be loaded can be specified as an option.
2045   // <pre>
2046   // AliSimulation sim;
2047   // sim.SetRunHLT("libAliHLTSample.so");
2048   // </pre>
2049   // will only load <tt>libAliHLTSample.so</tt>
2050
2051   // Other available options:
2052   // \li loglevel=<i>level</i> <br>
2053   //     logging level for this processing
2054   // \li alilog=off
2055   //     disable redirection of log messages to AliLog class
2056   // \li config=<i>macro</i>
2057   //     configuration macro
2058   // \li chains=<i>configuration</i>
2059   //     comma separated list of configurations to be run during simulation
2060   // \li rawfile=<i>file</i>
2061   //     source for the RawReader to be created, the default is <i>./</i> if
2062   //     raw data is simulated
2063
2064   int iResult=0;
2065
2066   if (!fpHLT && !CreateHLT()) {
2067     return kFALSE;
2068   }
2069   AliHLTSimulation* pHLT=fpHLT;
2070
2071   AliRunLoader* pRunLoader = LoadRun("READ");
2072   if (!pRunLoader) return kFALSE;
2073
2074   // initialize CDB storage, run number, set CDB lock
2075   // thats for the case of running HLT simulation without all the other steps
2076   // multiple calls are handled by the function, so we can just call
2077   InitCDB();
2078   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
2079   SetCDBLock();
2080   
2081   // init the HLT simulation
2082   TString options;
2083   if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2084   TString detStr = fWriteRawData;
2085   if (!IsSelected("HLT", detStr)) {
2086     options+=" writerawfiles=";
2087   } else {
2088     options+=" writerawfiles=HLT";
2089   }
2090
2091   if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2092     // as a matter of fact, HLT will run reconstruction and needs the RawReader
2093     // in order to get detector data. By default, RawReaderFile is used to read
2094     // the already simulated ddl files. Date and Root files from the raw data
2095     // are generated after the HLT simulation.
2096     options+=" rawfile=./";
2097   }
2098
2099   AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2100   if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2101     AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2102   } else {
2103     // run the HLT simulation
2104     AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2105     if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2106       AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2107     }
2108   }
2109
2110   // delete the instance
2111   AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2112   if (fctDelete==NULL || fctDelete(pHLT)<0) {
2113     AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2114   }
2115   pHLT=NULL;
2116
2117   return iResult>=0?kTRUE:kFALSE;
2118 }
2119
2120 //_____________________________________________________________________________
2121 Bool_t AliSimulation::RunQA()
2122 {
2123         // run the QA on summable hits, digits or digits
2124         
2125   if(!gAlice) return kFALSE;
2126         AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2127
2128         TString detectorsw("") ;  
2129         Bool_t rv = kTRUE ; 
2130   AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2131         detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ; 
2132         if ( detectorsw.IsNull() ) 
2133                 rv = kFALSE ; 
2134         return rv ; 
2135 }
2136
2137 //_____________________________________________________________________________
2138 Bool_t AliSimulation::SetRunQA(TString detAndAction) 
2139 {
2140         // Allows to run QA for a selected set of detectors
2141         // and a selected set of tasks among HITS, SDIGITS and DIGITS
2142         // all selected detectors run the same selected tasks
2143         
2144         if (!detAndAction.Contains(":")) {
2145                 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2146                 fRunQA = kFALSE ;
2147                 return kFALSE ;                 
2148         }
2149         Int_t colon = detAndAction.Index(":") ; 
2150         fQADetectors = detAndAction(0, colon) ; 
2151         if (fQADetectors.Contains("ALL") ){
2152     TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ; 
2153     Int_t minus = fQADetectors.Last('-') ; 
2154     TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ; 
2155     TString toRemove("") ;
2156     while (minus >= 0) {
2157       toRemove = fQADetectors(minus+1, fQADetectors.Length()) ; 
2158       toRemove = toRemove.Strip() ; 
2159       toKeep.ReplaceAll(toRemove, "") ; 
2160       fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ; 
2161       minus = fQADetectors.Last('-') ; 
2162     }
2163     fQADetectors = toKeep ; 
2164   }
2165   fQATasks   = detAndAction(colon+1, detAndAction.Sizeof() ) ; 
2166         if (fQATasks.Contains("ALL") ) {
2167                 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ; 
2168         } else {
2169                 fQATasks.ToUpper() ; 
2170                 TString tempo("") ; 
2171                 if ( fQATasks.Contains("HIT") ) 
2172                         tempo = Form("%d ", AliQAv1::kHITS) ; 
2173                 if ( fQATasks.Contains("SDIGIT") ) 
2174                         tempo += Form("%d ", AliQAv1::kSDIGITS) ; 
2175                 if ( fQATasks.Contains("DIGIT") ) 
2176                         tempo += Form("%d ", AliQAv1::kDIGITS) ; 
2177                 fQATasks = tempo ; 
2178                 if (fQATasks.IsNull()) {
2179                         AliInfo("No QA requested\n")  ;
2180                         fRunQA = kFALSE ;
2181                         return kTRUE ; 
2182                 }
2183         }       
2184         TString tempo(fQATasks) ; 
2185     tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS))  ;
2186     tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;  
2187     tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;    
2188         AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;  
2189         fRunQA = kTRUE ;
2190         AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ; 
2191         AliQAManager::QAManager()->SetTasks(fQATasks) ; 
2192   for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) 
2193     AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2194   
2195         return kTRUE; 
2196
2197
2198 //_____________________________________________________________________________
2199 void AliSimulation::ProcessEnvironmentVars()
2200 {
2201 // Extract run number and random generator seed from env variables
2202
2203     AliInfo("Processing environment variables");
2204     
2205     // Random Number seed
2206     
2207     // first check that seed is not already set
2208     if (fSeed == 0) {
2209         if (gSystem->Getenv("CONFIG_SEED")) {
2210                 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2211         }
2212     } else {
2213         if (gSystem->Getenv("CONFIG_SEED")) {
2214                 AliInfo(Form("Seed for random number generation already set (%d)"
2215                              ": CONFIG_SEED variable ignored!", fSeed));
2216         }
2217     }
2218    
2219     AliInfo(Form("Seed for random number generation = %d ", fSeed)); 
2220
2221     // Run Number
2222     
2223     // first check that run number is not already set
2224     if(fRun < 0) {    
2225         if (gSystem->Getenv("DC_RUN")) {
2226                 fRun = atoi(gSystem->Getenv("DC_RUN"));
2227         }
2228     } else {
2229         if (gSystem->Getenv("DC_RUN")) {
2230                 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2231         }
2232     }
2233     
2234     AliInfo(Form("Run number = %d", fRun)); 
2235 }
2236
2237 //---------------------------------------------------------------------
2238 void AliSimulation::WriteGRPEntry()
2239 {
2240   // Get the necessary information from galice (generator, trigger etc) and
2241   // write a GRP entry corresponding to the settings in the Config.C used
2242   // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2243
2244
2245   AliInfo("Writing global run parameters entry into the OCDB");
2246
2247   AliGRPObject* grpObj = new AliGRPObject();
2248
2249   grpObj->SetRunType("PHYSICS");
2250   grpObj->SetTimeStart(0);
2251   TDatime curtime;
2252   grpObj->SetTimeStart(0);
2253   grpObj->SetTimeEnd(curtime.Convert()); 
2254   grpObj->SetBeamEnergyIsSqrtSHalfGeV(); // new format of GRP: store sqrt(s)/2 in GeV
2255
2256   const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2257   if (gen) {
2258     grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
2259     TString projectile;
2260     Int_t a,z;
2261     gen->GetProjectile(projectile,a,z);
2262     TString target;
2263     gen->GetTarget(target,a,z);
2264     TString beamType = projectile + "-" + target;
2265     beamType.ReplaceAll(" ","");
2266     if (!beamType.CompareTo("-")) {
2267       grpObj->SetBeamType("UNKNOWN");
2268     }
2269     else {
2270       grpObj->SetBeamType(beamType);
2271       // Heavy ion run, the event specie is set to kHighMult
2272       fEventSpecie = AliRecoParam::kHighMult;
2273       if ((strcmp(beamType,"p-p") == 0) ||
2274           (strcmp(beamType,"p-")  == 0) ||
2275           (strcmp(beamType,"-p")  == 0) ||
2276           (strcmp(beamType,"P-P") == 0) ||
2277           (strcmp(beamType,"P-")  == 0) ||
2278           (strcmp(beamType,"-P")  == 0)) {
2279         // Proton run, the event specie is set to kLowMult
2280         fEventSpecie = AliRecoParam::kLowMult;
2281       } 
2282     }
2283   } else {
2284     AliWarning("Unknown beam type and energy! Setting energy to 0");
2285     grpObj->SetBeamEnergy(0);
2286     grpObj->SetBeamType("UNKNOWN");
2287   }
2288
2289   UInt_t detectorPattern  = 0;
2290   Int_t nDets = 0;
2291   TObjArray *detArray = gAlice->Detectors();
2292   for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2293     if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2294       detectorPattern |= (1 << iDet);
2295       nDets++;
2296     }
2297   }
2298   // CTP
2299   if (!fTriggerConfig.IsNull())
2300     detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2301
2302   // HLT
2303   if (!fRunHLT.IsNull())
2304     detectorPattern |= (1 << AliDAQ::kHLTId);
2305
2306   grpObj->SetNumberOfDetectors((Char_t)nDets);
2307   grpObj->SetDetectorMask((Int_t)detectorPattern);
2308   grpObj->SetLHCPeriod("LHC08c");
2309   grpObj->SetLHCState("STABLE_BEAMS");
2310   //
2311   AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2312   Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2313
2314   Float_t factorSol     = field ? field->GetFactorSol() : 0;
2315   Float_t currentSol    = TMath::Abs(factorSol)>1E-6 ? 
2316     TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
2317   //
2318   Float_t factorDip     = field ? field->GetFactorDip() : 0;
2319   Float_t currentDip    = 6000.*TMath::Abs(factorDip);
2320   //
2321   grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
2322   grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);  
2323   grpObj->SetL3Polarity(factorSol>0 ? 0:1);  
2324   grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
2325   grpObj->SetUniformBMap(field->IsUniform());            // for special MC with k5kGUniform map
2326   grpObj->SetPolarityConventionLHC();                    // LHC convention +/+ current -> -/- field main components
2327   //
2328   grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2329   
2330   //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2331
2332   // Now store the entry in OCDB
2333   AliCDBManager* man = AliCDBManager::Instance();
2334   
2335   man->SetLock(0, fKey);
2336   
2337   AliCDBStorage* sto = man->GetStorage(fGRPWriteLocation.Data());
2338   
2339
2340   AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2341   AliCDBMetaData *metadata= new AliCDBMetaData();
2342
2343   metadata->SetResponsible("alice-off@cern.ch");
2344   metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2345  
2346   sto->Put(grpObj,id,metadata);
2347   man->SetLock(1, fKey);
2348 }
2349
2350