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 ///////////////////////////////////////////////////////////////////////////////
108 #include <TVirtualMCApplication.h>
109 #include <TGeoManager.h>
110 #include <TObjString.h>
115 #include "AliCodeTimer.h"
116 #include "AliCDBStorage.h"
117 #include "AliCDBEntry.h"
118 #include "AliCDBManager.h"
119 #include "AliGeomManager.h"
120 #include "AliAlignObj.h"
121 #include "AliCentralTrigger.h"
123 #include "AliDigitizer.h"
124 #include "AliGenerator.h"
126 #include "AliModule.h"
128 #include "AliRunDigitizer.h"
129 #include "AliRunLoader.h"
130 #include "AliSimulation.h"
131 #include "AliVertexGenFile.h"
132 #include "AliCentralTrigger.h"
133 #include "AliCTPRawData.h"
134 #include "AliRawReaderFile.h"
135 #include "AliRawReaderRoot.h"
136 #include "AliRawReaderDate.h"
138 #include "AliHeader.h"
139 #include "AliGenEventHeader.h"
141 #include "AliHLTSimulation.h"
142 #include "AliSysInfo.h"
144 #include "AliGRPObject.h"
146 ClassImp(AliSimulation)
148 AliSimulation *AliSimulation::fgInstance = 0;
149 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
151 //_____________________________________________________________________________
152 AliSimulation::AliSimulation(const char* configFileName,
153 const char* name, const char* title) :
156 fRunGeneration(kTRUE),
157 fRunSimulation(kTRUE),
158 fLoadAlignFromCDB(kTRUE),
159 fLoadAlObjsListOfDets("ALL"),
163 fMakeDigitsFromHits(""),
165 fRawDataFileName(""),
166 fDeleteIntermediateFiles(kFALSE),
167 fWriteSelRawData(kFALSE),
168 fStopOnError(kFALSE),
170 fConfigFileName(configFileName),
171 fGAliceFileName("galice.root"),
173 fBkgrdFileNames(NULL),
174 fAlignObjArray(NULL),
175 fUseBkgrdVertex(kTRUE),
176 fRegionOfInterest(kFALSE),
181 fInitCDBCalled(kFALSE),
182 fInitRunNumberCalled(kFALSE),
183 fSetRunNumberFromDataCalled(kFALSE),
184 fEmbeddingFlag(kFALSE),
190 fWriteGRPEntry(kTRUE)
192 // create simulation object with default parameters
194 SetGAliceFile("galice.root");
197 fQASteer = new AliQADataMakerSteer("sim") ;
198 fQASteer->SetActiveDetectors(fQADetectors) ;
199 fQATasks = Form("%d %d %d", AliQA::kHITS, AliQA::kSDIGITS, AliQA::kDIGITS) ;
200 fQASteer->SetTasks(fQATasks) ;
203 //_____________________________________________________________________________
204 AliSimulation::~AliSimulation()
208 fEventsPerFile.Delete();
209 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
210 // delete fAlignObjArray; fAlignObjArray=0;
212 if (fBkgrdFileNames) {
213 fBkgrdFileNames->Delete();
214 delete fBkgrdFileNames;
217 fSpecCDBUri.Delete();
218 if (fgInstance==this) fgInstance = 0;
222 AliCodeTimer::Instance()->Print();
226 //_____________________________________________________________________________
227 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
229 // set the number of events for one run
234 //_____________________________________________________________________________
235 void AliSimulation::InitCDB()
237 // activate a default CDB storage
238 // First check if we have any CDB storage set, because it is used
239 // to retrieve the calibration and alignment constants
241 if (fInitCDBCalled) return;
242 fInitCDBCalled = kTRUE;
244 AliCDBManager* man = AliCDBManager::Instance();
245 if (man->IsDefaultStorageSet())
247 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
248 AliWarning("Default CDB storage has been already set !");
249 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
250 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
251 fCDBUri = man->GetDefaultStorage()->GetURI();
254 if (fCDBUri.Length() > 0)
256 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
257 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
258 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
260 fCDBUri="local://$ALICE_ROOT";
261 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
262 AliWarning("Default CDB storage not yet set !!!!");
263 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
264 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
267 man->SetDefaultStorage(fCDBUri);
270 // Now activate the detector specific CDB storage locations
271 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
272 TObject* obj = fSpecCDBUri[i];
274 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
275 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
276 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
277 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
282 //_____________________________________________________________________________
283 void AliSimulation::InitRunNumber(){
284 // check run number. If not set, set it to 0 !!!!
286 if (fInitRunNumberCalled) return;
287 fInitRunNumberCalled = kTRUE;
289 AliCDBManager* man = AliCDBManager::Instance();
290 if (man->GetRun() >= 0)
292 AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
293 "Use external variable DC_RUN or AliSimulation::SetRun()!"));
297 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
298 AliDebug(2, Form("Setting CDB run number to: %d",fRun));
299 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
302 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
303 AliWarning("Run number not yet set !!!!");
304 AliWarning(Form("Setting it now to: %d", fRun));
305 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
314 //_____________________________________________________________________________
315 void AliSimulation::SetCDBLock() {
316 // Set CDB lock: from now on it is forbidden to reset the run number
317 // or the default storage or to activate any further storage!
319 AliCDBManager::Instance()->SetLock(1);
322 //_____________________________________________________________________________
323 void AliSimulation::SetDefaultStorage(const char* uri) {
324 // Store the desired default CDB storage location
325 // Activate it later within the Run() method
331 //_____________________________________________________________________________
332 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
333 // Store a detector-specific CDB storage location
334 // Activate it later within the Run() method
336 AliCDBPath aPath(calibType);
337 if(!aPath.IsValid()){
338 AliError(Form("Not a valid path: %s", calibType));
342 TObject* obj = fSpecCDBUri.FindObject(calibType);
343 if (obj) fSpecCDBUri.Remove(obj);
344 fSpecCDBUri.Add(new TNamed(calibType, uri));
348 //_____________________________________________________________________________
349 void AliSimulation::SetRunNumber(Int_t run)
352 // Activate it later within the Run() method
357 //_____________________________________________________________________________
358 void AliSimulation::SetSeed(Int_t seed)
361 // Activate it later within the Run() method
366 //_____________________________________________________________________________
367 Bool_t AliSimulation::SetRunNumberFromData()
369 // Set the CDB manager run number
370 // The run number is retrieved from gAlice
372 if (fSetRunNumberFromDataCalled) return kTRUE;
373 fSetRunNumberFromDataCalled = kTRUE;
375 AliCDBManager* man = AliCDBManager::Instance();
376 Int_t runData = -1, runCDB = -1;
378 AliRunLoader* runLoader = LoadRun("READ");
379 if (!runLoader) return kFALSE;
381 runData = runLoader->GetHeader()->GetRun();
385 runCDB = man->GetRun();
387 if (runCDB != runData) {
388 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
389 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
390 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
391 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
396 man->SetRun(runData);
399 if(man->GetRun() < 0) {
400 AliError("Run number not properly initalized!");
409 //_____________________________________________________________________________
410 void AliSimulation::SetConfigFile(const char* fileName)
412 // set the name of the config file
414 fConfigFileName = fileName;
417 //_____________________________________________________________________________
418 void AliSimulation::SetGAliceFile(const char* fileName)
420 // set the name of the galice file
421 // the path is converted to an absolute one if it is relative
423 fGAliceFileName = fileName;
424 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
425 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
427 fGAliceFileName = absFileName;
428 delete[] absFileName;
431 AliDebug(2, Form("galice file name set to %s", fileName));
434 //_____________________________________________________________________________
435 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
438 // set the number of events per file for the given detector and data type
439 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
441 TNamed* obj = new TNamed(detector, type);
442 obj->SetUniqueID(nEvents);
443 fEventsPerFile.Add(obj);
446 //_____________________________________________________________________________
447 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
449 // Read the alignment objects from CDB.
450 // Each detector is supposed to have the
451 // alignment objects in DET/Align/Data CDB path.
452 // All the detector objects are then collected,
453 // sorted by geometry level (starting from ALIC) and
454 // then applied to the TGeo geometry.
455 // Finally an overlaps check is performed.
457 if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
458 AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
462 // initialize CDB storage, run number, set CDB lock
464 // if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
467 Bool_t delRunLoader = kFALSE;
469 runLoader = LoadRun("READ");
470 if (!runLoader) return kFALSE;
471 delRunLoader = kTRUE;
474 // Export ideal geometry
475 if(!gAlice->IsRootGeometry()) AliGeomManager::GetGeometry()->Export("geometry.root");
477 // Load alignment data from CDB and apply to geometry through AliGeomManager
478 if(fLoadAlignFromCDB){
480 TString detStr = fLoadAlObjsListOfDets;
481 TString loadAlObjsListOfDets = "";
483 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
484 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
485 AliModule* det = (AliModule*) detArray->At(iDet);
486 if (!det || !det->IsActive()) continue;
487 if (IsSelected(det->GetName(), detStr)) {
488 //add det to list of dets to be aligned from CDB
489 loadAlObjsListOfDets += det->GetName();
490 loadAlObjsListOfDets += " ";
492 } // end loop over detectors
493 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
494 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
496 // Check if the array with alignment objects was
497 // provided by the user. If yes, apply the objects
498 // to the present TGeo geometry
499 if (fAlignObjArray) {
500 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
501 AliError("The misalignment of one or more volumes failed!"
502 "Compare the list of simulated detectors and the list of detector alignment data!");
503 if (delRunLoader) delete runLoader;
509 // Update the internal geometry of modules (ITS needs it)
510 TString detStr = fLoadAlObjsListOfDets;
511 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
512 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
514 AliModule* det = (AliModule*) detArray->At(iDet);
515 if (!det || !det->IsActive()) continue;
516 if (IsSelected(det->GetName(), detStr)) {
517 det->UpdateInternalGeometry();
519 } // end loop over detectors
522 if (delRunLoader) delete runLoader;
527 //_____________________________________________________________________________
528 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
530 // add a file with background events for merging
532 TObjString* fileNameStr = new TObjString(fileName);
533 fileNameStr->SetUniqueID(nSignalPerBkgrd);
534 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
535 fBkgrdFileNames->Add(fileNameStr);
538 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
540 // add a file with background events for embeddin
541 MergeWith(fileName, nSignalPerBkgrd);
542 fEmbeddingFlag = kTRUE;
545 //_____________________________________________________________________________
546 Bool_t AliSimulation::Run(Int_t nEvents)
548 // run the generation, simulation and digitization
553 // Load run number and seed from environmental vars
554 ProcessEnvironmentVars();
556 gRandom->SetSeed(fSeed);
558 if (nEvents > 0) fNEvents = nEvents;
560 // generation and simulation -> hits
561 if (fRunGeneration) {
562 if (!RunSimulation()) if (fStopOnError) return kFALSE;
565 // initialize CDB storage from external environment
566 // (either CDB manager or AliSimulation setters),
567 // if not already done in RunSimulation()
570 // Set run number in CDBManager from data
571 // From this point on the run number must be always loaded from data!
572 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
574 // Set CDB lock: from now on it is forbidden to reset the run number
575 // or the default storage or to activate any further storage!
578 // If RunSimulation was not called, load the geometry and misalign it
579 if (!AliGeomManager::GetGeometry()) {
580 // Initialize the geometry manager
581 AliGeomManager::LoadGeometry("geometry.root");
583 // // Check that the consistency of symbolic names for the activated subdetectors
584 // // in the geometry loaded by AliGeomManager
585 // AliRunLoader* runLoader = LoadRun("READ");
586 // if (!runLoader) return kFALSE;
588 // TString detsToBeChecked = "";
589 // TObjArray* detArray = runLoader->GetAliRun()->Detectors();
590 // for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
591 // AliModule* det = (AliModule*) detArray->At(iDet);
592 // if (!det || !det->IsActive()) continue;
593 // detsToBeChecked += det->GetName();
594 // detsToBeChecked += " ";
595 // } // end loop over detectors
596 // if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
597 if(!AliGeomManager::CheckSymNamesLUT("ALL"))
598 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
600 if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
602 if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
606 // hits -> summable digits
607 AliSysInfo::AddStamp("Start_sdigitization");
608 if (!fMakeSDigits.IsNull()) {
609 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
612 AliSysInfo::AddStamp("Stop_sdigitization");
614 AliSysInfo::AddStamp("Start_digitization");
615 // summable digits -> digits
616 if (!fMakeDigits.IsNull()) {
617 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
618 if (fStopOnError) return kFALSE;
621 AliSysInfo::AddStamp("Stop_digitization");
626 if (!fMakeDigitsFromHits.IsNull()) {
627 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
628 AliWarning(Form("Merging and direct creation of digits from hits "
629 "was selected for some detectors. "
630 "No merging will be done for the following detectors: %s",
631 fMakeDigitsFromHits.Data()));
633 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
634 if (fStopOnError) return kFALSE;
641 if (!RunTrigger(fMakeTrigger,fMakeDigits)) {
642 if (fStopOnError) return kFALSE;
647 // digits -> raw data
648 if (!fWriteRawData.IsNull()) {
649 if (!WriteRawData(fWriteRawData, fRawDataFileName,
650 fDeleteIntermediateFiles,fWriteSelRawData)) {
651 if (fStopOnError) return kFALSE;
656 // run HLT simulation on simulated digit data if raw data is not
657 // simulated, otherwise its called as part of WriteRawData
658 if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
660 if (fStopOnError) return kFALSE;
666 Bool_t rv = RunQA() ;
672 // Cleanup of CDB manager: cache and active storages!
673 AliCDBManager::Instance()->ClearCache();
678 //_____________________________________________________________________________
679 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
685 // initialize CDB storage from external environment
686 // (either CDB manager or AliSimulation setters),
687 // if not already done in RunSimulation()
690 // Set run number in CDBManager from data
691 // From this point on the run number must be always loaded from data!
692 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
694 // Set CDB lock: from now on it is forbidden to reset the run number
695 // or the default storage or to activate any further storage!
698 AliRunLoader* runLoader = LoadRun("READ");
699 if (!runLoader) return kFALSE;
700 TString trconfiguration = config;
702 if (trconfiguration.IsNull()) {
703 if (strcmp(gAlice->GetTriggerDescriptor(),"")) {
704 trconfiguration = gAlice->GetTriggerDescriptor();
707 AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
710 runLoader->MakeTree( "GG" );
711 AliCentralTrigger* aCTP = runLoader->GetTrigger();
712 // Load Configuration
713 if (!aCTP->LoadConfiguration( trconfiguration ))
717 if( !aCTP->RunTrigger( runLoader , detectors ) ) {
729 //_____________________________________________________________________________
730 Bool_t AliSimulation::WriteTriggerRawData()
732 // Writes the CTP (trigger) DDL raw data
733 // Details of the format are given in the
734 // trigger TDR - pages 134 and 135.
735 AliCTPRawData writer;
741 //_____________________________________________________________________________
742 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
744 // run the generation and simulation
748 // initialize CDB storage and run number from external environment
749 // (either CDB manager or AliSimulation setters)
755 AliError("no gAlice object. Restart aliroot and try again.");
758 if (gAlice->Modules()->GetEntries() > 0) {
759 AliError("gAlice was already run. Restart aliroot and try again.");
763 AliInfo(Form("initializing gAlice with config file %s",
764 fConfigFileName.Data()));
765 StdoutToAliInfo(StderrToAliError(
766 gAlice->Init(fConfigFileName.Data());
769 // Get the trigger descriptor string
770 // Either from AliSimulation or from
772 if (fMakeTrigger.IsNull()) {
773 if (strcmp(gAlice->GetTriggerDescriptor(),""))
774 fMakeTrigger = gAlice->GetTriggerDescriptor();
777 gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
779 // Set run number in CDBManager
780 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
782 AliRunLoader* runLoader = AliRunLoader::GetRunLoader();
784 AliError(Form("gAlice has no run loader object. "
785 "Check your config file: %s", fConfigFileName.Data()));
788 SetGAliceFile(runLoader->GetFileName());
791 #if ROOT_VERSION_CODE < 331527
792 AliGeomManager::SetGeometry(gGeoManager);
794 // Check that the consistency of symbolic names for the activated subdetectors
795 // in the geometry loaded by AliGeomManager
796 TString detsToBeChecked = "";
797 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
798 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
799 AliModule* det = (AliModule*) detArray->At(iDet);
800 if (!det || !det->IsActive()) continue;
801 detsToBeChecked += det->GetName();
802 detsToBeChecked += " ";
803 } // end loop over detectors
804 if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
805 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
806 MisalignGeometry(runLoader);
809 // AliRunLoader* runLoader = AliRunLoader::GetRunLoader();
811 // AliError(Form("gAlice has no run loader object. "
812 // "Check your config file: %s", fConfigFileName.Data()));
815 // SetGAliceFile(runLoader->GetFileName());
817 if (!gAlice->Generator()) {
818 AliError(Form("gAlice has no generator object. "
819 "Check your config file: %s", fConfigFileName.Data()));
823 // Write GRP entry corresponding to the setting found in Cofig.C
827 if (nEvents <= 0) nEvents = fNEvents;
829 // get vertex from background file in case of merging
830 if (fUseBkgrdVertex &&
831 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
832 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
833 const char* fileName = ((TObjString*)
834 (fBkgrdFileNames->At(0)))->GetName();
835 AliInfo(Form("The vertex will be taken from the background "
836 "file %s with nSignalPerBackground = %d",
837 fileName, signalPerBkgrd));
838 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
839 gAlice->Generator()->SetVertexGenerator(vtxGen);
842 if (!fRunSimulation) {
843 gAlice->Generator()->SetTrackingFlag(0);
846 // set the number of events per file for given detectors and data types
847 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
848 if (!fEventsPerFile[i]) continue;
849 const char* detName = fEventsPerFile[i]->GetName();
850 const char* typeName = fEventsPerFile[i]->GetTitle();
851 TString loaderName(detName);
852 loaderName += "Loader";
853 AliLoader* loader = runLoader->GetLoader(loaderName);
855 AliError(Form("RunSimulation", "no loader for %s found\n"
856 "Number of events per file not set for %s %s",
857 detName, typeName, detName));
860 AliDataLoader* dataLoader =
861 loader->GetDataLoader(typeName);
863 AliError(Form("no data loader for %s found\n"
864 "Number of events per file not set for %s %s",
865 typeName, detName, typeName));
868 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
869 AliDebug(1, Form("number of events per file set to %d for %s %s",
870 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
873 AliInfo("running gAlice");
874 AliSysInfo::AddStamp("Start_simulation");
875 StdoutToAliInfo(StderrToAliError(
876 gAlice->RunMC(nEvents);
878 AliSysInfo::AddStamp("Stop_simulation");
884 //_____________________________________________________________________________
885 Bool_t AliSimulation::RunSDigitization(const char* detectors)
887 // run the digitization and produce summable digits
888 static Int_t eventNr=0;
891 // initialize CDB storage, run number, set CDB lock
893 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
896 AliRunLoader* runLoader = LoadRun();
897 if (!runLoader) return kFALSE;
899 TString detStr = detectors;
900 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
901 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
902 AliModule* det = (AliModule*) detArray->At(iDet);
903 if (!det || !det->IsActive()) continue;
904 if (IsSelected(det->GetName(), detStr)) {
905 AliInfo(Form("creating summable digits for %s", det->GetName()));
906 AliCodeTimerAuto(Form("creating summable digits for %s", det->GetName()));
908 AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
912 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
913 AliError(Form("the following detectors were not found: %s",
915 if (fStopOnError) return kFALSE;
924 //_____________________________________________________________________________
925 Bool_t AliSimulation::RunDigitization(const char* detectors,
926 const char* excludeDetectors)
928 // run the digitization and produce digits from sdigits
932 // initialize CDB storage, run number, set CDB lock
934 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
937 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
938 if (gAlice) delete gAlice;
942 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
943 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
944 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
945 // manager->SetEmbeddingFlag(fEmbeddingFlag);
946 manager->SetInputStream(0, fGAliceFileName.Data());
947 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
948 const char* fileName = ((TObjString*)
949 (fBkgrdFileNames->At(iStream-1)))->GetName();
950 manager->SetInputStream(iStream, fileName);
953 TString detStr = detectors;
954 TString detExcl = excludeDetectors;
955 manager->GetInputStream(0)->ImportgAlice();
956 AliRunLoader* runLoader =
957 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
958 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
959 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
960 AliModule* det = (AliModule*) detArray->At(iDet);
961 if (!det || !det->IsActive()) continue;
962 if (IsSelected(det->GetName(), detStr) &&
963 !IsSelected(det->GetName(), detExcl)) {
964 AliDigitizer* digitizer = det->CreateDigitizer(manager);
967 AliError(Form("no digitizer for %s", det->GetName()));
968 if (fStopOnError) return kFALSE;
970 digitizer->SetRegionOfInterest(fRegionOfInterest);
975 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
976 AliError(Form("the following detectors were not found: %s",
978 if (fStopOnError) return kFALSE;
981 if (!manager->GetListOfTasks()->IsEmpty()) {
982 AliInfo("executing digitization");
991 //_____________________________________________________________________________
992 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
994 // run the digitization and produce digits from hits
998 // initialize CDB storage, run number, set CDB lock
1000 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1003 AliRunLoader* runLoader = LoadRun("READ");
1004 if (!runLoader) return kFALSE;
1006 TString detStr = detectors;
1007 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1008 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1009 AliModule* det = (AliModule*) detArray->At(iDet);
1010 if (!det || !det->IsActive()) continue;
1011 if (IsSelected(det->GetName(), detStr)) {
1012 AliInfo(Form("creating digits from hits for %s", det->GetName()));
1017 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1018 AliError(Form("the following detectors were not found: %s",
1020 if (fStopOnError) return kFALSE;
1026 //_____________________________________________________________________________
1027 Bool_t AliSimulation::WriteRawData(const char* detectors,
1028 const char* fileName,
1029 Bool_t deleteIntermediateFiles,
1032 // convert the digits to raw data
1033 // First DDL raw data files for the given detectors are created.
1034 // If a file name is given, the DDL files are then converted to a DATE file.
1035 // If deleteIntermediateFiles is true, the DDL raw files are deleted
1037 // If the file name has the extension ".root", the DATE file is converted
1039 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1040 // 'selrawdata' flag can be used to enable writing of detectors raw data
1041 // accoring to the trigger cluster.
1043 AliCodeTimerAuto("")
1045 TString detStr = detectors;
1046 if (!WriteRawFiles(detStr.Data())) {
1047 if (fStopOnError) return kFALSE;
1050 // run HLT simulation on simulated DDL raw files
1051 // and produce HLT ddl raw files to be included in date/root file
1052 if (IsSelected("HLT", detStr) && !fRunHLT.IsNull()) {
1054 if (fStopOnError) return kFALSE;
1058 TString dateFileName(fileName);
1059 if (!dateFileName.IsNull()) {
1060 Bool_t rootOutput = dateFileName.EndsWith(".root");
1061 if (rootOutput) dateFileName += ".date";
1062 TString selDateFileName;
1064 selDateFileName = "selected.";
1065 selDateFileName+= dateFileName;
1067 if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1068 if (fStopOnError) return kFALSE;
1070 if (deleteIntermediateFiles) {
1071 AliRunLoader* runLoader = LoadRun("READ");
1072 if (runLoader) for (Int_t iEvent = 0;
1073 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1075 sprintf(command, "rm -r raw%d", iEvent);
1076 gSystem->Exec(command);
1081 if (!ConvertDateToRoot(dateFileName, fileName)) {
1082 if (fStopOnError) return kFALSE;
1084 if (deleteIntermediateFiles) {
1085 gSystem->Unlink(dateFileName);
1088 TString selFileName = "selected.";
1089 selFileName += fileName;
1090 if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1091 if (fStopOnError) return kFALSE;
1093 if (deleteIntermediateFiles) {
1094 gSystem->Unlink(selDateFileName);
1103 //_____________________________________________________________________________
1104 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1106 // convert the digits to raw data DDL files
1108 AliCodeTimerAuto("")
1110 AliRunLoader* runLoader = LoadRun("READ");
1111 if (!runLoader) return kFALSE;
1113 // write raw data to DDL files
1114 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1115 AliInfo(Form("processing event %d", iEvent));
1116 runLoader->GetEvent(iEvent);
1117 TString baseDir = gSystem->WorkingDirectory();
1119 sprintf(dirName, "raw%d", iEvent);
1120 gSystem->MakeDirectory(dirName);
1121 if (!gSystem->ChangeDirectory(dirName)) {
1122 AliError(Form("couldn't change to directory %s", dirName));
1123 if (fStopOnError) return kFALSE; else continue;
1126 ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1129 TString detStr = detectors;
1130 if (IsSelected("HLT", detStr)) {
1131 // Do nothing. "HLT" will be removed from detStr and HLT raw
1132 // data files are generated in RunHLT.
1135 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1136 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1137 AliModule* det = (AliModule*) detArray->At(iDet);
1138 if (!det || !det->IsActive()) continue;
1139 if (IsSelected(det->GetName(), detStr)) {
1140 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1145 if (!WriteTriggerRawData())
1146 if (fStopOnError) return kFALSE;
1148 gSystem->ChangeDirectory(baseDir);
1149 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1150 AliError(Form("the following detectors were not found: %s",
1152 if (fStopOnError) return kFALSE;
1161 //_____________________________________________________________________________
1162 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1163 const char* selDateFileName)
1165 // convert raw data DDL files to a DATE file with the program "dateStream"
1166 // The second argument is not empty when the user decides to write
1167 // the detectors raw data according to the trigger cluster.
1169 AliCodeTimerAuto("")
1171 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1173 AliError("the program dateStream was not found");
1174 if (fStopOnError) return kFALSE;
1179 AliRunLoader* runLoader = LoadRun("READ");
1180 if (!runLoader) return kFALSE;
1182 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1183 Bool_t selrawdata = kFALSE;
1184 if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1187 // Note the option -s. It is used in order to avoid
1188 // the generation of SOR/EOR events.
1189 sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d",
1190 dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1191 FILE* pipe = gSystem->OpenPipe(command, "w");
1194 AliError(Form("Cannot execute command: %s",command));
1198 Int_t selEvents = 0;
1199 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1201 UInt_t detectorPattern = 0;
1202 runLoader->GetEvent(iEvent);
1203 if (!runLoader->LoadTrigger()) {
1204 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1205 detectorPattern = aCTP->GetClusterMask();
1206 // Check if the event was triggered by CTP
1208 if (aCTP->GetClassMask()) selEvents++;
1212 AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1214 AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1215 selrawdata = kFALSE;
1219 fprintf(pipe, "GDC DetectorPattern %u\n", detectorPattern);
1223 // loop over detectors and DDLs
1224 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1225 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1227 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1228 Int_t ldcID = Int_t(ldc + 0.0001);
1229 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1231 char rawFileName[256];
1232 sprintf(rawFileName, "raw%d/%s",
1233 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1235 // check existence and size of raw data file
1236 FILE* file = fopen(rawFileName, "rb");
1237 if (!file) continue;
1238 fseek(file, 0, SEEK_END);
1239 unsigned long size = ftell(file);
1241 if (!size) continue;
1243 if (ldcID != prevLDC) {
1244 fprintf(pipe, " LDC Id %d\n", ldcID);
1247 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1252 Int_t result = gSystem->ClosePipe(pipe);
1254 if (!(selrawdata && selEvents > 0)) {
1256 return (result == 0);
1259 AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1261 sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d",
1262 selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1263 FILE* pipe2 = gSystem->OpenPipe(command, "w");
1265 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1267 // Get the trigger decision and cluster
1268 UInt_t detectorPattern = 0;
1270 runLoader->GetEvent(iEvent);
1271 if (!runLoader->LoadTrigger()) {
1272 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1273 if (aCTP->GetClassMask() == 0) continue;
1274 detectorPattern = aCTP->GetClusterMask();
1275 detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1276 AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1279 fprintf(pipe2, "GDC DetectorPattern %u\n", detectorPattern);
1283 // loop over detectors and DDLs
1284 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1285 // Write only raw data from detectors that
1286 // are contained in the trigger cluster(s)
1287 if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1289 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1291 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1292 Int_t ldcID = Int_t(ldc + 0.0001);
1293 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1295 char rawFileName[256];
1296 sprintf(rawFileName, "raw%d/%s",
1297 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1299 // check existence and size of raw data file
1300 FILE* file = fopen(rawFileName, "rb");
1301 if (!file) continue;
1302 fseek(file, 0, SEEK_END);
1303 unsigned long size = ftell(file);
1305 if (!size) continue;
1307 if (ldcID != prevLDC) {
1308 fprintf(pipe2, " LDC Id %d\n", ldcID);
1311 fprintf(pipe2, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1316 Int_t result2 = gSystem->ClosePipe(pipe2);
1319 return ((result == 0) && (result2 == 0));
1322 //_____________________________________________________________________________
1323 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1324 const char* rootFileName)
1326 // convert a DATE file to a root file with the program "alimdc"
1329 const Int_t kDBSize = 2000000000;
1330 const Int_t kTagDBSize = 1000000000;
1331 const Bool_t kFilter = kFALSE;
1332 const Int_t kCompression = 1;
1334 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1336 AliError("the program alimdc was not found");
1337 if (fStopOnError) return kFALSE;
1342 AliInfo(Form("converting DATE file %s to root file %s",
1343 dateFileName, rootFileName));
1345 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1346 const char* tagDBFS = "/tmp/mdc1/tags";
1348 // User defined file system locations
1349 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1350 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1351 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1352 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1353 if (gSystem->Getenv("ALIMDC_TAGDB"))
1354 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1356 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1357 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1358 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1360 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1361 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1362 gSystem->Exec(Form("mkdir %s",tagDBFS));
1364 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1365 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1366 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1368 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1369 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1370 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1372 return (result == 0);
1376 //_____________________________________________________________________________
1377 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1379 // delete existing run loaders, open a new one and load gAlice
1381 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1382 AliRunLoader* runLoader =
1383 AliRunLoader::Open(fGAliceFileName.Data(),
1384 AliConfig::GetDefaultEventFolderName(), mode);
1386 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1389 runLoader->LoadgAlice();
1390 runLoader->LoadHeader();
1391 gAlice = runLoader->GetAliRun();
1393 AliError(Form("no gAlice object found in file %s",
1394 fGAliceFileName.Data()));
1400 //_____________________________________________________________________________
1401 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1403 // get or calculate the number of signal events per background event
1405 if (!fBkgrdFileNames) return 1;
1406 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1407 if (nBkgrdFiles == 0) return 1;
1409 // get the number of signal events
1411 AliRunLoader* runLoader =
1412 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1413 if (!runLoader) return 1;
1415 nEvents = runLoader->GetNumberOfEvents();
1420 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1421 // get the number of background events
1422 const char* fileName = ((TObjString*)
1423 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1424 AliRunLoader* runLoader =
1425 AliRunLoader::Open(fileName, "BKGRD");
1426 if (!runLoader) continue;
1427 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1430 // get or calculate the number of signal per background events
1431 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1432 if (nSignalPerBkgrd <= 0) {
1433 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1434 } else if (result && (result != nSignalPerBkgrd)) {
1435 AliInfo(Form("the number of signal events per background event "
1436 "will be changed from %d to %d for stream %d",
1437 nSignalPerBkgrd, result, iBkgrdFile+1));
1438 nSignalPerBkgrd = result;
1441 if (!result) result = nSignalPerBkgrd;
1442 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1443 AliWarning(Form("not enough background events (%d) for %d signal events "
1444 "using %d signal per background events for stream %d",
1445 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1452 //_____________________________________________________________________________
1453 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1455 // check whether detName is contained in detectors
1456 // if yes, it is removed from detectors
1458 // check if all detectors are selected
1459 if ((detectors.CompareTo("ALL") == 0) ||
1460 detectors.BeginsWith("ALL ") ||
1461 detectors.EndsWith(" ALL") ||
1462 detectors.Contains(" ALL ")) {
1467 // search for the given detector
1468 Bool_t result = kFALSE;
1469 if ((detectors.CompareTo(detName) == 0) ||
1470 detectors.BeginsWith(detName+" ") ||
1471 detectors.EndsWith(" "+detName) ||
1472 detectors.Contains(" "+detName+" ")) {
1473 detectors.ReplaceAll(detName, "");
1477 // clean up the detectors string
1478 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1479 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1480 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1485 //_____________________________________________________________________________
1486 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName)
1489 // Steering routine to convert raw data in directory rawDirectory/ to fake SDigits.
1490 // These can be used for embedding of MC tracks into RAW data using the standard
1491 // merging procedure.
1493 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1496 AliError("no gAlice object. Restart aliroot and try again.");
1499 if (gAlice->Modules()->GetEntries() > 0) {
1500 AliError("gAlice was already run. Restart aliroot and try again.");
1504 AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1505 StdoutToAliInfo(StderrToAliError(gAlice->Init(fConfigFileName.Data());););
1509 //AliCDBManager* man = AliCDBManager::Instance();
1510 //man->SetRun(0); // Should this come from rawdata header ?
1514 // Get the runloader
1515 AliRunLoader* runLoader = AliRunLoader::GetRunLoader();
1517 // Open esd file if available
1518 TFile* esdFile = TFile::Open(esdFileName);
1519 Bool_t esdOK = (esdFile != 0);
1520 AliESD* esd = new AliESD;
1523 treeESD = (TTree*) esdFile->Get("esdTree");
1525 AliWarning("No ESD tree found");
1528 treeESD->SetBranchAddress("ESD", &esd);
1532 // Create the RawReader
1533 TString fileName(rawDirectory);
1534 AliRawReader* rawReader = 0x0;
1535 if (fileName.EndsWith("/")) {
1536 rawReader = new AliRawReaderFile(fileName);
1537 } else if (fileName.EndsWith(".root")) {
1538 rawReader = new AliRawReaderRoot(fileName);
1539 } else if (!fileName.IsNull()) {
1540 rawReader = new AliRawReaderDate(fileName);
1542 // if (!fEquipIdMap.IsNull() && fRawReader)
1543 // fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1545 // Get list of detectors
1546 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1549 AliHeader* header = runLoader->GetHeader();
1551 TString detStr = fMakeSDigits;
1555 if (!(rawReader->NextEvent())) break;
1558 for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1559 AliModule* det = (AliModule*) detArray->At(iDet);
1560 if (!det || !det->IsActive()) continue;
1561 if (IsSelected(det->GetName(), detStr)) {
1562 AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1563 det->Raw2SDigits(rawReader);
1570 // If ESD information available obtain reconstructed vertex and store in header.
1572 treeESD->GetEvent(nev);
1573 const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1574 Double_t position[3];
1575 esdVertex->GetXYZ(position);
1576 AliGenEventHeader* mcHeader = new AliGenEventHeader("ESD");
1579 for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1580 mcHeader->SetPrimaryVertex(mcV);
1581 header->Reset(0,nev);
1582 header->SetGenEventHeader(mcHeader);
1583 printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1588 runLoader->TreeE()->Fill();
1589 runLoader->SetNextEvent();
1595 runLoader->CdGAFile();
1596 runLoader->WriteHeader("OVERWRITE");
1597 runLoader->WriteRunLoader();
1602 //_____________________________________________________________________________
1603 Int_t AliSimulation::GetDetIndex(const char* detector)
1605 // return the detector index corresponding to detector
1607 for (index = 0; index < fgkNDetectors ; index++) {
1608 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1614 //_____________________________________________________________________________
1615 Bool_t AliSimulation::RunHLT()
1617 // Run the HLT simulation
1618 // HLT simulation is implemented in HLT/sim/AliHLTSimulation
1619 // Disabled if fRunHLT is empty, default vaule is "default".
1620 // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
1621 // The default simulation depends on the HLT component libraries and their
1622 // corresponding agents which define components and chains to run. See
1623 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
1624 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
1626 // The libraries to be loaded can be specified as an option.
1628 // AliSimulation sim;
1629 // sim.SetRunHLT("libAliHLTSample.so");
1631 // will only load <tt>libAliHLTSample.so</tt>
1633 // Other available options:
1634 // \li loglevel=<i>level</i> <br>
1635 // logging level for this processing
1637 // disable redirection of log messages to AliLog class
1638 // \li config=<i>macro</i>
1639 // configuration macro
1640 // \li chains=<i>configuration</i>
1641 // comma separated list of configurations to be run during simulation
1642 // \li rawfile=<i>file</i>
1643 // source for the RawReader to be created, the default is <i>./</i> if
1644 // raw data is simulated
1647 AliRunLoader* pRunLoader = LoadRun("READ");
1648 if (!pRunLoader) return kFALSE;
1650 // initialize CDB storage, run number, set CDB lock
1652 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1655 // load the library dynamically
1656 gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1658 // check for the library version
1659 AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
1661 AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
1664 if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
1665 AliError(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
1669 // print compile info
1670 typedef void (*CompileInfo)( char*& date, char*& time);
1671 CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
1675 (*fctInfo)(date,time);
1676 if (!date) {date=new Char_t[8]; strcpy(date,"unknown");}
1677 if (!time) {time=new Char_t[8]; strcpy(time,"unknown");}
1678 AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
1682 AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
1685 // create instance of the HLT simulation
1686 AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
1687 AliHLTSimulation* pHLT=NULL;
1688 if (fctCreate==NULL || (pHLT=(fctCreate()))==NULL) {
1689 AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
1693 // init the HLT simulation
1695 if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
1696 TString detStr = fWriteRawData;
1697 if (!IsSelected("HLT", detStr)) {
1698 options+=" writerawfiles=";
1700 options+=" writerawfiles=HLT";
1703 if (!detStr.IsNull() && !options.Contains("rawfile=")) {
1704 // as a matter of fact, HLT will run reconstruction and needs the RawReader
1705 // in order to get detector data. By default, RawReaderFile is used to read
1706 // the already simulated ddl files. Date and Root files from the raw data
1707 // are generated after the HLT simulation.
1708 options+=" rawfile=./";
1711 AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
1712 if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
1713 AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
1715 // run the HLT simulation
1716 AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
1717 if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
1718 AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
1722 // delete the instance
1723 AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
1724 if (fctDelete==NULL || fctDelete(pHLT)<0) {
1725 AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
1729 return iResult>=0?kTRUE:kFALSE;
1732 //_____________________________________________________________________________
1733 Bool_t AliSimulation::RunQA()
1735 // run the QA on summable hits, digits or digits
1737 if(!gAlice) return kFALSE;
1738 fQASteer->SetRunLoader(AliRunLoader::GetRunLoader()) ;
1740 TString detectorsw("") ;
1742 detectorsw = fQASteer->Run(fQADetectors.Data()) ;
1743 if ( detectorsw.IsNull() )
1746 fQASteer->EndOfCycle(detectorsw) ;
1750 //_____________________________________________________________________________
1751 Bool_t AliSimulation::SetRunQA(TString detAndAction)
1753 // Allows to run QA for a selected set of detectors
1754 // and a selected set of tasks among HITS, SDIGITS and DIGITS
1755 // all selected detectors run the same selected tasks
1757 if (!detAndAction.Contains(":")) {
1758 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
1762 Int_t colon = detAndAction.Index(":") ;
1763 fQADetectors = detAndAction(0, colon) ;
1764 if (fQADetectors.Contains("ALL") )
1765 fQADetectors = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
1766 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
1767 if (fQATasks.Contains("ALL") ) {
1768 fQATasks = Form("%d %d %d", AliQA::kHITS, AliQA::kSDIGITS, AliQA::kDIGITS) ;
1770 fQATasks.ToUpper() ;
1772 if ( fQATasks.Contains("HIT") )
1773 tempo = Form("%d ", AliQA::kHITS) ;
1774 if ( fQATasks.Contains("SDIGIT") )
1775 tempo += Form("%d ", AliQA::kSDIGITS) ;
1776 if ( fQATasks.Contains("DIGIT") )
1777 tempo += Form("%d ", AliQA::kDIGITS) ;
1779 if (fQATasks.IsNull()) {
1780 AliInfo("No QA requested\n") ;
1785 TString tempo(fQATasks) ;
1786 tempo.ReplaceAll(Form("%d", AliQA::kHITS), AliQA::GetTaskName(AliQA::kHITS)) ;
1787 tempo.ReplaceAll(Form("%d", AliQA::kSDIGITS), AliQA::GetTaskName(AliQA::kSDIGITS)) ;
1788 tempo.ReplaceAll(Form("%d", AliQA::kDIGITS), AliQA::GetTaskName(AliQA::kDIGITS)) ;
1789 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
1791 fQASteer->SetActiveDetectors(fQADetectors) ;
1792 fQASteer->SetTasks(fQATasks) ;
1796 //_____________________________________________________________________________
1797 void AliSimulation::ProcessEnvironmentVars()
1799 // Extract run number and random generator seed from env variables
1801 AliInfo("Processing environment variables");
1803 // Random Number seed
1805 // first check that seed is not already set
1807 if (gSystem->Getenv("CONFIG_SEED")) {
1808 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
1811 if (gSystem->Getenv("CONFIG_SEED")) {
1812 AliInfo(Form("Seed for random number generation already set (%d)"
1813 ": CONFIG_SEED variable ignored!", fSeed));
1817 AliInfo(Form("Seed for random number generation = %d ", fSeed));
1821 // first check that run number is not already set
1823 if (gSystem->Getenv("DC_RUN")) {
1824 fRun = atoi(gSystem->Getenv("DC_RUN"));
1827 if (gSystem->Getenv("DC_RUN")) {
1828 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
1832 AliInfo(Form("Run number = %d", fRun));
1835 //---------------------------------------------------------------------
1837 void AliSimulation::WriteGRPEntry()
1839 // Get the necessary information from galice (generator, trigger etc) and
1840 // write a GRP entry corresponding to the settings in the Config.C used
1841 // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
1844 AliInfo("Writing global run parameters entry into the OCDB");
1846 AliGRPObject* grpObj = new AliGRPObject();
1848 grpObj->SetRunType("PHYSICS");
1849 grpObj->SetTimeStart(0);
1850 grpObj->SetTimeEnd(9999);
1852 const AliGenerator *gen = gAlice->Generator();
1854 grpObj->SetBeamEnergy(gen->GetEnergyCMS());
1857 gen->GetProjectile(projectile,a,z);
1859 gen->GetTarget(target,a,z);
1860 TString beamType = projectile + "-" + target;
1861 beamType.ReplaceAll(" ","");
1862 if (!beamType.CompareTo("-")) {
1864 grpObj->SetBeamType("UNKNOWN");
1867 grpObj->SetBeamType(beamType);
1871 AliWarning("Unknown beam type and energy! Setting energy to 0");
1872 grpObj->SetBeamEnergy(0);
1873 grpObj->SetBeamType("UNKNOWN");
1876 UInt_t detectorPattern = 0;
1878 TObjArray *detArray = gAlice->Detectors();
1879 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
1880 if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
1881 detectorPattern |= (1 << iDet);
1886 if (!fMakeTrigger.IsNull() || strcmp(gAlice->GetTriggerDescriptor(),""))
1887 detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
1890 if (!fRunHLT.IsNull())
1891 detectorPattern |= (1 << AliDAQ::kHLTId);
1893 grpObj->SetNumberOfDetectors((Char_t)nDets);
1894 grpObj->SetDetectorMask((Int_t)detectorPattern);
1895 grpObj->SetLHCPeriod("LHC08c");
1896 grpObj->SetLHCState("STABLE_BEAMS");
1897 grpObj->SetLHCLuminosity(0,(AliGRPObject::Stats)0);
1898 grpObj->SetBeamIntensity(0,(AliGRPObject::Stats)0);
1900 AliMagF *field = gAlice->Field();
1901 Float_t solenoidField = TMath::Abs(field->SolenoidField());
1902 Float_t factor = field->Factor();
1903 Float_t l3current = TMath::Abs(factor)*solenoidField*30000./5.;
1904 grpObj->SetL3Current(l3current,(AliGRPObject::Stats)0);
1907 grpObj->SetL3Polarity(0);
1908 grpObj->SetDipolePolarity(0);
1911 grpObj->SetL3Polarity(1);
1912 grpObj->SetDipolePolarity(1);
1915 if (TMath::Abs(factor) != 0)
1916 grpObj->SetDipoleCurrent(6000,(AliGRPObject::Stats)0);
1918 grpObj->SetDipoleCurrent(0,(AliGRPObject::Stats)0);
1920 grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
1922 //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
1924 // Now store the entry in OCDB
1925 AliCDBManager* man = AliCDBManager::Instance();
1927 AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
1928 AliCDBMetaData *metadata= new AliCDBMetaData();
1930 metadata->SetResponsible("alice-off@cern.ch");
1931 metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
1933 man->Put(grpObj,id,metadata);