]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliSimulation.cxx
Adding class to store QA thresholds (Barthelemy von Haller).
[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 //    // Check that the consistency of symbolic names for the activated subdetectors
656 //    // in the geometry loaded by AliGeomManager
657 //    AliRunLoader* runLoader = LoadRun("READ");
658 //    if (!runLoader) return kFALSE;
659 //
660 //    TString detsToBeChecked = "";
661 //    TObjArray* detArray = runLoader->GetAliRun()->Detectors();
662 //    for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
663 //      AliModule* det = (AliModule*) detArray->At(iDet);
664 //      if (!det || !det->IsActive()) continue;
665 //      detsToBeChecked += det->GetName();
666 //      detsToBeChecked += " ";
667 //    } // end loop over detectors
668 //    if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
669     if(!AliGeomManager::CheckSymNamesLUT("ALL"))
670         AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
671         
672     if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
673     // Misalign geometry
674     if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
675   }
676   AliSysInfo::AddStamp("MissalignGeometry");
677
678
679   // hits -> summable digits
680   AliSysInfo::AddStamp("Start_sdigitization");
681   if (!fMakeSDigits.IsNull()) {
682     if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
683  
684   }
685   AliSysInfo::AddStamp("Stop_sdigitization");
686   
687   AliSysInfo::AddStamp("Start_digitization");  
688   // summable digits -> digits  
689   if (!fMakeDigits.IsNull()) {
690     if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
691       if (fStopOnError) return kFALSE;
692     }
693    }
694   AliSysInfo::AddStamp("Stop_digitization");
695
696   
697   
698   // hits -> digits
699   if (!fMakeDigitsFromHits.IsNull()) {
700     if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
701       AliWarning(Form("Merging and direct creation of digits from hits " 
702                  "was selected for some detectors. "
703                  "No merging will be done for the following detectors: %s",
704                  fMakeDigitsFromHits.Data()));
705     }
706     if (!RunHitsDigitization(fMakeDigitsFromHits)) {
707       if (fStopOnError) return kFALSE;
708     }
709   }
710
711   AliSysInfo::AddStamp("Hits2Digits");
712   
713   
714   // digits -> trigger
715   if (!RunTrigger(fTriggerConfig,fMakeDigits)) {
716     if (fStopOnError) return kFALSE;
717   }
718
719   AliSysInfo::AddStamp("RunTrigger");
720   
721   
722   // digits -> raw data
723   if (!fWriteRawData.IsNull()) {
724     if (!WriteRawData(fWriteRawData, fRawDataFileName, 
725                       fDeleteIntermediateFiles,fWriteSelRawData)) {
726       if (fStopOnError) return kFALSE;
727     }
728   }
729
730   AliSysInfo::AddStamp("WriteRaw");
731   
732   // run HLT simulation on simulated digit data if raw data is not
733   // simulated, otherwise its called as part of WriteRawData
734   if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
735     if (!RunHLT()) {
736       if (fStopOnError) return kFALSE;
737     }
738   }
739
740   AliSysInfo::AddStamp("RunHLT");
741   
742   //QA
743         if (fRunQA) {
744                 Bool_t rv = RunQA() ; 
745                 if (!rv)
746                         if (fStopOnError) 
747                                 return kFALSE ;         
748         }
749
750   AliSysInfo::AddStamp("RunQA");
751
752   // Cleanup of CDB manager: cache and active storages!
753   AliCDBManager::Instance()->ClearCache();
754
755   return kTRUE;
756 }
757
758 //_______________________________________________________________________
759 Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
760                      Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
761                      Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
762 {
763   //
764   // Generates lego plots of:
765   //    - radiation length map phi vs theta
766   //    - radiation length map phi vs eta
767   //    - interaction length map
768   //    - g/cm2 length map
769   //
770   //  ntheta    bins in theta, eta
771   //  themin    minimum angle in theta (degrees)
772   //  themax    maximum angle in theta (degrees)
773   //  nphi      bins in phi
774   //  phimin    minimum angle in phi (degrees)
775   //  phimax    maximum angle in phi (degrees)
776   //  rmin      minimum radius
777   //  rmax      maximum radius
778   //  
779   //
780   //  The number of events generated = ntheta*nphi
781   //  run input parameters in macro setup (default="Config.C")
782   //
783   //  Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
784   //Begin_Html
785   /*
786     <img src="picts/AliRunLego1.gif">
787   */
788   //End_Html
789   //Begin_Html
790   /*
791     <img src="picts/AliRunLego2.gif">
792   */
793   //End_Html
794   //Begin_Html
795   /*
796     <img src="picts/AliRunLego3.gif">
797   */
798   //End_Html
799   //
800
801 // run the generation and simulation
802
803   AliCodeTimerAuto("",0)
804
805   // initialize CDB storage and run number from external environment
806   // (either CDB manager or AliSimulation setters)
807   InitCDB();
808   InitRunNumber();
809   SetCDBLock();
810   
811   if (!gAlice) {
812     AliError("no gAlice object. Restart aliroot and try again.");
813     return kFALSE;
814   }
815   if (gAlice->Modules()->GetEntries() > 0) {
816     AliError("gAlice was already run. Restart aliroot and try again.");
817     return kFALSE;
818   }
819
820   AliInfo(Form("initializing gAlice with config file %s",
821           fConfigFileName.Data()));
822
823   // Number of events 
824     if (nev == -1) nev  = nc1 * nc2;
825     
826   // check if initialisation has been done
827   // If runloader has been initialized, set the number of events per file to nc1 * nc2
828     
829   // Set new generator
830   if (!gener) gener  = new AliLegoGenerator();
831   //
832   // Configure Generator
833
834   gener->SetRadiusRange(rmin, rmax);
835   gener->SetZMax(zmax);
836   gener->SetCoor1Range(nc1, c1min, c1max);
837   gener->SetCoor2Range(nc2, c2min, c2max);
838   
839   
840   //Create Lego object  
841   fLego = new AliLego("lego",gener);
842
843   //__________________________________________________________________________
844
845   gAlice->Announce();
846
847   gROOT->LoadMacro(setup);
848   gInterpreter->ProcessLine(gAlice->GetConfigFunction());
849
850   if(AliCDBManager::Instance()->GetRun() >= 0) { 
851     SetRunNumber(AliCDBManager::Instance()->GetRun());
852   } else {
853     AliWarning("Run number not initialized!!");
854   }
855
856   AliRunLoader::Instance()->CdGAFile();
857   
858   AliPDG::AddParticlesToPdgDataBase();  
859   
860   gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
861
862   gAlice->GetMCApp()->Init();
863   
864   
865   //Must be here because some MCs (G4) adds detectors here and not in Config.C
866   gAlice->InitLoaders();
867   AliRunLoader::Instance()->MakeTree("E");
868   
869   //
870   // Save stuff at the beginning of the file to avoid file corruption
871   AliRunLoader::Instance()->CdGAFile();
872   gAlice->Write();
873
874   //Save current generator
875   AliGenerator *gen=gAlice->GetMCApp()->Generator();
876   gAlice->GetMCApp()->ResetGenerator(gener);
877   //Prepare MC for Lego Run
878   gMC->InitLego();
879   
880   //Run Lego Object
881   
882   
883   AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
884   gMC->ProcessRun(nev);
885   
886   // End of this run, close files
887   FinishRun();
888   // Restore current generator
889   gAlice->GetMCApp()->ResetGenerator(gen);
890   // Delete Lego Object
891   delete fLego;
892
893   return kTRUE;
894 }
895
896 //_____________________________________________________________________________
897 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
898 {
899   // run the trigger
900
901   AliCodeTimerAuto("",0)
902
903   // initialize CDB storage from external environment
904   // (either CDB manager or AliSimulation setters),
905   // if not already done in RunSimulation()
906   InitCDB();
907   
908   // Set run number in CDBManager from data 
909   // From this point on the run number must be always loaded from data!
910   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
911   
912   // Set CDB lock: from now on it is forbidden to reset the run number
913   // or the default storage or to activate any further storage!
914   SetCDBLock();
915    
916    AliRunLoader* runLoader = LoadRun("READ");
917    if (!runLoader) return kFALSE;
918    TString trconfiguration = config;
919
920    if (trconfiguration.IsNull()) {
921      if(!fTriggerConfig.IsNull()) {
922        trconfiguration = fTriggerConfig;
923      }
924      else
925        AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
926    }
927
928    runLoader->MakeTree( "GG" );
929    AliCentralTrigger* aCTP = runLoader->GetTrigger();
930    // Load Configuration
931    if (!aCTP->LoadConfiguration( trconfiguration ))
932      return kFALSE;
933
934    // digits -> trigger
935    if( !aCTP->RunTrigger( runLoader , detectors ) ) {
936       if (fStopOnError) {
937         //  delete aCTP;
938         return kFALSE;
939       }
940    }
941
942    delete runLoader;
943
944    return kTRUE;
945 }
946
947 //_____________________________________________________________________________
948 Bool_t AliSimulation::WriteTriggerRawData()
949 {
950   // Writes the CTP (trigger) DDL raw data
951   // Details of the format are given in the
952   // trigger TDR - pages 134 and 135.
953   AliCTPRawData writer;
954   writer.RawData();
955
956   return kTRUE;
957 }
958
959 //_____________________________________________________________________________
960 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
961 {
962 // run the generation and simulation
963
964   AliCodeTimerAuto("",0)
965
966   // initialize CDB storage and run number from external environment
967   // (either CDB manager or AliSimulation setters)
968   AliSysInfo::AddStamp("RunSimulation_Begin");
969   InitCDB();
970   AliSysInfo::AddStamp("RunSimulation_InitCDB");
971   InitRunNumber();
972   SetCDBLock();
973   AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
974   
975   if (!gAlice) {
976     AliError("no gAlice object. Restart aliroot and try again.");
977     return kFALSE;
978   }
979   if (gAlice->Modules()->GetEntries() > 0) {
980     AliError("gAlice was already run. Restart aliroot and try again.");
981     return kFALSE;
982   }
983
984   AliInfo(Form("initializing gAlice with config file %s",
985           fConfigFileName.Data()));
986
987   //
988   // Initialize ALICE Simulation run
989   //
990   gAlice->Announce();
991
992   //
993   // If requested set the mag. field from the GRP entry.
994   // After this the field is loccked and cannot be changed by Config.C
995   if (fUseMagFieldFromGRP) {
996       AliGRPManager grpM;
997       grpM.ReadGRPEntry();
998       grpM.SetMagField();
999       AliInfo("Field is locked now. It cannot be changed in Config.C");
1000   }
1001 //
1002 // Execute Config.C
1003   gROOT->LoadMacro(fConfigFileName.Data());
1004   gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1005   AliSysInfo::AddStamp("RunSimulation_Config");
1006
1007 //
1008 // If requested obtain the vertex position and vertex sigma_z from the CDB
1009 // This overwrites the settings from the Config.C  
1010   if (fUseVertexFromCDB) {
1011       Double_t vtxPos[3] = {0., 0., 0.}; 
1012       Double_t vtxSig[3] = {0., 0., 0.};
1013       AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1014       AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1015       if(vertex) {
1016           if(vertex->GetXRes()>2.8) { // > pipe radius --> it's a dummy object, don't use it 
1017               entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1018               vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1019           }
1020       }
1021       if (vertex) {
1022           vertex->GetXYZ(vtxPos);
1023           vertex->GetSigmaXYZ(vtxSig);
1024           AliInfo("Overwriting Config.C vertex settings !");
1025           AliInfo(Form("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
1026                  vtxPos[0], vtxPos[1], vtxPos[2], vtxSig[2]));
1027
1028           AliGenerator *gen = gAlice->GetMCApp()->Generator();
1029           gen->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]);   // vertex position
1030           gen->SetSigmaZ(vtxSig[2]);
1031       }
1032   }
1033   
1034   if(AliCDBManager::Instance()->GetRun() >= 0) { 
1035     AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
1036     AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
1037   } else {
1038         AliWarning("Run number not initialized!!");
1039   }
1040   
1041    AliRunLoader::Instance()->CdGAFile();
1042    
1043
1044    AliPDG::AddParticlesToPdgDataBase();  
1045
1046    gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1047    AliSysInfo::AddStamp("RunSimulation_GetField");
1048    
1049    gAlice->GetMCApp()->Init();
1050    AliSysInfo::AddStamp("RunSimulation_InitMCApp");
1051
1052    //Must be here because some MCs (G4) adds detectors here and not in Config.C
1053    gAlice->InitLoaders();
1054    AliRunLoader::Instance()->MakeTree("E");
1055    AliRunLoader::Instance()->LoadKinematics("RECREATE");
1056    AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1057    AliRunLoader::Instance()->LoadHits("all","RECREATE");
1058    //
1059    // Save stuff at the beginning of the file to avoid file corruption
1060    AliRunLoader::Instance()->CdGAFile();
1061    gAlice->Write();
1062    gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1063    AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1064   //___________________________________________________________________________________________
1065   
1066   AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
1067
1068   // Set run number in CDBManager
1069   AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
1070
1071   AliRunLoader* runLoader = AliRunLoader::Instance();
1072   if (!runLoader) {
1073              AliError(Form("gAlice has no run loader object. "
1074                              "Check your config file: %s", fConfigFileName.Data()));
1075              return kFALSE;
1076   }
1077   SetGAliceFile(runLoader->GetFileName());
1078       
1079   // Misalign geometry
1080 #if ROOT_VERSION_CODE < 331527
1081   AliGeomManager::SetGeometry(gGeoManager);
1082   
1083   // Check that the consistency of symbolic names for the activated subdetectors
1084   // in the geometry loaded by AliGeomManager
1085   TString detsToBeChecked = "";
1086   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1087   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1088     AliModule* det = (AliModule*) detArray->At(iDet);
1089     if (!det || !det->IsActive()) continue;
1090     detsToBeChecked += det->GetName();
1091     detsToBeChecked += " ";
1092   } // end loop over detectors
1093   if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
1094     AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
1095   MisalignGeometry(runLoader);
1096   AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
1097 #endif
1098
1099 //   AliRunLoader* runLoader = AliRunLoader::Instance();
1100 //   if (!runLoader) {
1101 //     AliError(Form("gAlice has no run loader object. "
1102 //                   "Check your config file: %s", fConfigFileName.Data()));
1103 //     return kFALSE;
1104 //   }
1105 //   SetGAliceFile(runLoader->GetFileName());
1106
1107   if (!gAlice->GetMCApp()->Generator()) {
1108     AliError(Form("gAlice has no generator object. "
1109                   "Check your config file: %s", fConfigFileName.Data()));
1110     return kFALSE;
1111   }
1112
1113   // Write GRP entry corresponding to the setting found in Cofig.C
1114   if (fWriteGRPEntry)
1115     WriteGRPEntry();
1116   AliSysInfo::AddStamp("RunSimulation_WriteGRP");
1117
1118   if (nEvents <= 0) nEvents = fNEvents;
1119
1120   // get vertex from background file in case of merging
1121   if (fUseBkgrdVertex &&
1122       fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
1123     Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
1124     const char* fileName = ((TObjString*)
1125                             (fBkgrdFileNames->At(0)))->GetName();
1126     AliInfo(Form("The vertex will be taken from the background "
1127                  "file %s with nSignalPerBackground = %d", 
1128                  fileName, signalPerBkgrd));
1129     AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
1130     gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
1131   }
1132
1133   if (!fRunSimulation) {
1134     gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
1135   }
1136
1137   // set the number of events per file for given detectors and data types
1138   for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
1139     if (!fEventsPerFile[i]) continue;
1140     const char* detName = fEventsPerFile[i]->GetName();
1141     const char* typeName = fEventsPerFile[i]->GetTitle();
1142     TString loaderName(detName);
1143     loaderName += "Loader";
1144     AliLoader* loader = runLoader->GetLoader(loaderName);
1145     if (!loader) {
1146       AliError(Form("RunSimulation no loader for %s found\n Number of events per file not set for %s %s", 
1147                     detName, typeName, detName));
1148       continue;
1149     }
1150     AliDataLoader* dataLoader = 
1151       loader->GetDataLoader(typeName);
1152     if (!dataLoader) {
1153       AliError(Form("no data loader for %s found\n"
1154                     "Number of events per file not set for %s %s", 
1155                     typeName, detName, typeName));
1156       continue;
1157     }
1158     dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
1159     AliDebug(1, Form("number of events per file set to %d for %s %s",
1160                      fEventsPerFile[i]->GetUniqueID(), detName, typeName));
1161   }
1162
1163   AliInfo("running gAlice");
1164   AliSysInfo::AddStamp("Start_ProcessRun");
1165
1166   // Create the Root Tree with one branch per detector
1167   //Hits moved to begin event -> now we are crating separate tree for each event
1168
1169   gMC->ProcessRun(nEvents);
1170
1171   // End of this run, close files
1172   if(nEvents>0) FinishRun();
1173
1174   AliSysInfo::AddStamp("Stop_ProcessRun");
1175   delete runLoader;
1176
1177   return kTRUE;
1178 }
1179
1180 //_____________________________________________________________________________
1181 Bool_t AliSimulation::RunSDigitization(const char* detectors)
1182 {
1183 // run the digitization and produce summable digits
1184   static Int_t eventNr=0;
1185   AliCodeTimerAuto("",0) ;
1186
1187   // initialize CDB storage, run number, set CDB lock
1188   InitCDB();
1189   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1190   SetCDBLock();
1191   
1192   AliRunLoader* runLoader = LoadRun();
1193   if (!runLoader) return kFALSE;
1194
1195   TString detStr = detectors;
1196   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1197   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1198     AliModule* det = (AliModule*) detArray->At(iDet);
1199     if (!det || !det->IsActive()) continue;
1200     if (IsSelected(det->GetName(), detStr)) {
1201       AliInfo(Form("creating summable digits for %s", det->GetName()));
1202       AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
1203       det->Hits2SDigits();
1204       AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
1205       AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
1206     }
1207   }
1208
1209   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1210     AliError(Form("the following detectors were not found: %s",
1211                   detStr.Data()));
1212     if (fStopOnError) return kFALSE;
1213   }
1214   eventNr++;
1215   delete runLoader;
1216
1217   return kTRUE;
1218 }
1219
1220
1221 //_____________________________________________________________________________
1222 Bool_t AliSimulation::RunDigitization(const char* detectors, 
1223                                       const char* excludeDetectors)
1224 {
1225 // run the digitization and produce digits from sdigits
1226
1227   AliCodeTimerAuto("",0)
1228
1229   // initialize CDB storage, run number, set CDB lock
1230   InitCDB();
1231   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1232   SetCDBLock();
1233   
1234   delete AliRunLoader::Instance();
1235   delete gAlice;
1236   gAlice = NULL;
1237
1238   Int_t nStreams = 1;
1239   if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1240   Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1241   AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
1242   // manager->SetEmbeddingFlag(fEmbeddingFlag);
1243   manager->SetInputStream(0, fGAliceFileName.Data());
1244   for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1245     const char* fileName = ((TObjString*)
1246                             (fBkgrdFileNames->At(iStream-1)))->GetName();
1247     manager->SetInputStream(iStream, fileName);
1248   }
1249
1250   TString detStr = detectors;
1251   TString detExcl = excludeDetectors;
1252   manager->GetInputStream(0)->ImportgAlice();
1253   AliRunLoader* runLoader = 
1254     AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
1255   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1256   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1257     AliModule* det = (AliModule*) detArray->At(iDet);
1258     if (!det || !det->IsActive()) continue;
1259     if (IsSelected(det->GetName(), detStr) && 
1260         !IsSelected(det->GetName(), detExcl)) {
1261       AliDigitizer* digitizer = det->CreateDigitizer(manager);
1262       
1263       if (!digitizer) {
1264         AliError(Form("no digitizer for %s", det->GetName()));
1265         if (fStopOnError) return kFALSE;
1266       } else {
1267         digitizer->SetRegionOfInterest(fRegionOfInterest);
1268       }
1269     }
1270   }
1271
1272   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1273     AliError(Form("the following detectors were not found: %s", 
1274                   detStr.Data()));
1275     if (fStopOnError) return kFALSE;
1276   }
1277
1278   if (!manager->GetListOfTasks()->IsEmpty()) {
1279     AliInfo("executing digitization");
1280     manager->Exec("");
1281   }
1282
1283   delete manager;
1284
1285   return kTRUE;
1286 }
1287
1288 //_____________________________________________________________________________
1289 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1290 {
1291 // run the digitization and produce digits from hits
1292
1293   AliCodeTimerAuto("",0)
1294
1295   // initialize CDB storage, run number, set CDB lock
1296   InitCDB();
1297   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1298   SetCDBLock();
1299   
1300   AliRunLoader* runLoader = LoadRun("READ");
1301   if (!runLoader) return kFALSE;
1302
1303   TString detStr = detectors;
1304   TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1305   for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1306     AliModule* det = (AliModule*) detArray->At(iDet);
1307     if (!det || !det->IsActive()) continue;
1308     if (IsSelected(det->GetName(), detStr)) {
1309       AliInfo(Form("creating digits from hits for %s", det->GetName()));
1310       det->Hits2Digits();
1311     }
1312   }
1313
1314   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1315     AliError(Form("the following detectors were not found: %s", 
1316                   detStr.Data()));
1317     if (fStopOnError) return kFALSE;
1318   }
1319
1320   return kTRUE;
1321 }
1322
1323 //_____________________________________________________________________________
1324 Bool_t AliSimulation::WriteRawData(const char* detectors, 
1325                                    const char* fileName,
1326                                    Bool_t deleteIntermediateFiles,
1327                                    Bool_t selrawdata)
1328 {
1329 // convert the digits to raw data
1330 // First DDL raw data files for the given detectors are created.
1331 // If a file name is given, the DDL files are then converted to a DATE file.
1332 // If deleteIntermediateFiles is true, the DDL raw files are deleted 
1333 // afterwards.
1334 // If the file name has the extension ".root", the DATE file is converted
1335 // to a root file.
1336 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1337 // 'selrawdata' flag can be used to enable writing of detectors raw data
1338 // accoring to the trigger cluster.
1339
1340   AliCodeTimerAuto("",0)
1341   AliSysInfo::AddStamp("WriteRawData_Start");
1342   
1343   TString detStr = detectors;
1344   if (!WriteRawFiles(detStr.Data())) {
1345     if (fStopOnError) return kFALSE;
1346   }
1347   AliSysInfo::AddStamp("WriteRawFiles");
1348
1349   // run HLT simulation on simulated DDL raw files
1350   // and produce HLT ddl raw files to be included in date/root file
1351   // bugfix 2009-06-26: the decision whether to write HLT raw data
1352   // is taken in RunHLT. Here HLT always needs to be run in order to
1353   // create HLT digits, unless its switched off. This is due to the
1354   // special placement of the HLT between the generation of DDL files
1355   // and conversion to DATE/Root file.
1356   detStr.ReplaceAll("HLT", "");
1357   if (!fRunHLT.IsNull()) {
1358     if (!RunHLT()) {
1359       if (fStopOnError) return kFALSE;
1360     }
1361   }
1362   AliSysInfo::AddStamp("WriteRawData_RunHLT");
1363
1364   TString dateFileName(fileName);
1365   if (!dateFileName.IsNull()) {
1366     Bool_t rootOutput = dateFileName.EndsWith(".root");
1367     if (rootOutput) dateFileName += ".date";
1368     TString selDateFileName;
1369     if (selrawdata) {
1370       selDateFileName = "selected.";
1371       selDateFileName+= dateFileName;
1372     }
1373     if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1374       if (fStopOnError) return kFALSE;
1375     }
1376     AliSysInfo::AddStamp("ConvertRawFilesToDate");
1377     if (deleteIntermediateFiles) {
1378       AliRunLoader* runLoader = LoadRun("READ");
1379       if (runLoader) for (Int_t iEvent = 0; 
1380                           iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1381         char command[256];
1382         sprintf(command, "rm -r raw%d", iEvent);
1383         gSystem->Exec(command);
1384       }
1385       delete runLoader;
1386     }
1387
1388     if (rootOutput) {
1389       if (!ConvertDateToRoot(dateFileName, fileName)) {
1390         if (fStopOnError) return kFALSE;
1391       }
1392       AliSysInfo::AddStamp("ConvertDateToRoot");
1393       if (deleteIntermediateFiles) {
1394         gSystem->Unlink(dateFileName);
1395       }
1396       if (selrawdata) {
1397         TString selFileName = "selected.";
1398         selFileName        += fileName;
1399         if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1400           if (fStopOnError) return kFALSE;
1401         }
1402         if (deleteIntermediateFiles) {
1403           gSystem->Unlink(selDateFileName);
1404         }
1405       }
1406     }
1407   }
1408
1409   return kTRUE;
1410 }
1411
1412 //_____________________________________________________________________________
1413 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1414 {
1415 // convert the digits to raw data DDL files
1416
1417   AliCodeTimerAuto("",0)
1418   
1419   AliRunLoader* runLoader = LoadRun("READ");
1420   if (!runLoader) return kFALSE;
1421
1422   // write raw data to DDL files
1423   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1424     AliInfo(Form("processing event %d", iEvent));
1425     runLoader->GetEvent(iEvent);
1426     TString baseDir = gSystem->WorkingDirectory();
1427     char dirName[256];
1428     sprintf(dirName, "raw%d", iEvent);
1429     gSystem->MakeDirectory(dirName);
1430     if (!gSystem->ChangeDirectory(dirName)) {
1431       AliError(Form("couldn't change to directory %s", dirName));
1432       if (fStopOnError) return kFALSE; else continue;
1433     }
1434
1435     ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1436     runNbFile.close();
1437
1438     TString detStr = detectors;
1439     if (IsSelected("HLT", detStr)) {
1440       // Do nothing. "HLT" will be removed from detStr and HLT raw
1441       // data files are generated in RunHLT.
1442     }
1443
1444     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1445     for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1446       AliModule* det = (AliModule*) detArray->At(iDet);
1447       if (!det || !det->IsActive()) continue;
1448       if (IsSelected(det->GetName(), detStr)) {
1449         AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1450         det->Digits2Raw();
1451       }
1452     }
1453
1454     if (!WriteTriggerRawData())
1455       if (fStopOnError) return kFALSE;
1456
1457     gSystem->ChangeDirectory(baseDir);
1458     if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1459       AliError(Form("the following detectors were not found: %s", 
1460                     detStr.Data()));
1461       if (fStopOnError) return kFALSE;
1462     }
1463   }
1464
1465   delete runLoader;
1466   
1467   return kTRUE;
1468 }
1469
1470 //_____________________________________________________________________________
1471 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1472                                             const char* selDateFileName)
1473 {
1474 // convert raw data DDL files to a DATE file with the program "dateStream"
1475 // The second argument is not empty when the user decides to write
1476 // the detectors raw data according to the trigger cluster.
1477
1478   AliCodeTimerAuto("",0)
1479   
1480   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1481   if (!path) {
1482     AliError("the program dateStream was not found");
1483     if (fStopOnError) return kFALSE;
1484   } else {
1485     delete[] path;
1486   }
1487
1488   AliRunLoader* runLoader = LoadRun("READ");
1489   if (!runLoader) return kFALSE;
1490
1491   AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1492   Bool_t selrawdata = kFALSE;
1493   if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1494
1495   char command[256];
1496   // Note the option -s. It is used in order to avoid
1497   // the generation of SOR/EOR events.
1498   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1499           dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1500   FILE* pipe = gSystem->OpenPipe(command, "w");
1501
1502   if (!pipe) {
1503     AliError(Form("Cannot execute command: %s",command));
1504     return kFALSE;
1505   }
1506
1507   Int_t selEvents = 0;
1508   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1509
1510     UInt_t detectorPattern = 0;
1511     runLoader->GetEvent(iEvent);
1512     if (!runLoader->LoadTrigger()) {
1513       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1514       detectorPattern = aCTP->GetClusterMask();
1515       // Check if the event was triggered by CTP
1516       if (selrawdata) {
1517         if (aCTP->GetClassMask()) selEvents++;
1518       }
1519     }
1520     else {
1521       AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1522       if (selrawdata) {
1523         AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1524         selrawdata = kFALSE;
1525       }
1526     }
1527
1528     fprintf(pipe, "GDC DetectorPattern %u\n", detectorPattern);
1529     Float_t ldc = 0;
1530     Int_t prevLDC = -1;
1531
1532     // loop over detectors and DDLs
1533     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1534       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1535
1536         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1537         Int_t ldcID = Int_t(ldc + 0.0001);
1538         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1539
1540         char rawFileName[256];
1541         sprintf(rawFileName, "raw%d/%s", 
1542                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1543
1544         // check existence and size of raw data file
1545         FILE* file = fopen(rawFileName, "rb");
1546         if (!file) continue;
1547         fseek(file, 0, SEEK_END);
1548         unsigned long size = ftell(file);
1549         fclose(file);
1550         if (!size) continue;
1551
1552         if (ldcID != prevLDC) {
1553           fprintf(pipe, " LDC Id %d\n", ldcID);
1554           prevLDC = ldcID;
1555         }
1556         fprintf(pipe, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1557       }
1558     }
1559   }
1560
1561   Int_t result = gSystem->ClosePipe(pipe);
1562
1563   if (!(selrawdata && selEvents > 0)) {
1564     delete runLoader;
1565     return (result == 0);
1566   }
1567
1568   AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1569   
1570   sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d", 
1571           selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1572   FILE* pipe2 = gSystem->OpenPipe(command, "w");
1573
1574   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1575
1576     // Get the trigger decision and cluster
1577     UInt_t detectorPattern = 0;
1578     TString detClust;
1579     runLoader->GetEvent(iEvent);
1580     if (!runLoader->LoadTrigger()) {
1581       AliCentralTrigger *aCTP = runLoader->GetTrigger();
1582       if (aCTP->GetClassMask() == 0) continue;
1583       detectorPattern = aCTP->GetClusterMask();
1584       detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1585       AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1586     }
1587
1588     fprintf(pipe2, "GDC DetectorPattern %u\n", detectorPattern);
1589     Float_t ldc = 0;
1590     Int_t prevLDC = -1;
1591
1592     // loop over detectors and DDLs
1593     for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1594       // Write only raw data from detectors that
1595       // are contained in the trigger cluster(s)
1596       if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1597
1598       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1599
1600         Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1601         Int_t ldcID = Int_t(ldc + 0.0001);
1602         ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1603
1604         char rawFileName[256];
1605         sprintf(rawFileName, "raw%d/%s", 
1606                 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1607
1608         // check existence and size of raw data file
1609         FILE* file = fopen(rawFileName, "rb");
1610         if (!file) continue;
1611         fseek(file, 0, SEEK_END);
1612         unsigned long size = ftell(file);
1613         fclose(file);
1614         if (!size) continue;
1615
1616         if (ldcID != prevLDC) {
1617           fprintf(pipe2, " LDC Id %d\n", ldcID);
1618           prevLDC = ldcID;
1619         }
1620         fprintf(pipe2, "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
1621       }
1622     }
1623   }
1624
1625   Int_t result2 = gSystem->ClosePipe(pipe2);
1626
1627   delete runLoader;
1628   return ((result == 0) && (result2 == 0));
1629 }
1630
1631 //_____________________________________________________________________________
1632 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1633                                         const char* rootFileName)
1634 {
1635 // convert a DATE file to a root file with the program "alimdc"
1636
1637   // ALIMDC setup
1638   const Int_t kDBSize = 2000000000;
1639   const Int_t kTagDBSize = 1000000000;
1640   const Bool_t kFilter = kFALSE;
1641   const Int_t kCompression = 1;
1642
1643   char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1644   if (!path) {
1645     AliError("the program alimdc was not found");
1646     if (fStopOnError) return kFALSE;
1647   } else {
1648     delete[] path;
1649   }
1650
1651   AliInfo(Form("converting DATE file %s to root file %s", 
1652                dateFileName, rootFileName));
1653
1654   const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1655   const char* tagDBFS    = "/tmp/mdc1/tags";
1656
1657   // User defined file system locations
1658   if (gSystem->Getenv("ALIMDC_RAWDB1")) 
1659     rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1660   if (gSystem->Getenv("ALIMDC_RAWDB2")) 
1661     rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1662   if (gSystem->Getenv("ALIMDC_TAGDB")) 
1663     tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1664
1665   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1666   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1667   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1668
1669   gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1670   gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1671   gSystem->Exec(Form("mkdir %s",tagDBFS));
1672
1673   Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s", 
1674                                     kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1675   gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1676
1677   gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1678   gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1679   gSystem->Exec(Form("rm -rf %s",tagDBFS));
1680
1681   return (result == 0);
1682 }
1683
1684
1685 //_____________________________________________________________________________
1686 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1687 {
1688 // delete existing run loaders, open a new one and load gAlice
1689
1690   delete AliRunLoader::Instance();
1691   AliRunLoader* runLoader = 
1692     AliRunLoader::Open(fGAliceFileName.Data(), 
1693                        AliConfig::GetDefaultEventFolderName(), mode);
1694   if (!runLoader) {
1695     AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1696     return NULL;
1697   }
1698   runLoader->LoadgAlice();
1699   runLoader->LoadHeader();
1700   gAlice = runLoader->GetAliRun();
1701   if (!gAlice) {
1702     AliError(Form("no gAlice object found in file %s", 
1703                   fGAliceFileName.Data()));
1704     return NULL;
1705   }
1706   return runLoader;
1707 }
1708
1709 //_____________________________________________________________________________
1710 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1711 {
1712 // get or calculate the number of signal events per background event
1713
1714   if (!fBkgrdFileNames) return 1;
1715   Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1716   if (nBkgrdFiles == 0) return 1;
1717
1718   // get the number of signal events
1719   if (nEvents <= 0) {
1720     AliRunLoader* runLoader = 
1721         AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1722     if (!runLoader) return 1;
1723     
1724     nEvents = runLoader->GetNumberOfEvents();
1725     delete runLoader;
1726   }
1727
1728   Int_t result = 0;
1729   for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1730     // get the number of background events
1731     const char* fileName = ((TObjString*)
1732                             (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1733     AliRunLoader* runLoader =
1734       AliRunLoader::Open(fileName, "BKGRD");
1735     if (!runLoader) continue;
1736     Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1737     delete runLoader;
1738   
1739     // get or calculate the number of signal per background events
1740     Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1741     if (nSignalPerBkgrd <= 0) {
1742       nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1743     } else if (result && (result != nSignalPerBkgrd)) {
1744       AliInfo(Form("the number of signal events per background event "
1745                    "will be changed from %d to %d for stream %d", 
1746                    nSignalPerBkgrd, result, iBkgrdFile+1));
1747       nSignalPerBkgrd = result;
1748     }
1749
1750     if (!result) result = nSignalPerBkgrd;
1751     if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1752       AliWarning(Form("not enough background events (%d) for %d signal events "
1753                       "using %d signal per background events for stream %d",
1754                       nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1755     }
1756   }
1757
1758   return result;
1759 }
1760
1761 //_____________________________________________________________________________
1762 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1763 {
1764 // check whether detName is contained in detectors
1765 // if yes, it is removed from detectors
1766
1767   // check if all detectors are selected
1768   if ((detectors.CompareTo("ALL") == 0) ||
1769       detectors.BeginsWith("ALL ") ||
1770       detectors.EndsWith(" ALL") ||
1771       detectors.Contains(" ALL ")) {
1772     detectors = "ALL";
1773     return kTRUE;
1774   }
1775
1776   // search for the given detector
1777   Bool_t result = kFALSE;
1778   if ((detectors.CompareTo(detName) == 0) ||
1779       detectors.BeginsWith(detName+" ") ||
1780       detectors.EndsWith(" "+detName) ||
1781       detectors.Contains(" "+detName+" ")) {
1782     detectors.ReplaceAll(detName, "");
1783     result = kTRUE;
1784   }
1785
1786   // clean up the detectors string
1787   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
1788   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1789   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1790
1791   return result;
1792 }
1793
1794 //_____________________________________________________________________________
1795 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName, Int_t N) 
1796 {
1797 //
1798 // Steering routine  to convert raw data in directory rawDirectory/ to fake SDigits. 
1799 // These can be used for embedding of MC tracks into RAW data using the standard 
1800 // merging procedure.
1801 //
1802 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1803 //
1804   if (!gAlice) {
1805     AliError("no gAlice object. Restart aliroot and try again.");
1806     return kFALSE;
1807   }
1808   if (gAlice->Modules()->GetEntries() > 0) {
1809     AliError("gAlice was already run. Restart aliroot and try again.");
1810     return kFALSE;
1811   }
1812   
1813   AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1814   
1815   gAlice->Announce();
1816   
1817   gROOT->LoadMacro(fConfigFileName.Data());
1818   gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1819   
1820   if(AliCDBManager::Instance()->GetRun() >= 0) { 
1821         SetRunNumber(AliCDBManager::Instance()->GetRun());
1822   } else {
1823         AliWarning("Run number not initialized!!");
1824   }
1825   
1826    AliRunLoader::Instance()->CdGAFile();
1827     
1828    AliPDG::AddParticlesToPdgDataBase();  
1829
1830    gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1831    
1832    gAlice->GetMCApp()->Init();
1833
1834    //Must be here because some MCs (G4) adds detectors here and not in Config.C
1835    gAlice->InitLoaders();
1836    AliRunLoader::Instance()->MakeTree("E");
1837    AliRunLoader::Instance()->LoadKinematics("RECREATE");
1838    AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1839    AliRunLoader::Instance()->LoadHits("all","RECREATE");
1840
1841    //
1842    // Save stuff at the beginning of the file to avoid file corruption
1843    AliRunLoader::Instance()->CdGAFile();
1844    gAlice->Write();
1845 //
1846 //  Initialize CDB     
1847     InitCDB();
1848     //AliCDBManager* man = AliCDBManager::Instance();
1849     //man->SetRun(0); // Should this come from rawdata header ?
1850     
1851     Int_t iDet;
1852     //
1853     // Get the runloader
1854     AliRunLoader* runLoader = AliRunLoader::Instance();
1855     //
1856     // Open esd file if available
1857     TFile* esdFile = 0;
1858     TTree* treeESD = 0;
1859     AliESDEvent* esd = 0;
1860     if (esdFileName && (strlen(esdFileName)>0)) {
1861       esdFile = TFile::Open(esdFileName);
1862       if (esdFile) {
1863         esd = new AliESDEvent();
1864         esdFile->GetObject("esdTree", treeESD);
1865         if (treeESD) esd->ReadFromTree(treeESD);
1866       }
1867     }
1868
1869     //
1870     // Create the RawReader
1871     TString fileName(rawDirectory);
1872     AliRawReader* rawReader = 0x0;
1873     if (fileName.EndsWith("/")) {
1874       rawReader = new AliRawReaderFile(fileName);
1875     } else if (fileName.EndsWith(".root")) {
1876       rawReader = new AliRawReaderRoot(fileName);
1877     } else if (!fileName.IsNull()) {
1878       rawReader = new AliRawReaderDate(fileName);
1879     }
1880 //     if (!fEquipIdMap.IsNull() && fRawReader)
1881 //       fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1882     //
1883     // Get list of detectors
1884     TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1885     //
1886     // Get Header
1887     AliHeader* header = runLoader->GetHeader();
1888     // Event loop
1889     Int_t nev = 0;
1890     while(kTRUE) {
1891         if (!(rawReader->NextEvent())) break;
1892         runLoader->SetEventNumber(nev);
1893         runLoader->GetHeader()->Reset(rawReader->GetRunNumber(), 
1894                                       nev, nev);
1895         runLoader->TreeE()->Fill();
1896         runLoader->GetEvent(nev);
1897         AliInfo(Form("We are at event %d",nev));
1898         //
1899         // Detector loop
1900         TString detStr = fMakeSDigits;
1901         for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1902             AliModule* det = (AliModule*) detArray->At(iDet);
1903             if (!det || !det->IsActive()) continue;
1904             if (IsSelected(det->GetName(), detStr)) {
1905               AliInfo(Form("Calling Raw2SDigits for %s", det->GetName()));
1906               det->Raw2SDigits(rawReader);
1907               rawReader->Reset();
1908             }
1909         } // detectors
1910         
1911
1912         //
1913         //  If ESD information available obtain reconstructed vertex and store in header.
1914         if (treeESD) {
1915             treeESD->GetEvent(nev);
1916             const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1917             Double_t position[3];
1918             esdVertex->GetXYZ(position);
1919             AliGenEventHeader* mcHeader = new  AliGenEventHeader("ESD");
1920             TArrayF mcV;
1921             mcV.Set(3);
1922             for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1923             mcHeader->SetPrimaryVertex(mcV);
1924             header->Reset(0,nev);
1925             header->SetGenEventHeader(mcHeader);
1926             AliInfo(Form("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]));
1927         }
1928 //
1929 //      Finish the event
1930         AliInfo(Form("Finished event %d",nev));
1931         nev++;
1932         if (N>0&&nev>=N)
1933           break;
1934     } // events
1935  
1936     delete rawReader;
1937 //
1938 //  Finish the run 
1939     runLoader->CdGAFile();
1940     runLoader->WriteHeader("OVERWRITE");
1941     runLoader->WriteRunLoader();
1942
1943     return kTRUE;
1944 }
1945
1946 //_____________________________________________________________________________
1947 void AliSimulation::FinishRun()
1948 {
1949   //
1950   // Called at the end of the run.
1951   //
1952
1953   if(IsLegoRun()) 
1954    {
1955     AliDebug(1, "Finish Lego");
1956     AliRunLoader::Instance()->CdGAFile();
1957     fLego->FinishRun();
1958    }
1959   
1960   // Clean detector information
1961   TIter next(gAlice->Modules());
1962   AliModule *detector;
1963   while((detector = dynamic_cast<AliModule*>(next()))) {
1964     AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
1965     detector->FinishRun();
1966   }
1967   
1968   AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
1969   AliRunLoader::Instance()->WriteHeader("OVERWRITE");
1970
1971   // Write AliRun info and all detectors parameters
1972   AliRunLoader::Instance()->CdGAFile();
1973   gAlice->Write(0,TObject::kOverwrite);//write AliRun
1974   AliRunLoader::Instance()->Write(0,TObject::kOverwrite);//write RunLoader itself
1975   
1976   if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();  
1977   AliRunLoader::Instance()->Synchronize();
1978 }
1979
1980 //_____________________________________________________________________________
1981 Int_t AliSimulation::GetDetIndex(const char* detector)
1982 {
1983   // return the detector index corresponding to detector
1984   Int_t index = -1 ; 
1985   for (index = 0; index < fgkNDetectors ; index++) {
1986     if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1987           break ; 
1988   }     
1989   return index ; 
1990 }
1991
1992 //_____________________________________________________________________________
1993 Bool_t AliSimulation::CreateHLT()
1994 {
1995   // Init the HLT simulation.
1996   // The function  loads the library and creates the instance of AliHLTSimulation.
1997   // the main reason for the decoupled creation is to set the transient OCDB
1998   // objects before the OCDB is locked
1999
2000   // load the library dynamically
2001   gSystem->Load(ALIHLTSIMULATION_LIBRARY);
2002
2003   // check for the library version
2004   AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
2005   if (!fctVersion) {
2006     AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
2007     return kFALSE;
2008   }
2009   if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
2010     AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
2011   }
2012
2013   // print compile info
2014   typedef void (*CompileInfo)( const char*& date, const char*& time);
2015   CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
2016   if (fctInfo) {
2017     const char* date="";
2018     const char* time="";
2019     (*fctInfo)(date, time);
2020     if (!date) date="unknown";
2021     if (!time) time="unknown";
2022     AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
2023   } else {
2024     AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
2025   }
2026
2027   // create instance of the HLT simulation
2028   AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
2029   if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
2030     AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
2031     return kFALSE;    
2032   }
2033
2034   TString specObjects;
2035   for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
2036     if (specObjects.Length()>0) specObjects+=" ";
2037     specObjects+=fSpecCDBUri[i]->GetName();
2038   }
2039
2040   AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
2041   if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
2042     AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
2043   }
2044
2045   return kTRUE;
2046 }
2047
2048 //_____________________________________________________________________________
2049 Bool_t AliSimulation::RunHLT()
2050 {
2051   // Run the HLT simulation
2052   // HLT simulation is implemented in HLT/sim/AliHLTSimulation
2053   // Disabled if fRunHLT is empty, default vaule is "default".
2054   // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
2055   // The default simulation depends on the HLT component libraries and their
2056   // corresponding agents which define components and chains to run. See
2057   // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
2058   // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
2059   //
2060   // The libraries to be loaded can be specified as an option.
2061   // <pre>
2062   // AliSimulation sim;
2063   // sim.SetRunHLT("libAliHLTSample.so");
2064   // </pre>
2065   // will only load <tt>libAliHLTSample.so</tt>
2066
2067   // Other available options:
2068   // \li loglevel=<i>level</i> <br>
2069   //     logging level for this processing
2070   // \li alilog=off
2071   //     disable redirection of log messages to AliLog class
2072   // \li config=<i>macro</i>
2073   //     configuration macro
2074   // \li chains=<i>configuration</i>
2075   //     comma separated list of configurations to be run during simulation
2076   // \li rawfile=<i>file</i>
2077   //     source for the RawReader to be created, the default is <i>./</i> if
2078   //     raw data is simulated
2079
2080   int iResult=0;
2081
2082   if (!fpHLT && !CreateHLT()) {
2083     return kFALSE;
2084   }
2085   AliHLTSimulation* pHLT=fpHLT;
2086
2087   AliRunLoader* pRunLoader = LoadRun("READ");
2088   if (!pRunLoader) return kFALSE;
2089
2090   // initialize CDB storage, run number, set CDB lock
2091   // thats for the case of running HLT simulation without all the other steps
2092   // multiple calls are handled by the function, so we can just call
2093   InitCDB();
2094   if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
2095   SetCDBLock();
2096   
2097   // init the HLT simulation
2098   TString options;
2099   if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2100   TString detStr = fWriteRawData;
2101   if (!IsSelected("HLT", detStr)) {
2102     options+=" writerawfiles=";
2103   } else {
2104     options+=" writerawfiles=HLT";
2105   }
2106
2107   if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2108     // as a matter of fact, HLT will run reconstruction and needs the RawReader
2109     // in order to get detector data. By default, RawReaderFile is used to read
2110     // the already simulated ddl files. Date and Root files from the raw data
2111     // are generated after the HLT simulation.
2112     options+=" rawfile=./";
2113   }
2114
2115   AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2116   if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2117     AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2118   } else {
2119     // run the HLT simulation
2120     AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2121     if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2122       AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2123     }
2124   }
2125
2126   // delete the instance
2127   AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2128   if (fctDelete==NULL || fctDelete(pHLT)<0) {
2129     AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2130   }
2131   pHLT=NULL;
2132
2133   return iResult>=0?kTRUE:kFALSE;
2134 }
2135
2136 //_____________________________________________________________________________
2137 Bool_t AliSimulation::RunQA()
2138 {
2139         // run the QA on summable hits, digits or digits
2140         
2141     //if(!gAlice) return kFALSE;
2142         AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2143
2144         TString detectorsw("") ;  
2145         Bool_t rv = kTRUE ; 
2146   AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2147         detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ; 
2148         if ( detectorsw.IsNull() ) 
2149                 rv = kFALSE ; 
2150         return rv ; 
2151 }
2152
2153 //_____________________________________________________________________________
2154 Bool_t AliSimulation::SetRunQA(TString detAndAction) 
2155 {
2156         // Allows to run QA for a selected set of detectors
2157         // and a selected set of tasks among HITS, SDIGITS and DIGITS
2158         // all selected detectors run the same selected tasks
2159         
2160         if (!detAndAction.Contains(":")) {
2161                 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2162                 fRunQA = kFALSE ;
2163                 return kFALSE ;                 
2164         }
2165         Int_t colon = detAndAction.Index(":") ; 
2166         fQADetectors = detAndAction(0, colon) ; 
2167         if (fQADetectors.Contains("ALL") ){
2168     TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ; 
2169     Int_t minus = fQADetectors.Last('-') ; 
2170     TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ; 
2171     TString toRemove("") ;
2172     while (minus >= 0) {
2173       toRemove = fQADetectors(minus+1, fQADetectors.Length()) ; 
2174       toRemove = toRemove.Strip() ; 
2175       toKeep.ReplaceAll(toRemove, "") ; 
2176       fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ; 
2177       minus = fQADetectors.Last('-') ; 
2178     }
2179     fQADetectors = toKeep ; 
2180   }
2181   fQATasks   = detAndAction(colon+1, detAndAction.Sizeof() ) ; 
2182         if (fQATasks.Contains("ALL") ) {
2183                 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ; 
2184         } else {
2185                 fQATasks.ToUpper() ; 
2186                 TString tempo("") ; 
2187                 if ( fQATasks.Contains("HIT") ) 
2188                         tempo = Form("%d ", AliQAv1::kHITS) ; 
2189                 if ( fQATasks.Contains("SDIGIT") ) 
2190                         tempo += Form("%d ", AliQAv1::kSDIGITS) ; 
2191                 if ( fQATasks.Contains("DIGIT") ) 
2192                         tempo += Form("%d ", AliQAv1::kDIGITS) ; 
2193                 fQATasks = tempo ; 
2194                 if (fQATasks.IsNull()) {
2195                         AliInfo("No QA requested\n")  ;
2196                         fRunQA = kFALSE ;
2197                         return kTRUE ; 
2198                 }
2199         }       
2200         TString tempo(fQATasks) ; 
2201     tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS))  ;
2202     tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;  
2203     tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;    
2204         AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;  
2205         fRunQA = kTRUE ;
2206         AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ; 
2207         AliQAManager::QAManager()->SetTasks(fQATasks) ; 
2208   for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) 
2209     AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2210   
2211         return kTRUE; 
2212
2213
2214 //_____________________________________________________________________________
2215 void AliSimulation::ProcessEnvironmentVars()
2216 {
2217 // Extract run number and random generator seed from env variables
2218
2219     AliInfo("Processing environment variables");
2220     
2221     // Random Number seed
2222     
2223     // first check that seed is not already set
2224     if (fSeed == 0) {
2225         if (gSystem->Getenv("CONFIG_SEED")) {
2226                 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2227         }
2228     } else {
2229         if (gSystem->Getenv("CONFIG_SEED")) {
2230                 AliInfo(Form("Seed for random number generation already set (%d)"
2231                              ": CONFIG_SEED variable ignored!", fSeed));
2232         }
2233     }
2234    
2235     AliInfo(Form("Seed for random number generation = %d ", fSeed)); 
2236
2237     // Run Number
2238     
2239     // first check that run number is not already set
2240     if(fRun < 0) {    
2241         if (gSystem->Getenv("DC_RUN")) {
2242                 fRun = atoi(gSystem->Getenv("DC_RUN"));
2243         }
2244     } else {
2245         if (gSystem->Getenv("DC_RUN")) {
2246                 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2247         }
2248     }
2249     
2250     AliInfo(Form("Run number = %d", fRun)); 
2251 }
2252
2253 //---------------------------------------------------------------------
2254 void AliSimulation::WriteGRPEntry()
2255 {
2256   // Get the necessary information from galice (generator, trigger etc) and
2257   // write a GRP entry corresponding to the settings in the Config.C used
2258   // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2259
2260
2261   AliInfo("Writing global run parameters entry into the OCDB");
2262
2263   AliGRPObject* grpObj = new AliGRPObject();
2264
2265   grpObj->SetRunType("PHYSICS");
2266   grpObj->SetTimeStart(0);
2267   TDatime curtime;
2268   grpObj->SetTimeStart(0);
2269   grpObj->SetTimeEnd(curtime.Convert()); 
2270   grpObj->SetBeamEnergyIsSqrtSHalfGeV(); // new format of GRP: store sqrt(s)/2 in GeV
2271
2272   const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2273   Int_t a = 0;
2274   Int_t z = 0;
2275
2276   if (gen) {
2277     TString projectile;
2278     gen->GetProjectile(projectile,a,z);
2279     TString target;
2280     gen->GetTarget(target,a,z);
2281     TString beamType = projectile + "-" + target;
2282     beamType.ReplaceAll(" ","");
2283     if (!beamType.CompareTo("-")) {
2284       grpObj->SetBeamType("UNKNOWN");
2285       grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
2286     }
2287     else {
2288       grpObj->SetBeamType(beamType);
2289       if (z != 0) {
2290           grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 * a / z);
2291       } else {
2292           grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 );
2293       }
2294       // Heavy ion run, the event specie is set to kHighMult
2295       fEventSpecie = AliRecoParam::kHighMult;
2296       if ((strcmp(beamType,"p-p") == 0) ||
2297           (strcmp(beamType,"p-")  == 0) ||
2298           (strcmp(beamType,"-p")  == 0) ||
2299           (strcmp(beamType,"P-P") == 0) ||
2300           (strcmp(beamType,"P-")  == 0) ||
2301           (strcmp(beamType,"-P")  == 0)) {
2302         // Proton run, the event specie is set to kLowMult
2303         fEventSpecie = AliRecoParam::kLowMult;
2304       } 
2305     }
2306   } else {
2307     AliWarning("Unknown beam type and energy! Setting energy to 0");
2308     grpObj->SetBeamEnergy(0);
2309     grpObj->SetBeamType("UNKNOWN");
2310   }
2311
2312   UInt_t detectorPattern  = 0;
2313   Int_t nDets = 0;
2314   TObjArray *detArray = gAlice->Detectors();
2315   for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2316     if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2317       detectorPattern |= (1 << iDet);
2318       nDets++;
2319     }
2320   }
2321   // CTP
2322   if (!fTriggerConfig.IsNull())
2323     detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2324
2325   // HLT
2326   if (!fRunHLT.IsNull())
2327     detectorPattern |= (1 << AliDAQ::kHLTId);
2328
2329   grpObj->SetNumberOfDetectors((Char_t)nDets);
2330   grpObj->SetDetectorMask((Int_t)detectorPattern);
2331   grpObj->SetLHCPeriod("LHC08c");
2332   grpObj->SetLHCState("STABLE_BEAMS");
2333   //
2334   AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2335   Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2336
2337   Float_t factorSol     = field ? field->GetFactorSol() : 0;
2338   Float_t currentSol    = TMath::Abs(factorSol)>1E-6 ? 
2339     TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
2340   //
2341   Float_t factorDip     = field ? field->GetFactorDip() : 0;
2342   Float_t currentDip    = 6000.*TMath::Abs(factorDip);
2343   //
2344   grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
2345   grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);  
2346   grpObj->SetL3Polarity(factorSol>0 ? 0:1);  
2347   grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
2348   grpObj->SetUniformBMap(field->IsUniform());            // for special MC with k5kGUniform map
2349   grpObj->SetPolarityConventionLHC();                    // LHC convention +/+ current -> -/- field main components
2350   //
2351   grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2352   
2353   //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2354
2355   // Now store the entry in OCDB
2356   AliCDBManager* man = AliCDBManager::Instance();
2357   
2358   man->SetLock(0, fKey);
2359   
2360   AliCDBStorage* sto = man->GetStorage(fGRPWriteLocation.Data());
2361   
2362
2363   AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2364   AliCDBMetaData *metadata= new AliCDBMetaData();
2365
2366   metadata->SetResponsible("alice-off@cern.ch");
2367   metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2368  
2369   sto->Put(grpObj,id,metadata);
2370   man->SetLock(1, fKey);
2371 }
2372
2373