1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running generation, simulation and digitization //
22 // Hits, sdigits and digits are created for all detectors by typing: //
24 // AliSimulation sim; //
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 //
31 // sim.SetNumberOfEvents(n); //
33 // The name of the configuration file can be passed as argument to the //
34 // AliSimulation constructor or can be specified by //
36 // sim.SetConfigFile("..."); //
38 // The generation of particles and the simulation of detector hits can be //
39 // switched on or off by //
41 // sim.SetRunGeneration(kTRUE); // generation of primary particles //
42 // sim.SetRunSimulation(kFALSE); // but no tracking //
44 // For which detectors sdigits and digits will be created, can be steered //
47 // sim.SetMakeSDigits("ALL"); // make sdigits for all detectors //
48 // sim.SetMakeDigits("ITS TPC"); // make digits only for ITS and TPC //
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. //
55 // The creation of digits from hits instead of from sdigits can be selected //
58 // sim.SetMakeDigitsFromHits("TRD"); //
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. //
64 // Background events can be merged by calling //
66 // sim.MergeWith("background/galice.root", 2); //
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. //
75 // The output of raw data can be switched on by calling //
77 // sim.SetWriteRawData("MUON"); // write raw data for MUON //
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 //
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. //
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 //
98 // sim.SetEventsPerFile("PHOS", "Reconstructed Points", 3); //
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. //
106 ///////////////////////////////////////////////////////////////////////////////
110 #include <TGeoGlobalMagField.h>
111 #include <TGeoManager.h>
112 #include <TObjString.h>
115 #include <TVirtualMC.h>
116 #include <TVirtualMCApplication.h>
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"
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"
138 #include "AliLegoGenerator.h"
142 #include "AliModule.h"
144 #include "AliRawReaderDate.h"
145 #include "AliRawReaderFile.h"
146 #include "AliRawReaderRoot.h"
148 #include "AliRunDigitizer.h"
149 #include "AliRunLoader.h"
150 #include "AliSimulation.h"
151 #include "AliSysInfo.h"
152 #include "AliVertexGenFile.h"
154 ClassImp(AliSimulation)
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"};
159 //_____________________________________________________________________________
160 AliSimulation::AliSimulation(const char* configFileName,
161 const char* name, const char* title) :
164 fRunGeneration(kTRUE),
165 fRunSimulation(kTRUE),
166 fLoadAlignFromCDB(kTRUE),
167 fLoadAlObjsListOfDets("ALL"),
171 fMakeDigitsFromHits(""),
173 fRawDataFileName(""),
174 fDeleteIntermediateFiles(kFALSE),
175 fWriteSelRawData(kFALSE),
176 fStopOnError(kFALSE),
178 fConfigFileName(configFileName),
179 fGAliceFileName("galice.root"),
181 fBkgrdFileNames(NULL),
182 fAlignObjArray(NULL),
183 fUseBkgrdVertex(kTRUE),
184 fRegionOfInterest(kFALSE),
190 fInitCDBCalled(kFALSE),
191 fInitRunNumberCalled(kFALSE),
192 fSetRunNumberFromDataCalled(kFALSE),
193 fEmbeddingFlag(kFALSE),
196 fUseVertexFromCDB(0),
197 fUseMagFieldFromGRP(0),
198 fGRPWriteLocation(Form("local://%s", gSystem->pwd())),
202 fEventSpecie(AliRecoParam::kDefault),
203 fWriteQAExpertData(kTRUE),
207 fWriteGRPEntry(kTRUE)
209 // create simulation object with default parameters
211 SetGAliceFile("galice.root");
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) ;
220 //_____________________________________________________________________________
221 AliSimulation::~AliSimulation()
225 fEventsPerFile.Delete();
226 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
227 // delete fAlignObjArray; fAlignObjArray=0;
229 if (fBkgrdFileNames) {
230 fBkgrdFileNames->Delete();
231 delete fBkgrdFileNames;
234 fSpecCDBUri.Delete();
235 if (fgInstance==this) fgInstance = 0;
237 AliQAManager::QAManager()->ShowQA() ;
238 AliQAManager::Destroy() ;
239 AliCodeTimer::Instance()->Print();
243 //_____________________________________________________________________________
244 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
246 // set the number of events for one run
251 //_____________________________________________________________________________
252 void AliSimulation::InitQA()
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
258 if (fInitCDBCalled) return;
259 fInitCDBCalled = kTRUE;
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() ;
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();
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, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
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("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
286 qam->SetDefaultStorage(fQARefUri);
290 //_____________________________________________________________________________
291 void AliSimulation::InitCDB()
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
297 if (fInitCDBCalled) return;
298 fInitCDBCalled = kTRUE;
300 AliCDBManager* man = AliCDBManager::Instance();
301 if (man->IsDefaultStorageSet())
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();
310 if (fCDBUri.Length() > 0)
312 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
313 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
314 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
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("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
323 man->SetDefaultStorage(fCDBUri);
326 // Now activate the detector specific CDB storage locations
327 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
328 TObject* obj = fSpecCDBUri[i];
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());
338 //_____________________________________________________________________________
339 void AliSimulation::InitRunNumber(){
340 // check run number. If not set, set it to 0 !!!!
342 if (fInitRunNumberCalled) return;
343 fInitRunNumberCalled = kTRUE;
345 AliCDBManager* man = AliCDBManager::Instance();
346 if (man->GetRun() >= 0)
348 AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
349 "Use external variable DC_RUN or AliSimulation::SetRun()!"));
353 AliDebug(2, Form("Setting CDB run number to: %d",fRun));
356 AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
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!
370 ULong_t key = AliCDBManager::Instance()->SetLock(1);
374 //_____________________________________________________________________________
375 void AliSimulation::SetDefaultStorage(const char* uri) {
376 // Store the desired default CDB storage location
377 // Activate it later within the Run() method
383 //_____________________________________________________________________________
384 void AliSimulation::SetQARefDefaultStorage(const char* uri) {
385 // Store the desired default CDB storage location
386 // Activate it later within the Run() method
389 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
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
397 AliCDBPath aPath(calibType);
398 if(!aPath.IsValid()){
399 AliError(Form("Not a valid path: %s", calibType));
403 TObject* obj = fSpecCDBUri.FindObject(calibType);
404 if (obj) fSpecCDBUri.Remove(obj);
405 fSpecCDBUri.Add(new TNamed(calibType, uri));
409 //_____________________________________________________________________________
410 void AliSimulation::SetRunNumber(Int_t run)
413 // Activate it later within the Run() method
418 //_____________________________________________________________________________
419 void AliSimulation::SetSeed(Int_t seed)
422 // Activate it later within the Run() method
427 //_____________________________________________________________________________
428 Bool_t AliSimulation::SetRunNumberFromData()
430 // Set the CDB manager run number
431 // The run number is retrieved from gAlice
433 if (fSetRunNumberFromDataCalled) return kTRUE;
434 fSetRunNumberFromDataCalled = kTRUE;
436 AliCDBManager* man = AliCDBManager::Instance();
437 Int_t runData = -1, runCDB = -1;
439 AliRunLoader* runLoader = LoadRun("READ");
440 if (!runLoader) return kFALSE;
442 runData = runLoader->GetHeader()->GetRun();
446 runCDB = man->GetRun();
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("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
457 man->SetRun(runData);
460 if(man->GetRun() < 0) {
461 AliError("Run number not properly initalized!");
470 //_____________________________________________________________________________
471 void AliSimulation::SetConfigFile(const char* fileName)
473 // set the name of the config file
475 fConfigFileName = fileName;
478 //_____________________________________________________________________________
479 void AliSimulation::SetGAliceFile(const char* fileName)
481 // set the name of the galice file
482 // the path is converted to an absolute one if it is relative
484 fGAliceFileName = fileName;
485 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
486 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
488 fGAliceFileName = absFileName;
489 delete[] absFileName;
492 AliDebug(2, Form("galice file name set to %s", fileName));
495 //_____________________________________________________________________________
496 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
499 // set the number of events per file for the given detector and data type
500 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
502 TNamed* obj = new TNamed(detector, type);
503 obj->SetUniqueID(nEvents);
504 fEventsPerFile.Add(obj);
507 //_____________________________________________________________________________
508 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
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.
518 if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
519 AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
523 // initialize CDB storage, run number, set CDB lock
525 // if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
528 Bool_t delRunLoader = kFALSE;
530 runLoader = LoadRun("READ");
531 if (!runLoader) return kFALSE;
532 delRunLoader = kTRUE;
535 // Export ideal geometry
536 if(!IsGeometryFromFile()) AliGeomManager::GetGeometry()->Export("geometry.root");
538 // Load alignment data from CDB and apply to geometry through AliGeomManager
539 if(fLoadAlignFromCDB){
541 TString detStr = fLoadAlObjsListOfDets;
542 TString loadAlObjsListOfDets = "";
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 += " ";
553 } // end loop over detectors
554 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
555 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
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;
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++) {
575 AliModule* det = (AliModule*) detArray->At(iDet);
576 if (!det || !det->IsActive()) continue;
577 if (IsSelected(det->GetName(), detStr)) {
578 det->UpdateInternalGeometry();
580 } // end loop over detectors
583 if (delRunLoader) delete runLoader;
588 //_____________________________________________________________________________
589 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
591 // add a file with background events for merging
593 TObjString* fileNameStr = new TObjString(fileName);
594 fileNameStr->SetUniqueID(nSignalPerBkgrd);
595 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
596 fBkgrdFileNames->Add(fileNameStr);
599 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
601 // add a file with background events for embeddin
602 MergeWith(fileName, nSignalPerBkgrd);
603 fEmbeddingFlag = kTRUE;
606 //_____________________________________________________________________________
607 Bool_t AliSimulation::Run(Int_t nEvents)
609 // run the generation, simulation and digitization
612 AliCodeTimerAuto("",0)
613 AliSysInfo::AddStamp("Start_Run");
615 // Load run number and seed from environmental vars
616 ProcessEnvironmentVars();
617 AliSysInfo::AddStamp("ProcessEnvironmentVars");
619 gRandom->SetSeed(fSeed);
621 if (nEvents > 0) fNEvents = nEvents;
623 // create and setup the HLT instance
624 if (!fRunHLT.IsNull() && !CreateHLT()) {
625 if (fStopOnError) return kFALSE;
630 // generation and simulation -> hits
631 if (fRunGeneration) {
632 if (!RunSimulation()) if (fStopOnError) return kFALSE;
634 AliSysInfo::AddStamp("RunSimulation");
636 // initialize CDB storage from external environment
637 // (either CDB manager or AliSimulation setters),
638 // if not already done in RunSimulation()
640 AliSysInfo::AddStamp("InitCDB");
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;
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!
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;
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!");
672 if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
674 if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
676 AliSysInfo::AddStamp("MissalignGeometry");
679 // hits -> summable digits
680 AliSysInfo::AddStamp("Start_sdigitization");
681 if (!fMakeSDigits.IsNull()) {
682 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
685 AliSysInfo::AddStamp("Stop_sdigitization");
687 AliSysInfo::AddStamp("Start_digitization");
688 // summable digits -> digits
689 if (!fMakeDigits.IsNull()) {
690 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
691 if (fStopOnError) return kFALSE;
694 AliSysInfo::AddStamp("Stop_digitization");
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()));
706 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
707 if (fStopOnError) return kFALSE;
711 AliSysInfo::AddStamp("Hits2Digits");
715 if (!RunTrigger(fTriggerConfig,fMakeDigits)) {
716 if (fStopOnError) return kFALSE;
719 AliSysInfo::AddStamp("RunTrigger");
722 // digits -> raw data
723 if (!fWriteRawData.IsNull()) {
724 if (!WriteRawData(fWriteRawData, fRawDataFileName,
725 fDeleteIntermediateFiles,fWriteSelRawData)) {
726 if (fStopOnError) return kFALSE;
730 AliSysInfo::AddStamp("WriteRaw");
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()) {
736 if (fStopOnError) return kFALSE;
740 AliSysInfo::AddStamp("RunHLT");
744 Bool_t rv = RunQA() ;
750 AliSysInfo::AddStamp("RunQA");
752 // Cleanup of CDB manager: cache and active storages!
753 AliCDBManager::Instance()->ClearCache();
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)
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
770 // ntheta bins in theta, eta
771 // themin minimum angle in theta (degrees)
772 // themax maximum angle in theta (degrees)
774 // phimin minimum angle in phi (degrees)
775 // phimax maximum angle in phi (degrees)
776 // rmin minimum radius
777 // rmax maximum radius
780 // The number of events generated = ntheta*nphi
781 // run input parameters in macro setup (default="Config.C")
783 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
786 <img src="picts/AliRunLego1.gif">
791 <img src="picts/AliRunLego2.gif">
796 <img src="picts/AliRunLego3.gif">
801 // run the generation and simulation
803 AliCodeTimerAuto("",0)
805 // initialize CDB storage and run number from external environment
806 // (either CDB manager or AliSimulation setters)
812 AliError("no gAlice object. Restart aliroot and try again.");
815 if (gAlice->Modules()->GetEntries() > 0) {
816 AliError("gAlice was already run. Restart aliroot and try again.");
820 AliInfo(Form("initializing gAlice with config file %s",
821 fConfigFileName.Data()));
824 if (nev == -1) nev = nc1 * nc2;
826 // check if initialisation has been done
827 // If runloader has been initialized, set the number of events per file to nc1 * nc2
830 if (!gener) gener = new AliLegoGenerator();
832 // Configure Generator
834 gener->SetRadiusRange(rmin, rmax);
835 gener->SetZMax(zmax);
836 gener->SetCoor1Range(nc1, c1min, c1max);
837 gener->SetCoor2Range(nc2, c2min, c2max);
841 fLego = new AliLego("lego",gener);
843 //__________________________________________________________________________
847 gROOT->LoadMacro(setup);
848 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
850 if(AliCDBManager::Instance()->GetRun() >= 0) {
851 SetRunNumber(AliCDBManager::Instance()->GetRun());
853 AliWarning("Run number not initialized!!");
856 AliRunLoader::Instance()->CdGAFile();
858 AliPDG::AddParticlesToPdgDataBase();
860 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
862 gAlice->GetMCApp()->Init();
865 //Must be here because some MCs (G4) adds detectors here and not in Config.C
866 gAlice->InitLoaders();
867 AliRunLoader::Instance()->MakeTree("E");
870 // Save stuff at the beginning of the file to avoid file corruption
871 AliRunLoader::Instance()->CdGAFile();
874 //Save current generator
875 AliGenerator *gen=gAlice->GetMCApp()->Generator();
876 gAlice->GetMCApp()->ResetGenerator(gener);
877 //Prepare MC for Lego Run
883 AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
884 gMC->ProcessRun(nev);
886 // End of this run, close files
888 // Restore current generator
889 gAlice->GetMCApp()->ResetGenerator(gen);
890 // Delete Lego Object
896 //_____________________________________________________________________________
897 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
901 AliCodeTimerAuto("",0)
903 // initialize CDB storage from external environment
904 // (either CDB manager or AliSimulation setters),
905 // if not already done in RunSimulation()
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;
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!
916 AliRunLoader* runLoader = LoadRun("READ");
917 if (!runLoader) return kFALSE;
918 TString trconfiguration = config;
920 if (trconfiguration.IsNull()) {
921 if(!fTriggerConfig.IsNull()) {
922 trconfiguration = fTriggerConfig;
925 AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
928 runLoader->MakeTree( "GG" );
929 AliCentralTrigger* aCTP = runLoader->GetTrigger();
930 // Load Configuration
931 if (!aCTP->LoadConfiguration( trconfiguration ))
935 if( !aCTP->RunTrigger( runLoader , detectors ) ) {
947 //_____________________________________________________________________________
948 Bool_t AliSimulation::WriteTriggerRawData()
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;
959 //_____________________________________________________________________________
960 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
962 // run the generation and simulation
964 AliCodeTimerAuto("",0)
966 // initialize CDB storage and run number from external environment
967 // (either CDB manager or AliSimulation setters)
968 AliSysInfo::AddStamp("RunSimulation_Begin");
970 AliSysInfo::AddStamp("RunSimulation_InitCDB");
973 AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
976 AliError("no gAlice object. Restart aliroot and try again.");
979 if (gAlice->Modules()->GetEntries() > 0) {
980 AliError("gAlice was already run. Restart aliroot and try again.");
984 AliInfo(Form("initializing gAlice with config file %s",
985 fConfigFileName.Data()));
988 // Initialize ALICE Simulation run
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) {
999 AliInfo("Field is locked now. It cannot be changed in Config.C");
1003 gROOT->LoadMacro(fConfigFileName.Data());
1004 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1005 AliSysInfo::AddStamp("RunSimulation_Config");
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/MeanVertexSPD");
1014 AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1016 vertex->GetXYZ(vtxPos);
1017 vertex->GetSigmaXYZ(vtxSig);
1018 AliInfo("Overwriting Config.C vertex settings !");
1019 AliInfo(Form("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
1020 vtxPos[0], vtxPos[1], vtxPos[2], vtxSig[2]));
1022 AliGenerator *gen = gAlice->GetMCApp()->Generator();
1023 gen->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]); // vertex position
1024 gen->SetSigmaZ(vtxSig[2]);
1028 if(AliCDBManager::Instance()->GetRun() >= 0) {
1029 AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
1030 AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
1032 AliWarning("Run number not initialized!!");
1035 AliRunLoader::Instance()->CdGAFile();
1038 AliPDG::AddParticlesToPdgDataBase();
1040 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1041 AliSysInfo::AddStamp("RunSimulation_GetField");
1043 gAlice->GetMCApp()->Init();
1044 AliSysInfo::AddStamp("RunSimulation_InitMCApp");
1046 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1047 gAlice->InitLoaders();
1048 AliRunLoader::Instance()->MakeTree("E");
1049 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1050 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1051 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1053 // Save stuff at the beginning of the file to avoid file corruption
1054 AliRunLoader::Instance()->CdGAFile();
1056 gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1057 AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1058 //___________________________________________________________________________________________
1060 AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
1062 // Set run number in CDBManager
1063 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
1065 AliRunLoader* runLoader = AliRunLoader::Instance();
1067 AliError(Form("gAlice has no run loader object. "
1068 "Check your config file: %s", fConfigFileName.Data()));
1071 SetGAliceFile(runLoader->GetFileName());
1073 // Misalign geometry
1074 #if ROOT_VERSION_CODE < 331527
1075 AliGeomManager::SetGeometry(gGeoManager);
1077 // Check that the consistency of symbolic names for the activated subdetectors
1078 // in the geometry loaded by AliGeomManager
1079 TString detsToBeChecked = "";
1080 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1081 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1082 AliModule* det = (AliModule*) detArray->At(iDet);
1083 if (!det || !det->IsActive()) continue;
1084 detsToBeChecked += det->GetName();
1085 detsToBeChecked += " ";
1086 } // end loop over detectors
1087 if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
1088 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
1089 MisalignGeometry(runLoader);
1090 AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
1093 // AliRunLoader* runLoader = AliRunLoader::Instance();
1094 // if (!runLoader) {
1095 // AliError(Form("gAlice has no run loader object. "
1096 // "Check your config file: %s", fConfigFileName.Data()));
1099 // SetGAliceFile(runLoader->GetFileName());
1101 if (!gAlice->GetMCApp()->Generator()) {
1102 AliError(Form("gAlice has no generator object. "
1103 "Check your config file: %s", fConfigFileName.Data()));
1107 // Write GRP entry corresponding to the setting found in Cofig.C
1110 AliSysInfo::AddStamp("RunSimulation_WriteGRP");
1112 if (nEvents <= 0) nEvents = fNEvents;
1114 // get vertex from background file in case of merging
1115 if (fUseBkgrdVertex &&
1116 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
1117 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
1118 const char* fileName = ((TObjString*)
1119 (fBkgrdFileNames->At(0)))->GetName();
1120 AliInfo(Form("The vertex will be taken from the background "
1121 "file %s with nSignalPerBackground = %d",
1122 fileName, signalPerBkgrd));
1123 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
1124 gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
1127 if (!fRunSimulation) {
1128 gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
1131 // set the number of events per file for given detectors and data types
1132 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
1133 if (!fEventsPerFile[i]) continue;
1134 const char* detName = fEventsPerFile[i]->GetName();
1135 const char* typeName = fEventsPerFile[i]->GetTitle();
1136 TString loaderName(detName);
1137 loaderName += "Loader";
1138 AliLoader* loader = runLoader->GetLoader(loaderName);
1140 AliError(Form("RunSimulation no loader for %s found\n Number of events per file not set for %s %s",
1141 detName, typeName, detName));
1144 AliDataLoader* dataLoader =
1145 loader->GetDataLoader(typeName);
1147 AliError(Form("no data loader for %s found\n"
1148 "Number of events per file not set for %s %s",
1149 typeName, detName, typeName));
1152 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
1153 AliDebug(1, Form("number of events per file set to %d for %s %s",
1154 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
1157 AliInfo("running gAlice");
1158 AliSysInfo::AddStamp("Start_ProcessRun");
1160 // Create the Root Tree with one branch per detector
1161 //Hits moved to begin event -> now we are crating separate tree for each event
1163 gMC->ProcessRun(nEvents);
1165 // End of this run, close files
1166 if(nEvents>0) FinishRun();
1168 AliSysInfo::AddStamp("Stop_ProcessRun");
1174 //_____________________________________________________________________________
1175 Bool_t AliSimulation::RunSDigitization(const char* detectors)
1177 // run the digitization and produce summable digits
1178 static Int_t eventNr=0;
1179 AliCodeTimerAuto("",0) ;
1181 // initialize CDB storage, run number, set CDB lock
1183 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1186 AliRunLoader* runLoader = LoadRun();
1187 if (!runLoader) return kFALSE;
1189 TString detStr = detectors;
1190 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1191 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1192 AliModule* det = (AliModule*) detArray->At(iDet);
1193 if (!det || !det->IsActive()) continue;
1194 if (IsSelected(det->GetName(), detStr)) {
1195 AliInfo(Form("creating summable digits for %s", det->GetName()));
1196 AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
1197 det->Hits2SDigits();
1198 AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
1199 AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
1203 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1204 AliError(Form("the following detectors were not found: %s",
1206 if (fStopOnError) return kFALSE;
1215 //_____________________________________________________________________________
1216 Bool_t AliSimulation::RunDigitization(const char* detectors,
1217 const char* excludeDetectors)
1219 // run the digitization and produce digits from sdigits
1221 AliCodeTimerAuto("",0)
1223 // initialize CDB storage, run number, set CDB lock
1225 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1228 delete AliRunLoader::Instance();
1233 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1234 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1235 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
1236 // manager->SetEmbeddingFlag(fEmbeddingFlag);
1237 manager->SetInputStream(0, fGAliceFileName.Data());
1238 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1239 const char* fileName = ((TObjString*)
1240 (fBkgrdFileNames->At(iStream-1)))->GetName();
1241 manager->SetInputStream(iStream, fileName);
1244 TString detStr = detectors;
1245 TString detExcl = excludeDetectors;
1246 manager->GetInputStream(0)->ImportgAlice();
1247 AliRunLoader* runLoader =
1248 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
1249 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1250 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1251 AliModule* det = (AliModule*) detArray->At(iDet);
1252 if (!det || !det->IsActive()) continue;
1253 if (IsSelected(det->GetName(), detStr) &&
1254 !IsSelected(det->GetName(), detExcl)) {
1255 AliDigitizer* digitizer = det->CreateDigitizer(manager);
1258 AliError(Form("no digitizer for %s", det->GetName()));
1259 if (fStopOnError) return kFALSE;
1261 digitizer->SetRegionOfInterest(fRegionOfInterest);
1266 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1267 AliError(Form("the following detectors were not found: %s",
1269 if (fStopOnError) return kFALSE;
1272 if (!manager->GetListOfTasks()->IsEmpty()) {
1273 AliInfo("executing digitization");
1282 //_____________________________________________________________________________
1283 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1285 // run the digitization and produce digits from hits
1287 AliCodeTimerAuto("",0)
1289 // initialize CDB storage, run number, set CDB lock
1291 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1294 AliRunLoader* runLoader = LoadRun("READ");
1295 if (!runLoader) return kFALSE;
1297 TString detStr = detectors;
1298 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1299 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1300 AliModule* det = (AliModule*) detArray->At(iDet);
1301 if (!det || !det->IsActive()) continue;
1302 if (IsSelected(det->GetName(), detStr)) {
1303 AliInfo(Form("creating digits from hits for %s", det->GetName()));
1308 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1309 AliError(Form("the following detectors were not found: %s",
1311 if (fStopOnError) return kFALSE;
1317 //_____________________________________________________________________________
1318 Bool_t AliSimulation::WriteRawData(const char* detectors,
1319 const char* fileName,
1320 Bool_t deleteIntermediateFiles,
1323 // convert the digits to raw data
1324 // First DDL raw data files for the given detectors are created.
1325 // If a file name is given, the DDL files are then converted to a DATE file.
1326 // If deleteIntermediateFiles is true, the DDL raw files are deleted
1328 // If the file name has the extension ".root", the DATE file is converted
1330 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1331 // 'selrawdata' flag can be used to enable writing of detectors raw data
1332 // accoring to the trigger cluster.
1334 AliCodeTimerAuto("",0)
1335 AliSysInfo::AddStamp("WriteRawData_Start");
1337 TString detStr = detectors;
1338 if (!WriteRawFiles(detStr.Data())) {
1339 if (fStopOnError) return kFALSE;
1341 AliSysInfo::AddStamp("WriteRawFiles");
1343 // run HLT simulation on simulated DDL raw files
1344 // and produce HLT ddl raw files to be included in date/root file
1345 // bugfix 2009-06-26: the decision whether to write HLT raw data
1346 // is taken in RunHLT. Here HLT always needs to be run in order to
1347 // create HLT digits, unless its switched off. This is due to the
1348 // special placement of the HLT between the generation of DDL files
1349 // and conversion to DATE/Root file.
1350 detStr.ReplaceAll("HLT", "");
1351 if (!fRunHLT.IsNull()) {
1353 if (fStopOnError) return kFALSE;
1356 AliSysInfo::AddStamp("WriteRawData_RunHLT");
1358 TString dateFileName(fileName);
1359 if (!dateFileName.IsNull()) {
1360 Bool_t rootOutput = dateFileName.EndsWith(".root");
1361 if (rootOutput) dateFileName += ".date";
1362 TString selDateFileName;
1364 selDateFileName = "selected.";
1365 selDateFileName+= dateFileName;
1367 if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1368 if (fStopOnError) return kFALSE;
1370 AliSysInfo::AddStamp("ConvertRawFilesToDate");
1371 if (deleteIntermediateFiles) {
1372 AliRunLoader* runLoader = LoadRun("READ");
1373 if (runLoader) for (Int_t iEvent = 0;
1374 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1376 sprintf(command, "rm -r raw%d", iEvent);
1377 gSystem->Exec(command);
1383 if (!ConvertDateToRoot(dateFileName, fileName)) {
1384 if (fStopOnError) return kFALSE;
1386 AliSysInfo::AddStamp("ConvertDateToRoot");
1387 if (deleteIntermediateFiles) {
1388 gSystem->Unlink(dateFileName);
1391 TString selFileName = "selected.";
1392 selFileName += fileName;
1393 if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1394 if (fStopOnError) return kFALSE;
1396 if (deleteIntermediateFiles) {
1397 gSystem->Unlink(selDateFileName);
1406 //_____________________________________________________________________________
1407 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1409 // convert the digits to raw data DDL files
1411 AliCodeTimerAuto("",0)
1413 AliRunLoader* runLoader = LoadRun("READ");
1414 if (!runLoader) return kFALSE;
1416 // write raw data to DDL files
1417 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1418 AliInfo(Form("processing event %d", iEvent));
1419 runLoader->GetEvent(iEvent);
1420 TString baseDir = gSystem->WorkingDirectory();
1422 sprintf(dirName, "raw%d", iEvent);
1423 gSystem->MakeDirectory(dirName);
1424 if (!gSystem->ChangeDirectory(dirName)) {
1425 AliError(Form("couldn't change to directory %s", dirName));
1426 if (fStopOnError) return kFALSE; else continue;
1429 ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1432 TString detStr = detectors;
1433 if (IsSelected("HLT", detStr)) {
1434 // Do nothing. "HLT" will be removed from detStr and HLT raw
1435 // data files are generated in RunHLT.
1438 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1439 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1440 AliModule* det = (AliModule*) detArray->At(iDet);
1441 if (!det || !det->IsActive()) continue;
1442 if (IsSelected(det->GetName(), detStr)) {
1443 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1448 if (!WriteTriggerRawData())
1449 if (fStopOnError) return kFALSE;
1451 gSystem->ChangeDirectory(baseDir);
1452 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1453 AliError(Form("the following detectors were not found: %s",
1455 if (fStopOnError) return kFALSE;
1464 //_____________________________________________________________________________
1465 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1466 const char* selDateFileName)
1468 // convert raw data DDL files to a DATE file with the program "dateStream"
1469 // The second argument is not empty when the user decides to write
1470 // the detectors raw data according to the trigger cluster.
1472 AliCodeTimerAuto("",0)
1474 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1476 AliError("the program dateStream was not found");
1477 if (fStopOnError) return kFALSE;
1482 AliRunLoader* runLoader = LoadRun("READ");
1483 if (!runLoader) return kFALSE;
1485 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1486 Bool_t selrawdata = kFALSE;
1487 if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1490 // Note the option -s. It is used in order to avoid
1491 // the generation of SOR/EOR events.
1492 sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d",
1493 dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1494 FILE* pipe = gSystem->OpenPipe(command, "w");
1497 AliError(Form("Cannot execute command: %s",command));
1501 Int_t selEvents = 0;
1502 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1504 UInt_t detectorPattern = 0;
1505 runLoader->GetEvent(iEvent);
1506 if (!runLoader->LoadTrigger()) {
1507 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1508 detectorPattern = aCTP->GetClusterMask();
1509 // Check if the event was triggered by CTP
1511 if (aCTP->GetClassMask()) selEvents++;
1515 AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1517 AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1518 selrawdata = kFALSE;
1522 fprintf(pipe, "GDC DetectorPattern %u\n", detectorPattern);
1526 // loop over detectors and DDLs
1527 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1528 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1530 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1531 Int_t ldcID = Int_t(ldc + 0.0001);
1532 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1534 char rawFileName[256];
1535 sprintf(rawFileName, "raw%d/%s",
1536 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1538 // check existence and size of raw data file
1539 FILE* file = fopen(rawFileName, "rb");
1540 if (!file) continue;
1541 fseek(file, 0, SEEK_END);
1542 unsigned long size = ftell(file);
1544 if (!size) continue;
1546 if (ldcID != prevLDC) {
1547 fprintf(pipe, " LDC Id %d\n", ldcID);
1550 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1555 Int_t result = gSystem->ClosePipe(pipe);
1557 if (!(selrawdata && selEvents > 0)) {
1559 return (result == 0);
1562 AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1564 sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d",
1565 selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1566 FILE* pipe2 = gSystem->OpenPipe(command, "w");
1568 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1570 // Get the trigger decision and cluster
1571 UInt_t detectorPattern = 0;
1573 runLoader->GetEvent(iEvent);
1574 if (!runLoader->LoadTrigger()) {
1575 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1576 if (aCTP->GetClassMask() == 0) continue;
1577 detectorPattern = aCTP->GetClusterMask();
1578 detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1579 AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1582 fprintf(pipe2, "GDC DetectorPattern %u\n", detectorPattern);
1586 // loop over detectors and DDLs
1587 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1588 // Write only raw data from detectors that
1589 // are contained in the trigger cluster(s)
1590 if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1592 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1594 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1595 Int_t ldcID = Int_t(ldc + 0.0001);
1596 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1598 char rawFileName[256];
1599 sprintf(rawFileName, "raw%d/%s",
1600 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1602 // check existence and size of raw data file
1603 FILE* file = fopen(rawFileName, "rb");
1604 if (!file) continue;
1605 fseek(file, 0, SEEK_END);
1606 unsigned long size = ftell(file);
1608 if (!size) continue;
1610 if (ldcID != prevLDC) {
1611 fprintf(pipe2, " LDC Id %d\n", ldcID);
1614 fprintf(pipe2, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1619 Int_t result2 = gSystem->ClosePipe(pipe2);
1622 return ((result == 0) && (result2 == 0));
1625 //_____________________________________________________________________________
1626 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1627 const char* rootFileName)
1629 // convert a DATE file to a root file with the program "alimdc"
1632 const Int_t kDBSize = 2000000000;
1633 const Int_t kTagDBSize = 1000000000;
1634 const Bool_t kFilter = kFALSE;
1635 const Int_t kCompression = 1;
1637 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1639 AliError("the program alimdc was not found");
1640 if (fStopOnError) return kFALSE;
1645 AliInfo(Form("converting DATE file %s to root file %s",
1646 dateFileName, rootFileName));
1648 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1649 const char* tagDBFS = "/tmp/mdc1/tags";
1651 // User defined file system locations
1652 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1653 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1654 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1655 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1656 if (gSystem->Getenv("ALIMDC_TAGDB"))
1657 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1659 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1660 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1661 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1663 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1664 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1665 gSystem->Exec(Form("mkdir %s",tagDBFS));
1667 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1668 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1669 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1671 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1672 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1673 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1675 return (result == 0);
1679 //_____________________________________________________________________________
1680 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1682 // delete existing run loaders, open a new one and load gAlice
1684 delete AliRunLoader::Instance();
1685 AliRunLoader* runLoader =
1686 AliRunLoader::Open(fGAliceFileName.Data(),
1687 AliConfig::GetDefaultEventFolderName(), mode);
1689 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1692 runLoader->LoadgAlice();
1693 runLoader->LoadHeader();
1694 gAlice = runLoader->GetAliRun();
1696 AliError(Form("no gAlice object found in file %s",
1697 fGAliceFileName.Data()));
1703 //_____________________________________________________________________________
1704 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1706 // get or calculate the number of signal events per background event
1708 if (!fBkgrdFileNames) return 1;
1709 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1710 if (nBkgrdFiles == 0) return 1;
1712 // get the number of signal events
1714 AliRunLoader* runLoader =
1715 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1716 if (!runLoader) return 1;
1718 nEvents = runLoader->GetNumberOfEvents();
1723 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1724 // get the number of background events
1725 const char* fileName = ((TObjString*)
1726 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1727 AliRunLoader* runLoader =
1728 AliRunLoader::Open(fileName, "BKGRD");
1729 if (!runLoader) continue;
1730 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1733 // get or calculate the number of signal per background events
1734 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1735 if (nSignalPerBkgrd <= 0) {
1736 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1737 } else if (result && (result != nSignalPerBkgrd)) {
1738 AliInfo(Form("the number of signal events per background event "
1739 "will be changed from %d to %d for stream %d",
1740 nSignalPerBkgrd, result, iBkgrdFile+1));
1741 nSignalPerBkgrd = result;
1744 if (!result) result = nSignalPerBkgrd;
1745 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1746 AliWarning(Form("not enough background events (%d) for %d signal events "
1747 "using %d signal per background events for stream %d",
1748 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1755 //_____________________________________________________________________________
1756 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1758 // check whether detName is contained in detectors
1759 // if yes, it is removed from detectors
1761 // check if all detectors are selected
1762 if ((detectors.CompareTo("ALL") == 0) ||
1763 detectors.BeginsWith("ALL ") ||
1764 detectors.EndsWith(" ALL") ||
1765 detectors.Contains(" ALL ")) {
1770 // search for the given detector
1771 Bool_t result = kFALSE;
1772 if ((detectors.CompareTo(detName) == 0) ||
1773 detectors.BeginsWith(detName+" ") ||
1774 detectors.EndsWith(" "+detName) ||
1775 detectors.Contains(" "+detName+" ")) {
1776 detectors.ReplaceAll(detName, "");
1780 // clean up the detectors string
1781 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1782 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1783 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1788 //_____________________________________________________________________________
1789 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName, Int_t N)
1792 // Steering routine to convert raw data in directory rawDirectory/ to fake SDigits.
1793 // These can be used for embedding of MC tracks into RAW data using the standard
1794 // merging procedure.
1796 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1799 AliError("no gAlice object. Restart aliroot and try again.");
1802 if (gAlice->Modules()->GetEntries() > 0) {
1803 AliError("gAlice was already run. Restart aliroot and try again.");
1807 AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1811 gROOT->LoadMacro(fConfigFileName.Data());
1812 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1814 if(AliCDBManager::Instance()->GetRun() >= 0) {
1815 SetRunNumber(AliCDBManager::Instance()->GetRun());
1817 AliWarning("Run number not initialized!!");
1820 AliRunLoader::Instance()->CdGAFile();
1822 AliPDG::AddParticlesToPdgDataBase();
1824 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1826 gAlice->GetMCApp()->Init();
1828 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1829 gAlice->InitLoaders();
1830 AliRunLoader::Instance()->MakeTree("E");
1831 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1832 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1833 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1836 // Save stuff at the beginning of the file to avoid file corruption
1837 AliRunLoader::Instance()->CdGAFile();
1842 //AliCDBManager* man = AliCDBManager::Instance();
1843 //man->SetRun(0); // Should this come from rawdata header ?
1847 // Get the runloader
1848 AliRunLoader* runLoader = AliRunLoader::Instance();
1850 // Open esd file if available
1853 AliESDEvent* esd = 0;
1854 if (esdFileName && (strlen(esdFileName)>0)) {
1855 esdFile = TFile::Open(esdFileName);
1857 esd = new AliESDEvent();
1858 esdFile->GetObject("esdTree", treeESD);
1859 if (treeESD) esd->ReadFromTree(treeESD);
1864 // Create the RawReader
1865 TString fileName(rawDirectory);
1866 AliRawReader* rawReader = 0x0;
1867 if (fileName.EndsWith("/")) {
1868 rawReader = new AliRawReaderFile(fileName);
1869 } else if (fileName.EndsWith(".root")) {
1870 rawReader = new AliRawReaderRoot(fileName);
1871 } else if (!fileName.IsNull()) {
1872 rawReader = new AliRawReaderDate(fileName);
1874 // if (!fEquipIdMap.IsNull() && fRawReader)
1875 // fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1877 // Get list of detectors
1878 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1881 AliHeader* header = runLoader->GetHeader();
1885 if (!(rawReader->NextEvent())) break;
1886 runLoader->SetEventNumber(nev);
1887 runLoader->GetHeader()->Reset(rawReader->GetRunNumber(),
1889 runLoader->TreeE()->Fill();
1890 runLoader->GetEvent(nev);
1891 AliInfo(Form("We are at event %d",nev));
1894 TString detStr = fMakeSDigits;
1895 for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1896 AliModule* det = (AliModule*) detArray->At(iDet);
1897 if (!det || !det->IsActive()) continue;
1898 if (IsSelected(det->GetName(), detStr)) {
1899 AliInfo(Form("Calling Raw2SDigits for %s", det->GetName()));
1900 det->Raw2SDigits(rawReader);
1907 // If ESD information available obtain reconstructed vertex and store in header.
1909 treeESD->GetEvent(nev);
1910 const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1911 Double_t position[3];
1912 esdVertex->GetXYZ(position);
1913 AliGenEventHeader* mcHeader = new AliGenEventHeader("ESD");
1916 for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1917 mcHeader->SetPrimaryVertex(mcV);
1918 header->Reset(0,nev);
1919 header->SetGenEventHeader(mcHeader);
1920 AliInfo(Form("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]));
1924 AliInfo(Form("Finished event %d",nev));
1933 runLoader->CdGAFile();
1934 runLoader->WriteHeader("OVERWRITE");
1935 runLoader->WriteRunLoader();
1940 //_____________________________________________________________________________
1941 void AliSimulation::FinishRun()
1944 // Called at the end of the run.
1949 AliDebug(1, "Finish Lego");
1950 AliRunLoader::Instance()->CdGAFile();
1954 // Clean detector information
1955 TIter next(gAlice->Modules());
1956 AliModule *detector;
1957 while((detector = dynamic_cast<AliModule*>(next()))) {
1958 AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
1959 detector->FinishRun();
1962 AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
1963 AliRunLoader::Instance()->WriteHeader("OVERWRITE");
1965 // Write AliRun info and all detectors parameters
1966 AliRunLoader::Instance()->CdGAFile();
1967 gAlice->Write(0,TObject::kOverwrite);//write AliRun
1968 AliRunLoader::Instance()->Write(0,TObject::kOverwrite);//write RunLoader itself
1970 if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();
1971 AliRunLoader::Instance()->Synchronize();
1974 //_____________________________________________________________________________
1975 Int_t AliSimulation::GetDetIndex(const char* detector)
1977 // return the detector index corresponding to detector
1979 for (index = 0; index < fgkNDetectors ; index++) {
1980 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1986 //_____________________________________________________________________________
1987 Bool_t AliSimulation::CreateHLT()
1989 // Init the HLT simulation.
1990 // The function loads the library and creates the instance of AliHLTSimulation.
1991 // the main reason for the decoupled creation is to set the transient OCDB
1992 // objects before the OCDB is locked
1994 // load the library dynamically
1995 gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1997 // check for the library version
1998 AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
2000 AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
2003 if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
2004 AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
2007 // print compile info
2008 typedef void (*CompileInfo)( const char*& date, const char*& time);
2009 CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
2011 const char* date="";
2012 const char* time="";
2013 (*fctInfo)(date, time);
2014 if (!date) date="unknown";
2015 if (!time) time="unknown";
2016 AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
2018 AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
2021 // create instance of the HLT simulation
2022 AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
2023 if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
2024 AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
2028 TString specObjects;
2029 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
2030 if (specObjects.Length()>0) specObjects+=" ";
2031 specObjects+=fSpecCDBUri[i]->GetName();
2034 AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
2035 if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
2036 AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
2042 //_____________________________________________________________________________
2043 Bool_t AliSimulation::RunHLT()
2045 // Run the HLT simulation
2046 // HLT simulation is implemented in HLT/sim/AliHLTSimulation
2047 // Disabled if fRunHLT is empty, default vaule is "default".
2048 // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
2049 // The default simulation depends on the HLT component libraries and their
2050 // corresponding agents which define components and chains to run. See
2051 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
2052 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
2054 // The libraries to be loaded can be specified as an option.
2056 // AliSimulation sim;
2057 // sim.SetRunHLT("libAliHLTSample.so");
2059 // will only load <tt>libAliHLTSample.so</tt>
2061 // Other available options:
2062 // \li loglevel=<i>level</i> <br>
2063 // logging level for this processing
2065 // disable redirection of log messages to AliLog class
2066 // \li config=<i>macro</i>
2067 // configuration macro
2068 // \li chains=<i>configuration</i>
2069 // comma separated list of configurations to be run during simulation
2070 // \li rawfile=<i>file</i>
2071 // source for the RawReader to be created, the default is <i>./</i> if
2072 // raw data is simulated
2076 if (!fpHLT && !CreateHLT()) {
2079 AliHLTSimulation* pHLT=fpHLT;
2081 AliRunLoader* pRunLoader = LoadRun("READ");
2082 if (!pRunLoader) return kFALSE;
2084 // initialize CDB storage, run number, set CDB lock
2085 // thats for the case of running HLT simulation without all the other steps
2086 // multiple calls are handled by the function, so we can just call
2088 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
2091 // init the HLT simulation
2093 if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2094 TString detStr = fWriteRawData;
2095 if (!IsSelected("HLT", detStr)) {
2096 options+=" writerawfiles=";
2098 options+=" writerawfiles=HLT";
2101 if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2102 // as a matter of fact, HLT will run reconstruction and needs the RawReader
2103 // in order to get detector data. By default, RawReaderFile is used to read
2104 // the already simulated ddl files. Date and Root files from the raw data
2105 // are generated after the HLT simulation.
2106 options+=" rawfile=./";
2109 AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2110 if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2111 AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2113 // run the HLT simulation
2114 AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2115 if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2116 AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2120 // delete the instance
2121 AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2122 if (fctDelete==NULL || fctDelete(pHLT)<0) {
2123 AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2127 return iResult>=0?kTRUE:kFALSE;
2130 //_____________________________________________________________________________
2131 Bool_t AliSimulation::RunQA()
2133 // run the QA on summable hits, digits or digits
2135 if(!gAlice) return kFALSE;
2136 AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2138 TString detectorsw("") ;
2140 AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2141 detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ;
2142 if ( detectorsw.IsNull() )
2147 //_____________________________________________________________________________
2148 Bool_t AliSimulation::SetRunQA(TString detAndAction)
2150 // Allows to run QA for a selected set of detectors
2151 // and a selected set of tasks among HITS, SDIGITS and DIGITS
2152 // all selected detectors run the same selected tasks
2154 if (!detAndAction.Contains(":")) {
2155 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2159 Int_t colon = detAndAction.Index(":") ;
2160 fQADetectors = detAndAction(0, colon) ;
2161 if (fQADetectors.Contains("ALL") ){
2162 TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2163 Int_t minus = fQADetectors.Last('-') ;
2164 TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2165 TString toRemove("") ;
2166 while (minus >= 0) {
2167 toRemove = fQADetectors(minus+1, fQADetectors.Length()) ;
2168 toRemove = toRemove.Strip() ;
2169 toKeep.ReplaceAll(toRemove, "") ;
2170 fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ;
2171 minus = fQADetectors.Last('-') ;
2173 fQADetectors = toKeep ;
2175 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2176 if (fQATasks.Contains("ALL") ) {
2177 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
2179 fQATasks.ToUpper() ;
2181 if ( fQATasks.Contains("HIT") )
2182 tempo = Form("%d ", AliQAv1::kHITS) ;
2183 if ( fQATasks.Contains("SDIGIT") )
2184 tempo += Form("%d ", AliQAv1::kSDIGITS) ;
2185 if ( fQATasks.Contains("DIGIT") )
2186 tempo += Form("%d ", AliQAv1::kDIGITS) ;
2188 if (fQATasks.IsNull()) {
2189 AliInfo("No QA requested\n") ;
2194 TString tempo(fQATasks) ;
2195 tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS)) ;
2196 tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;
2197 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;
2198 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2200 AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ;
2201 AliQAManager::QAManager()->SetTasks(fQATasks) ;
2202 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++)
2203 AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2208 //_____________________________________________________________________________
2209 void AliSimulation::ProcessEnvironmentVars()
2211 // Extract run number and random generator seed from env variables
2213 AliInfo("Processing environment variables");
2215 // Random Number seed
2217 // first check that seed is not already set
2219 if (gSystem->Getenv("CONFIG_SEED")) {
2220 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2223 if (gSystem->Getenv("CONFIG_SEED")) {
2224 AliInfo(Form("Seed for random number generation already set (%d)"
2225 ": CONFIG_SEED variable ignored!", fSeed));
2229 AliInfo(Form("Seed for random number generation = %d ", fSeed));
2233 // first check that run number is not already set
2235 if (gSystem->Getenv("DC_RUN")) {
2236 fRun = atoi(gSystem->Getenv("DC_RUN"));
2239 if (gSystem->Getenv("DC_RUN")) {
2240 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2244 AliInfo(Form("Run number = %d", fRun));
2247 //---------------------------------------------------------------------
2248 void AliSimulation::WriteGRPEntry()
2250 // Get the necessary information from galice (generator, trigger etc) and
2251 // write a GRP entry corresponding to the settings in the Config.C used
2252 // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2255 AliInfo("Writing global run parameters entry into the OCDB");
2257 AliGRPObject* grpObj = new AliGRPObject();
2259 grpObj->SetRunType("PHYSICS");
2260 grpObj->SetTimeStart(0);
2262 grpObj->SetTimeStart(0);
2263 grpObj->SetTimeEnd(curtime.Convert());
2264 grpObj->SetBeamEnergyIsSqrtSHalfGeV(); // new format of GRP: store sqrt(s)/2 in GeV
2266 const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2268 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
2271 gen->GetProjectile(projectile,a,z);
2273 gen->GetTarget(target,a,z);
2274 TString beamType = projectile + "-" + target;
2275 beamType.ReplaceAll(" ","");
2276 if (!beamType.CompareTo("-")) {
2277 grpObj->SetBeamType("UNKNOWN");
2280 grpObj->SetBeamType(beamType);
2281 // Heavy ion run, the event specie is set to kHighMult
2282 fEventSpecie = AliRecoParam::kHighMult;
2283 if ((strcmp(beamType,"p-p") == 0) ||
2284 (strcmp(beamType,"p-") == 0) ||
2285 (strcmp(beamType,"-p") == 0) ||
2286 (strcmp(beamType,"P-P") == 0) ||
2287 (strcmp(beamType,"P-") == 0) ||
2288 (strcmp(beamType,"-P") == 0)) {
2289 // Proton run, the event specie is set to kLowMult
2290 fEventSpecie = AliRecoParam::kLowMult;
2294 AliWarning("Unknown beam type and energy! Setting energy to 0");
2295 grpObj->SetBeamEnergy(0);
2296 grpObj->SetBeamType("UNKNOWN");
2299 UInt_t detectorPattern = 0;
2301 TObjArray *detArray = gAlice->Detectors();
2302 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2303 if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2304 detectorPattern |= (1 << iDet);
2309 if (!fTriggerConfig.IsNull())
2310 detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2313 if (!fRunHLT.IsNull())
2314 detectorPattern |= (1 << AliDAQ::kHLTId);
2316 grpObj->SetNumberOfDetectors((Char_t)nDets);
2317 grpObj->SetDetectorMask((Int_t)detectorPattern);
2318 grpObj->SetLHCPeriod("LHC08c");
2319 grpObj->SetLHCState("STABLE_BEAMS");
2321 AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2322 Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2324 Float_t factorSol = field ? field->GetFactorSol() : 0;
2325 Float_t currentSol = TMath::Abs(factorSol)>1E-6 ?
2326 TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
2328 Float_t factorDip = field ? field->GetFactorDip() : 0;
2329 Float_t currentDip = 6000.*TMath::Abs(factorDip);
2331 grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
2332 grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);
2333 grpObj->SetL3Polarity(factorSol>0 ? 0:1);
2334 grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
2335 grpObj->SetUniformBMap(field->IsUniform()); // for special MC with k5kGUniform map
2336 grpObj->SetPolarityConventionLHC(); // LHC convention +/+ current -> -/- field main components
2338 grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2340 //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2342 // Now store the entry in OCDB
2343 AliCDBManager* man = AliCDBManager::Instance();
2345 man->SetLock(0, fKey);
2347 AliCDBStorage* sto = man->GetStorage(fGRPWriteLocation.Data());
2350 AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2351 AliCDBMetaData *metadata= new AliCDBMetaData();
2353 metadata->SetResponsible("alice-off@cern.ch");
2354 metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2356 sto->Put(grpObj,id,metadata);
2357 man->SetLock(1, fKey);