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 **************************************************************************/
16 /* $Id: AliSimulation.cxx 63204 2013-06-26 13:33:28Z rgrosso $ */
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 ///////////////////////////////////////////////////////////////////////////////
109 #include <TGeoGlobalMagField.h>
110 #include <TGeoManager.h>
111 #include <TObjString.h>
114 #include <TVirtualMC.h>
115 #include <TVirtualMCApplication.h>
117 #include <TInterpreter.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 "AliDigitizationInput.h"
149 #include "AliRunLoader.h"
150 #include "AliStack.h"
151 #include "AliSimulation.h"
152 #include "AliSysInfo.h"
153 #include "AliVertexGenFile.h"
156 ClassImp(AliSimulation)
158 AliSimulation *AliSimulation::fgInstance = 0;
159 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE","AD","MFT","HLT"};
161 //_____________________________________________________________________________
162 AliSimulation::AliSimulation(const char* configFileName,
163 const char* name, const char* title) :
166 fRunGeneratorOnly(kFALSE),
167 fRunGeneration(kTRUE),
168 fRunSimulation(kTRUE),
169 fLoadAlignFromCDB(kTRUE),
170 fLoadAlObjsListOfDets("ALL"),
174 fMakeDigitsFromHits(""),
176 fRawDataFileName(""),
177 fDeleteIntermediateFiles(kFALSE),
178 fWriteSelRawData(kFALSE),
179 fStopOnError(kFALSE),
180 fUseMonitoring(kFALSE),
182 fConfigFileName(configFileName),
183 fGAliceFileName("galice.root"),
185 fBkgrdFileNames(NULL),
186 fAlignObjArray(NULL),
187 fUseBkgrdVertex(kTRUE),
188 fRegionOfInterest(kFALSE),
194 fInitCDBCalled(kFALSE),
195 fInitRunNumberCalled(kFALSE),
196 fSetRunNumberFromDataCalled(kFALSE),
197 fEmbeddingFlag(kFALSE),
200 fUseVertexFromCDB(0),
201 fUseMagFieldFromGRP(0),
202 fGRPWriteLocation(Form("local://%s", gSystem->pwd())),
203 fUseTimeStampFromCDB(0),
209 fEventSpecie(AliRecoParam::kDefault),
210 fWriteQAExpertData(kTRUE),
214 fWriteGRPEntry(kTRUE)
216 // create simulation object with default parameters
218 SetGAliceFile("galice.root");
221 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
222 qam->SetActiveDetectors(fQADetectors) ;
223 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
224 qam->SetTasks(fQATasks) ;
227 //_____________________________________________________________________________
228 AliSimulation::~AliSimulation()
232 fEventsPerFile.Delete();
233 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
234 // delete fAlignObjArray; fAlignObjArray=0;
236 if (fBkgrdFileNames) {
237 fBkgrdFileNames->Delete();
238 delete fBkgrdFileNames;
241 fSpecCDBUri.Delete();
242 if (fgInstance==this) fgInstance = 0;
244 AliQAManager::QAManager()->ShowQA() ;
245 AliQAManager::Destroy() ;
246 AliCodeTimer::Instance()->Print();
250 //_____________________________________________________________________________
251 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
253 // set the number of events for one run
258 //_____________________________________________________________________________
259 void AliSimulation::InitQA()
261 // activate a default CDB storage
262 // First check if we have any CDB storage set, because it is used
263 // to retrieve the calibration and alignment constants
265 if (fInitCDBCalled) return;
266 fInitCDBCalled = kTRUE;
268 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
269 qam->SetActiveDetectors(fQADetectors) ;
270 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
271 qam->SetTasks(fQATasks) ;
272 if (fWriteQAExpertData)
273 qam->SetWriteExpert() ;
275 if (qam->IsDefaultStorageSet()) {
276 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
277 AliWarning("Default QA reference storage has been already set !");
278 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fQARefUri.Data()));
279 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
280 fQARefUri = qam->GetDefaultStorage()->GetURI();
282 if (fQARefUri.Length() > 0) {
283 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
284 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
285 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
287 fQARefUri="local://$ALICE_ROOT/QARef";
288 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
289 AliWarning("Default QA reference storage not yet set !!!!");
290 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
291 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
293 qam->SetDefaultStorage(fQARefUri);
297 //_____________________________________________________________________________
298 void AliSimulation::InitCDB()
300 // activate a default CDB storage
301 // First check if we have any CDB storage set, because it is used
302 // to retrieve the calibration and alignment constants
304 if (fInitCDBCalled) return;
305 fInitCDBCalled = kTRUE;
307 AliCDBManager* man = AliCDBManager::Instance();
308 if (man->IsDefaultStorageSet())
310 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
311 AliWarning("Default CDB storage has been already set !");
312 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
313 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
314 fCDBUri = man->GetDefaultStorage()->GetURI();
317 if (fCDBUri.Length() > 0)
319 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
320 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
321 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
323 fCDBUri="local://$ALICE_ROOT/OCDB";
324 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
325 AliWarning("Default CDB storage not yet set !!!!");
326 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
327 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
330 man->SetDefaultStorage(fCDBUri);
333 // Now activate the detector specific CDB storage locations
334 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
335 TObject* obj = fSpecCDBUri[i];
337 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
338 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
339 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
340 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
345 //_____________________________________________________________________________
346 void AliSimulation::InitRunNumber(){
347 // check run number. If not set, set it to 0 !!!!
349 if (fInitRunNumberCalled) return;
350 fInitRunNumberCalled = kTRUE;
352 AliCDBManager* man = AliCDBManager::Instance();
353 if (man->GetRun() >= 0)
355 AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
356 "Use external variable DC_RUN or AliSimulation::SetRun()!"));
360 AliDebug(2, Form("Setting CDB run number to: %d",fRun));
363 AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
372 //_____________________________________________________________________________
373 void AliSimulation::SetCDBLock() {
374 // Set CDB lock: from now on it is forbidden to reset the run number
375 // or the default storage or to activate any further storage!
377 ULong64_t key = AliCDBManager::Instance()->SetLock(1);
381 //_____________________________________________________________________________
382 void AliSimulation::SetDefaultStorage(const char* uri) {
383 // Store the desired default CDB storage location
384 // Activate it later within the Run() method
390 //_____________________________________________________________________________
391 void AliSimulation::SetQARefDefaultStorage(const char* uri) {
392 // Store the desired default CDB storage location
393 // Activate it later within the Run() method
396 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
399 //_____________________________________________________________________________
400 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
401 // Store a detector-specific CDB storage location
402 // Activate it later within the Run() method
404 AliCDBPath aPath(calibType);
405 if(!aPath.IsValid()){
406 AliError(Form("Not a valid path: %s", calibType));
410 TObject* obj = fSpecCDBUri.FindObject(calibType);
411 if (obj) fSpecCDBUri.Remove(obj);
412 fSpecCDBUri.Add(new TNamed(calibType, uri));
416 //_____________________________________________________________________________
417 void AliSimulation::SetRunNumber(Int_t run)
420 // Activate it later within the Run() method
425 //_____________________________________________________________________________
426 void AliSimulation::SetSeed(Int_t seed)
429 // Activate it later within the Run() method
434 //_____________________________________________________________________________
435 Bool_t AliSimulation::SetRunNumberFromData()
437 // Set the CDB manager run number
438 // The run number is retrieved from gAlice
440 if (fSetRunNumberFromDataCalled) return kTRUE;
441 fSetRunNumberFromDataCalled = kTRUE;
443 AliCDBManager* man = AliCDBManager::Instance();
444 Int_t runData = -1, runCDB = -1;
446 AliRunLoader* runLoader = LoadRun("READ");
447 if (!runLoader) return kFALSE;
449 runData = runLoader->GetHeader()->GetRun();
453 runCDB = man->GetRun();
455 if (runCDB != runData) {
456 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
457 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
458 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
459 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
464 man->SetRun(runData);
467 if(man->GetRun() < 0) {
468 AliError("Run number not properly initalized!");
477 //_____________________________________________________________________________
478 void AliSimulation::SetConfigFile(const char* fileName)
480 // set the name of the config file
482 fConfigFileName = fileName;
485 //_____________________________________________________________________________
486 void AliSimulation::SetGAliceFile(const char* fileName)
488 // set the name of the galice file
489 // the path is converted to an absolute one if it is relative
491 fGAliceFileName = fileName;
492 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
493 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
495 fGAliceFileName = absFileName;
496 delete[] absFileName;
499 AliDebug(2, Form("galice file name set to %s", fileName));
502 //_____________________________________________________________________________
503 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
506 // set the number of events per file for the given detector and data type
507 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
509 TNamed* obj = new TNamed(detector, type);
510 obj->SetUniqueID(nEvents);
511 fEventsPerFile.Add(obj);
514 //_____________________________________________________________________________
515 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
517 // Read the alignment objects from CDB.
518 // Each detector is supposed to have the
519 // alignment objects in DET/Align/Data CDB path.
520 // All the detector objects are then collected,
521 // sorted by geometry level (starting from ALIC) and
522 // then applied to the TGeo geometry.
523 // Finally an overlaps check is performed.
525 if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
526 AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
530 // initialize CDB storage, run number, set CDB lock
532 // if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
535 Bool_t delRunLoader = kFALSE;
537 runLoader = LoadRun("READ");
538 if (!runLoader) return kFALSE;
539 delRunLoader = kTRUE;
542 // Export ideal geometry
543 if(!IsGeometryFromFile()) AliGeomManager::GetGeometry()->Export("geometry.root");
545 // Load alignment data from CDB and apply to geometry through AliGeomManager
546 if(fLoadAlignFromCDB){
548 TString detStr = fLoadAlObjsListOfDets;
549 TString loadAlObjsListOfDets = "";
551 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
552 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
553 AliModule* det = (AliModule*) detArray->At(iDet);
554 if (!det || !det->IsActive()) continue;
555 if (IsSelected(det->GetName(), detStr)) {
556 //add det to list of dets to be aligned from CDB
557 loadAlObjsListOfDets += det->GetName();
558 loadAlObjsListOfDets += " ";
560 } // end loop over detectors
561 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
562 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
564 // Check if the array with alignment objects was
565 // provided by the user. If yes, apply the objects
566 // to the present TGeo geometry
567 if (fAlignObjArray) {
568 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
569 AliError("The misalignment of one or more volumes failed!"
570 "Compare the list of simulated detectors and the list of detector alignment data!");
571 if (delRunLoader) delete runLoader;
577 // Update the internal geometry of modules (ITS needs it)
578 TString detStr = fLoadAlObjsListOfDets;
579 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
580 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
582 AliModule* det = (AliModule*) detArray->At(iDet);
583 if (!det || !det->IsActive()) continue;
584 if (IsSelected(det->GetName(), detStr)) {
585 det->UpdateInternalGeometry();
587 } // end loop over detectors
590 if (delRunLoader) delete runLoader;
595 //_____________________________________________________________________________
596 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
598 // add a file with background events for merging
600 TObjString* fileNameStr = new TObjString(fileName);
601 fileNameStr->SetUniqueID(nSignalPerBkgrd);
602 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
603 fBkgrdFileNames->Add(fileNameStr);
606 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
608 // add a file with background events for embeddin
609 MergeWith(fileName, nSignalPerBkgrd);
610 fEmbeddingFlag = kTRUE;
613 //_____________________________________________________________________________
614 Bool_t AliSimulation::Run(Int_t nEvents)
616 // run the generation, simulation and digitization
619 AliCodeTimerAuto("",0)
620 AliSysInfo::AddStamp("Start_Run");
622 // Load run number and seed from environmental vars
623 ProcessEnvironmentVars();
624 AliSysInfo::AddStamp("ProcessEnvironmentVars");
626 gRandom->SetSeed(fSeed);
628 if (nEvents > 0) fNEvents = nEvents;
630 // Run generator-only code on demand
631 if (fRunGeneratorOnly)
633 if(!RunGeneratorOnly())
635 if (fStopOnError) return kFALSE;
641 // create and setup the HLT instance
642 if (!fRunHLT.IsNull() && !CreateHLT()) {
643 if (fStopOnError) return kFALSE;
648 // generation and simulation -> hits
649 if (fRunGeneration) {
650 if (!RunSimulation()) if (fStopOnError) return kFALSE;
652 AliSysInfo::AddStamp("RunSimulation");
654 // initialize CDB storage from external environment
655 // (either CDB manager or AliSimulation setters),
656 // if not already done in RunSimulation()
658 AliSysInfo::AddStamp("InitCDB");
660 // Set run number in CDBManager from data
661 // From this point on the run number must be always loaded from data!
662 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
664 // Set CDB lock: from now on it is forbidden to reset the run number
665 // or the default storage or to activate any further storage!
668 // If RunSimulation was not called, load the geometry and misalign it
669 if (!AliGeomManager::GetGeometry()) {
670 // Initialize the geometry manager
671 AliGeomManager::LoadGeometry("geometry.root");
672 AliSysInfo::AddStamp("GetGeometry");
673 // // Check that the consistency of symbolic names for the activated subdetectors
674 // // in the geometry loaded by AliGeomManager
675 // AliRunLoader* runLoader = LoadRun("READ");
676 // if (!runLoader) return kFALSE;
678 // TString detsToBeChecked = "";
679 // TObjArray* detArray = runLoader->GetAliRun()->Detectors();
680 // for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
681 // AliModule* det = (AliModule*) detArray->At(iDet);
682 // if (!det || !det->IsActive()) continue;
683 // detsToBeChecked += det->GetName();
684 // detsToBeChecked += " ";
685 // } // end loop over detectors
686 // if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
687 if(!AliGeomManager::CheckSymNamesLUT("ALL"))
688 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
690 if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
692 if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
694 AliSysInfo::AddStamp("MissalignGeometry");
697 // hits -> summable digits
698 AliSysInfo::AddStamp("Start_sdigitization");
699 if (!fMakeSDigits.IsNull()) {
700 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
703 AliSysInfo::AddStamp("Stop_sdigitization");
705 AliSysInfo::AddStamp("Start_digitization");
706 // summable digits -> digits
707 if (!fMakeDigits.IsNull()) {
708 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
709 if (fStopOnError) return kFALSE;
712 AliSysInfo::AddStamp("Stop_digitization");
717 if (!fMakeDigitsFromHits.IsNull()) {
718 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
719 AliWarning(Form("Merging and direct creation of digits from hits "
720 "was selected for some detectors. "
721 "No merging will be done for the following detectors: %s",
722 fMakeDigitsFromHits.Data()));
724 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
725 if (fStopOnError) return kFALSE;
729 AliSysInfo::AddStamp("Hits2Digits");
733 if (!fTriggerConfig.IsNull() && !RunTrigger(fTriggerConfig,fMakeDigits)) {
734 if (fStopOnError) return kFALSE;
737 AliSysInfo::AddStamp("RunTrigger");
740 // digits -> raw data
741 if (!fWriteRawData.IsNull()) {
742 if (!WriteRawData(fWriteRawData, fRawDataFileName,
743 fDeleteIntermediateFiles,fWriteSelRawData)) {
744 if (fStopOnError) return kFALSE;
748 AliSysInfo::AddStamp("WriteRaw");
750 // run HLT simulation on simulated digit data if raw data is not
751 // simulated, otherwise its called as part of WriteRawData
752 if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
754 if (fStopOnError) return kFALSE;
758 AliSysInfo::AddStamp("RunHLT");
762 Bool_t rv = RunQA() ;
768 AliSysInfo::AddStamp("RunQA");
770 TString snapshotFileOut("");
771 if(TString(gSystem->Getenv("OCDB_SNAPSHOT_CREATE")) == TString("kTRUE")){
772 AliInfo(" ******** Creating the snapshot! *********");
773 TString snapshotFile(gSystem->Getenv("OCDB_SNAPSHOT_FILENAME"));
774 if(!(snapshotFile.IsNull() || snapshotFile.IsWhitespace()))
775 snapshotFileOut = snapshotFile;
777 snapshotFileOut="OCDB.root";
778 AliCDBManager::Instance()->DumpToSnapshotFile(snapshotFileOut.Data(),kFALSE);
781 // Cleanup of CDB manager: cache and active storages!
782 AliCDBManager::Instance()->ClearCache();
787 //_______________________________________________________________________
788 Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
789 Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
790 Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
793 // Generates lego plots of:
794 // - radiation length map phi vs theta
795 // - radiation length map phi vs eta
796 // - interaction length map
797 // - g/cm2 length map
799 // ntheta bins in theta, eta
800 // themin minimum angle in theta (degrees)
801 // themax maximum angle in theta (degrees)
803 // phimin minimum angle in phi (degrees)
804 // phimax maximum angle in phi (degrees)
805 // rmin minimum radius
806 // rmax maximum radius
809 // The number of events generated = ntheta*nphi
810 // run input parameters in macro setup (default="Config.C")
812 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
815 <img src="picts/AliRunLego1.gif">
820 <img src="picts/AliRunLego2.gif">
825 <img src="picts/AliRunLego3.gif">
830 // run the generation and simulation
832 AliCodeTimerAuto("",0)
834 // initialize CDB storage and run number from external environment
835 // (either CDB manager or AliSimulation setters)
841 AliError("no gAlice object. Restart aliroot and try again.");
844 if (gAlice->Modules()->GetEntries() > 0) {
845 AliError("gAlice was already run. Restart aliroot and try again.");
849 AliInfo(Form("initializing gAlice with config file %s",
850 fConfigFileName.Data()));
853 if (nev == -1) nev = nc1 * nc2;
855 // check if initialisation has been done
856 // If runloader has been initialized, set the number of events per file to nc1 * nc2
859 if (!gener) gener = new AliLegoGenerator();
861 // Configure Generator
863 gener->SetRadiusRange(rmin, rmax);
864 gener->SetZMax(zmax);
865 gener->SetCoor1Range(nc1, c1min, c1max);
866 gener->SetCoor2Range(nc2, c2min, c2max);
870 fLego = new AliLego("lego",gener);
872 //__________________________________________________________________________
876 gROOT->LoadMacro(setup);
877 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
879 if(AliCDBManager::Instance()->GetRun() >= 0) {
880 SetRunNumber(AliCDBManager::Instance()->GetRun());
882 AliWarning("Run number not initialized!!");
885 AliRunLoader::Instance()->CdGAFile();
887 AliPDG::AddParticlesToPdgDataBase();
889 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
891 gAlice->GetMCApp()->Init();
894 //Must be here because some MCs (G4) adds detectors here and not in Config.C
895 gAlice->InitLoaders();
896 AliRunLoader::Instance()->MakeTree("E");
899 // Save stuff at the beginning of the file to avoid file corruption
900 AliRunLoader::Instance()->CdGAFile();
903 //Save current generator
904 AliGenerator *gen=gAlice->GetMCApp()->Generator();
905 gAlice->GetMCApp()->ResetGenerator(gener);
906 //Prepare MC for Lego Run
912 AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
913 gMC->ProcessRun(nev);
915 // End of this run, close files
917 // Restore current generator
918 gAlice->GetMCApp()->ResetGenerator(gen);
919 // Delete Lego Object
925 //_____________________________________________________________________________
926 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
930 AliCodeTimerAuto("",0)
932 // initialize CDB storage from external environment
933 // (either CDB manager or AliSimulation setters),
934 // if not already done in RunSimulation()
937 // Set run number in CDBManager from data
938 // From this point on the run number must be always loaded from data!
939 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
941 // Set CDB lock: from now on it is forbidden to reset the run number
942 // or the default storage or to activate any further storage!
945 AliRunLoader* runLoader = LoadRun("READ");
946 if (!runLoader) return kFALSE;
947 TString trconfiguration = config;
949 if (trconfiguration.IsNull()) {
950 if(!fTriggerConfig.IsNull()) {
951 trconfiguration = fTriggerConfig;
954 AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
957 runLoader->MakeTree( "GG" );
958 AliCentralTrigger* aCTP = runLoader->GetTrigger();
959 // Load Configuration
960 if (!aCTP->LoadConfiguration( trconfiguration ))
964 if( !aCTP->RunTrigger( runLoader , detectors ) ) {
976 //_____________________________________________________________________________
977 Bool_t AliSimulation::WriteTriggerRawData()
979 // Writes the CTP (trigger) DDL raw data
980 // Details of the format are given in the
981 // trigger TDR - pages 134 and 135.
982 AliCTPRawData writer;
988 //_____________________________________________________________________________
989 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
991 // run the generation and simulation
993 AliCodeTimerAuto("",0)
995 // initialize CDB storage and run number from external environment
996 // (either CDB manager or AliSimulation setters)
997 AliSysInfo::AddStamp("RunSimulation_Begin");
999 AliSysInfo::AddStamp("RunSimulation_InitCDB");
1003 AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
1006 AliError("no gAlice object. Restart aliroot and try again.");
1009 if (gAlice->Modules()->GetEntries() > 0) {
1010 AliError("gAlice was already run. Restart aliroot and try again.");
1014 // Setup monitoring if requested
1015 gAlice->GetMCApp()->SetUseMonitoring(fUseMonitoring);
1017 AliInfo(Form("initializing gAlice with config file %s",
1018 fConfigFileName.Data()));
1021 // Initialize ALICE Simulation run
1026 // If requested set the mag. field from the GRP entry.
1027 // After this the field is loccked and cannot be changed by Config.C
1028 if (fUseMagFieldFromGRP) {
1030 grpM.ReadGRPEntry();
1032 AliInfo("Field is locked now. It cannot be changed in Config.C");
1036 TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
1037 gROOT->LoadMacro(fConfigFileName.Data());
1038 Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
1039 if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
1040 AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
1042 AliSysInfo::AddStamp("RunSimulation_Config");
1045 // If requested obtain the vertex position and vertex sigma_z from the CDB
1046 // This overwrites the settings from the Config.C
1047 if (fUseVertexFromCDB) {
1048 Double_t vtxPos[3] = {0., 0., 0.};
1049 Double_t vtxSig[3] = {0., 0., 0.};
1050 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1052 AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1054 if(vertex->GetXRes()>2.8) { // > pipe radius --> it's a dummy object, don't use it
1055 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1056 if (entry) vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1060 vertex->GetXYZ(vtxPos);
1061 vertex->GetSigmaXYZ(vtxSig);
1062 AliInfo("Overwriting Config.C vertex settings !");
1063 AliInfo(Form("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
1064 vtxPos[0], vtxPos[1], vtxPos[2], vtxSig[2]));
1066 AliGenerator *gen = gAlice->GetMCApp()->Generator();
1067 gen->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]); // vertex position
1068 gen->SetSigmaZ(vtxSig[2]);
1073 // If requested we take the SOR and EOR time-stamps from the GRP and use them
1074 // in order to generate the event time-stamps
1075 if (fUseTimeStampFromCDB) {
1077 grpM.ReadGRPEntry();
1078 const AliGRPObject *grpObj = grpM.GetGRPData();
1079 if (!grpObj || (grpObj->GetTimeEnd() <= grpObj->GetTimeStart())) {
1080 AliError("Missing GRP or bad SOR/EOR time-stamps! Switching off the time-stamp generation from GRP!");
1081 fTimeStart = fTimeEnd = 0;
1082 fUseTimeStampFromCDB = kFALSE;
1085 fTimeStart = grpObj->GetTimeStart();
1086 fTimeEnd = grpObj->GetTimeEnd();
1090 if(AliCDBManager::Instance()->GetRun() >= 0) {
1091 AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
1092 AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
1094 AliWarning("Run number not initialized!!");
1097 AliRunLoader::Instance()->CdGAFile();
1100 AliPDG::AddParticlesToPdgDataBase();
1102 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1103 AliSysInfo::AddStamp("RunSimulation_GetField");
1104 gAlice->GetMCApp()->Init();
1105 AliSysInfo::AddStamp("RunSimulation_InitMCApp");
1107 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1108 gAlice->InitLoaders();
1109 AliRunLoader::Instance()->MakeTree("E");
1110 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1111 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1112 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1114 // Save stuff at the beginning of the file to avoid file corruption
1115 AliRunLoader::Instance()->CdGAFile();
1117 gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1118 AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1119 //___________________________________________________________________________________________
1121 AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
1123 // Set run number in CDBManager
1124 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
1126 AliRunLoader* runLoader = AliRunLoader::Instance();
1128 AliError(Form("gAlice has no run loader object. "
1129 "Check your config file: %s", fConfigFileName.Data()));
1132 SetGAliceFile(runLoader->GetFileName());
1134 // Misalign geometry
1135 #if ROOT_VERSION_CODE < 331527
1136 AliGeomManager::SetGeometry(gGeoManager);
1138 // Check that the consistency of symbolic names for the activated subdetectors
1139 // in the geometry loaded by AliGeomManager
1140 TString detsToBeChecked = "";
1141 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1142 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1143 AliModule* det = (AliModule*) detArray->At(iDet);
1144 if (!det || !det->IsActive()) continue;
1145 detsToBeChecked += det->GetName();
1146 detsToBeChecked += " ";
1147 } // end loop over detectors
1148 if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
1149 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
1150 MisalignGeometry(runLoader);
1151 AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
1154 // AliRunLoader* runLoader = AliRunLoader::Instance();
1155 // if (!runLoader) {
1156 // AliError(Form("gAlice has no run loader object. "
1157 // "Check your config file: %s", fConfigFileName.Data()));
1160 // SetGAliceFile(runLoader->GetFileName());
1162 if (!gAlice->GetMCApp()->Generator()) {
1163 AliError(Form("gAlice has no generator object. "
1164 "Check your config file: %s", fConfigFileName.Data()));
1168 // Write GRP entry corresponding to the setting found in Cofig.C
1171 AliSysInfo::AddStamp("RunSimulation_WriteGRP");
1173 if (nEvents <= 0) nEvents = fNEvents;
1175 // get vertex from background file in case of merging
1176 if (fUseBkgrdVertex &&
1177 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
1178 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
1179 const char* fileName = ((TObjString*)
1180 (fBkgrdFileNames->At(0)))->GetName();
1181 AliInfo(Form("The vertex will be taken from the background "
1182 "file %s with nSignalPerBackground = %d",
1183 fileName, signalPerBkgrd));
1184 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
1185 gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
1188 if (!fRunSimulation) {
1189 gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
1192 // set the number of events per file for given detectors and data types
1193 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
1194 if (!fEventsPerFile[i]) continue;
1195 const char* detName = fEventsPerFile[i]->GetName();
1196 const char* typeName = fEventsPerFile[i]->GetTitle();
1197 TString loaderName(detName);
1198 loaderName += "Loader";
1199 AliLoader* loader = runLoader->GetLoader(loaderName);
1201 AliError(Form("RunSimulation no loader for %s found\n Number of events per file not set for %s %s",
1202 detName, typeName, detName));
1205 AliDataLoader* dataLoader =
1206 loader->GetDataLoader(typeName);
1208 AliError(Form("no data loader for %s found\n"
1209 "Number of events per file not set for %s %s",
1210 typeName, detName, typeName));
1213 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
1214 AliDebug(1, Form("number of events per file set to %d for %s %s",
1215 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
1218 AliInfo("running gAlice");
1219 AliSysInfo::AddStamp("Start_ProcessRun");
1221 // Create the Root Tree with one branch per detector
1222 //Hits moved to begin event -> now we are crating separate tree for each event
1223 gMC->ProcessRun(nEvents);
1225 // End of this run, close files
1226 if(nEvents>0) FinishRun();
1228 AliSysInfo::AddStamp("Stop_ProcessRun");
1234 //_____________________________________________________________________________
1235 Bool_t AliSimulation::RunGeneratorOnly()
1238 TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
1239 gROOT->LoadMacro(fConfigFileName.Data());
1240 Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
1241 if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
1242 AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
1245 // Setup the runloader and generator, check if everything is OK
1246 AliRunLoader* runLoader = AliRunLoader::Instance();
1247 AliGenerator* generator = gAlice->GetMCApp()->Generator();
1249 AliError(Form("gAlice has no run loader object. "
1250 "Check your config file: %s", fConfigFileName.Data()));
1254 AliError(Form("gAlice has no generator object. "
1255 "Check your config file: %s", fConfigFileName.Data()));
1259 runLoader->LoadKinematics("RECREATE");
1260 runLoader->MakeTree("E");
1262 // Create stack and header
1263 runLoader->MakeStack();
1264 AliStack* stack = runLoader->Stack();
1265 AliHeader* header = runLoader->GetHeader();
1267 // Intialize generator
1269 generator->SetStack(stack);
1271 // Run main generator loop
1273 for (Int_t iev=0; iev<fNEvents; iev++)
1276 header->Reset(0,iev);
1277 runLoader->SetEventNumber(iev);
1279 runLoader->MakeTree("K");
1282 generator->Generate();
1285 header->SetNprimary(stack->GetNprimary());
1286 header->SetNtrack(stack->GetNtrack());
1287 stack->FinishEvent();
1288 header->SetStack(stack);
1289 runLoader->TreeE()->Fill();
1290 runLoader->WriteKinematics("OVERWRITE");
1294 generator->FinishRun();
1296 runLoader->WriteHeader("OVERWRITE");
1303 //_____________________________________________________________________________
1304 Bool_t AliSimulation::RunSDigitization(const char* detectors)
1306 // run the digitization and produce summable digits
1307 static Int_t eventNr=0;
1308 AliCodeTimerAuto("",0) ;
1310 // initialize CDB storage, run number, set CDB lock
1312 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1315 AliRunLoader* runLoader = LoadRun();
1316 if (!runLoader) return kFALSE;
1318 TString detStr = detectors;
1319 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1320 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1321 AliModule* det = (AliModule*) detArray->At(iDet);
1322 if (!det || !det->IsActive()) continue;
1323 if (IsSelected(det->GetName(), detStr)) {
1324 AliInfo(Form("creating summable digits for %s", det->GetName()));
1325 AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
1326 det->Hits2SDigits();
1327 AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
1328 AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
1332 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1333 AliError(Form("the following detectors were not found: %s",
1335 if (fStopOnError) return kFALSE;
1344 //_____________________________________________________________________________
1345 Bool_t AliSimulation::RunDigitization(const char* detectors,
1346 const char* excludeDetectors)
1348 // run the digitization and produce digits from sdigits
1350 AliCodeTimerAuto("",0)
1352 // initialize CDB storage, run number, set CDB lock
1354 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1357 delete AliRunLoader::Instance();
1362 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1363 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1364 AliDigitizationInput digInp(nStreams, signalPerBkgrd);
1365 // digInp.SetEmbeddingFlag(fEmbeddingFlag);
1366 digInp.SetRegionOfInterest(fRegionOfInterest);
1367 digInp.SetInputStream(0, fGAliceFileName.Data());
1368 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1369 const char* fileName = ((TObjString*)(fBkgrdFileNames->At(iStream-1)))->GetName();
1370 digInp.SetInputStream(iStream, fileName);
1373 detArr.SetOwner(kTRUE);
1374 TString detStr = detectors;
1375 TString detExcl = excludeDetectors;
1376 if (!static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice()) {
1377 AliError("Error occured while getting gAlice from Input 0");
1380 AliRunLoader* runLoader = AliRunLoader::GetRunLoader(digInp.GetInputStream(0)->GetFolderName());
1381 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1382 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1383 AliModule* det = (AliModule*) detArray->At(iDet);
1384 if (!det || !det->IsActive()) continue;
1385 if (!IsSelected(det->GetName(), detStr) || IsSelected(det->GetName(), detExcl)) continue;
1386 AliDigitizer* digitizer = det->CreateDigitizer(&digInp);
1387 if (!digitizer || !digitizer->Init()) {
1388 AliError(Form("no digitizer for %s", det->GetName()));
1389 if (fStopOnError) return kFALSE;
1392 detArr.AddLast(digitizer);
1393 AliInfo(Form("Created digitizer from SDigits -> Digits for %s", det->GetName()));
1396 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1397 AliError(Form("the following detectors were not found: %s", detStr.Data()));
1398 if (fStopOnError) return kFALSE;
1401 Int_t ndigs = detArr.GetEntriesFast();
1402 Int_t eventsCreated = 0;
1403 AliRunLoader* outRl = digInp.GetOutRunLoader();
1404 while ((eventsCreated++ < fNEvents) || (fNEvents < 0)) {
1405 if (!digInp.ConnectInputTrees()) break;
1406 digInp.InitEvent(); //this must be after call of Connect Input tress.
1407 if (outRl) outRl->SetEventNumber(eventsCreated-1);
1408 static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice(); // use gAlice of the first input stream
1409 for (int id=0;id<ndigs;id++) ((AliDigitizer*)detArr[id])->Digitize("");
1410 digInp.FinishEvent();
1412 digInp.FinishGlobal();
1417 //_____________________________________________________________________________
1418 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1420 // run the digitization and produce digits from hits
1422 AliCodeTimerAuto("",0)
1424 // initialize CDB storage, run number, set CDB lock
1426 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1429 AliRunLoader* runLoader = LoadRun("READ");
1430 if (!runLoader) return kFALSE;
1432 TString detStr = detectors;
1433 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1434 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1435 AliModule* det = (AliModule*) detArray->At(iDet);
1436 if (!det || !det->IsActive()) continue;
1437 if (IsSelected(det->GetName(), detStr)) {
1438 AliInfo(Form("creating digits from hits for %s", det->GetName()));
1443 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1444 AliError(Form("the following detectors were not found: %s",
1446 if (fStopOnError) return kFALSE;
1452 //_____________________________________________________________________________
1453 Bool_t AliSimulation::WriteRawData(const char* detectors,
1454 const char* fileName,
1455 Bool_t deleteIntermediateFiles,
1458 // convert the digits to raw data
1459 // First DDL raw data files for the given detectors are created.
1460 // If a file name is given, the DDL files are then converted to a DATE file.
1461 // If deleteIntermediateFiles is true, the DDL raw files are deleted
1463 // If the file name has the extension ".root", the DATE file is converted
1465 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1466 // 'selrawdata' flag can be used to enable writing of detectors raw data
1467 // accoring to the trigger cluster.
1469 AliCodeTimerAuto("",0)
1470 AliSysInfo::AddStamp("WriteRawData_Start");
1472 TString detStr = detectors;
1473 if (!WriteRawFiles(detStr.Data())) {
1474 if (fStopOnError) return kFALSE;
1476 AliSysInfo::AddStamp("WriteRawFiles");
1478 // run HLT simulation on simulated DDL raw files
1479 // and produce HLT ddl raw files to be included in date/root file
1480 // bugfix 2009-06-26: the decision whether to write HLT raw data
1481 // is taken in RunHLT. Here HLT always needs to be run in order to
1482 // create HLT digits, unless its switched off. This is due to the
1483 // special placement of the HLT between the generation of DDL files
1484 // and conversion to DATE/Root file.
1485 detStr.ReplaceAll("HLT", "");
1486 if (!fRunHLT.IsNull()) {
1488 if (fStopOnError) return kFALSE;
1491 AliSysInfo::AddStamp("WriteRawData_RunHLT");
1493 TString dateFileName(fileName);
1494 if (!dateFileName.IsNull()) {
1495 Bool_t rootOutput = dateFileName.EndsWith(".root");
1496 if (rootOutput) dateFileName += ".date";
1497 TString selDateFileName;
1499 selDateFileName = "selected.";
1500 selDateFileName+= dateFileName;
1502 if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1503 if (fStopOnError) return kFALSE;
1505 AliSysInfo::AddStamp("ConvertRawFilesToDate");
1506 if (deleteIntermediateFiles) {
1507 AliRunLoader* runLoader = LoadRun("READ");
1508 if (runLoader) for (Int_t iEvent = 0;
1509 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1511 snprintf(command, 256, "rm -r raw%d", iEvent);
1512 gSystem->Exec(command);
1518 if (!ConvertDateToRoot(dateFileName, fileName)) {
1519 if (fStopOnError) return kFALSE;
1521 AliSysInfo::AddStamp("ConvertDateToRoot");
1522 if (deleteIntermediateFiles) {
1523 gSystem->Unlink(dateFileName);
1526 TString selFileName = "selected.";
1527 selFileName += fileName;
1528 if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1529 if (fStopOnError) return kFALSE;
1531 if (deleteIntermediateFiles) {
1532 gSystem->Unlink(selDateFileName);
1541 //_____________________________________________________________________________
1542 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1544 // convert the digits to raw data DDL files
1546 AliCodeTimerAuto("",0)
1548 AliRunLoader* runLoader = LoadRun("READ");
1549 if (!runLoader) return kFALSE;
1551 // write raw data to DDL files
1552 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1553 AliInfo(Form("processing event %d", iEvent));
1554 runLoader->GetEvent(iEvent);
1555 TString baseDir = gSystem->WorkingDirectory();
1557 snprintf(dirName, 256, "raw%d", iEvent);
1558 gSystem->MakeDirectory(dirName);
1559 if (!gSystem->ChangeDirectory(dirName)) {
1560 AliError(Form("couldn't change to directory %s", dirName));
1561 if (fStopOnError) return kFALSE; else continue;
1564 ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1567 TString detStr = detectors;
1568 if (IsSelected("HLT", detStr)) {
1569 // Do nothing. "HLT" will be removed from detStr and HLT raw
1570 // data files are generated in RunHLT.
1573 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1574 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1575 AliModule* det = (AliModule*) detArray->At(iDet);
1576 if (!det || !det->IsActive()) continue;
1577 if (IsSelected(det->GetName(), detStr)) {
1578 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1583 if (!WriteTriggerRawData())
1584 if (fStopOnError) return kFALSE;
1586 gSystem->ChangeDirectory(baseDir);
1587 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1588 AliError(Form("the following detectors were not found: %s",
1590 if (fStopOnError) return kFALSE;
1599 //_____________________________________________________________________________
1600 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1601 const char* selDateFileName)
1603 // convert raw data DDL files to a DATE file with the program "dateStream"
1604 // The second argument is not empty when the user decides to write
1605 // the detectors raw data according to the trigger cluster.
1607 AliCodeTimerAuto("",0)
1609 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1611 AliError("the program dateStream was not found");
1612 if (fStopOnError) return kFALSE;
1617 AliRunLoader* runLoader = LoadRun("READ");
1618 if (!runLoader) return kFALSE;
1620 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1621 Bool_t selrawdata = kFALSE;
1622 if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1625 // Note the option -s. It is used in order to avoid
1626 // the generation of SOR/EOR events.
1627 snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d",
1628 dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1629 FILE* pipe = gSystem->OpenPipe(command, "w");
1632 AliError(Form("Cannot execute command: %s",command));
1636 Int_t selEvents = 0;
1637 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1639 UInt_t detectorPattern = 0;
1640 runLoader->GetEvent(iEvent);
1641 if (!runLoader->LoadTrigger()) {
1642 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1643 detectorPattern = aCTP->GetClusterMask();
1644 // Check if the event was triggered by CTP
1646 if (aCTP->GetClassMask()) selEvents++;
1650 AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1652 AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1653 selrawdata = kFALSE;
1657 fprintf(pipe, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
1661 // loop over detectors and DDLs
1662 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1663 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1665 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1666 Int_t ldcID = Int_t(ldc + 0.0001);
1667 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1669 char rawFileName[256];
1670 snprintf(rawFileName, 256, "raw%d/%s",
1671 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1673 // check existence and size of raw data file
1674 FILE* file = fopen(rawFileName, "rb");
1675 if (!file) continue;
1676 fseek(file, 0, SEEK_END);
1677 unsigned long size = ftell(file);
1679 if (!size) continue;
1681 if (ldcID != prevLDC) {
1682 fprintf(pipe, " LDC Id %d\n", ldcID);
1685 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1690 Int_t result = gSystem->ClosePipe(pipe);
1692 if (!(selrawdata && selEvents > 0)) {
1694 return (result == 0);
1697 AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1699 snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d",
1700 selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1701 FILE* pipe2 = gSystem->OpenPipe(command, "w");
1703 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1705 // Get the trigger decision and cluster
1706 UInt_t detectorPattern = 0;
1708 runLoader->GetEvent(iEvent);
1709 if (!runLoader->LoadTrigger()) {
1710 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1711 if (aCTP->GetClassMask() == 0) continue;
1712 detectorPattern = aCTP->GetClusterMask();
1713 detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1714 AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1717 fprintf(pipe2, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
1721 // loop over detectors and DDLs
1722 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1723 // Write only raw data from detectors that
1724 // are contained in the trigger cluster(s)
1725 if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1727 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1729 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1730 Int_t ldcID = Int_t(ldc + 0.0001);
1731 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1733 char rawFileName[256];
1734 snprintf(rawFileName, 256, "raw%d/%s",
1735 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1737 // check existence and size of raw data file
1738 FILE* file = fopen(rawFileName, "rb");
1739 if (!file) continue;
1740 fseek(file, 0, SEEK_END);
1741 unsigned long size = ftell(file);
1743 if (!size) continue;
1745 if (ldcID != prevLDC) {
1746 fprintf(pipe2, " LDC Id %d\n", ldcID);
1749 fprintf(pipe2, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1754 Int_t result2 = gSystem->ClosePipe(pipe2);
1757 return ((result == 0) && (result2 == 0));
1760 //_____________________________________________________________________________
1761 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1762 const char* rootFileName)
1764 // convert a DATE file to a root file with the program "alimdc"
1767 const Int_t kDBSize = 2000000000;
1768 const Int_t kTagDBSize = 1000000000;
1769 const Bool_t kFilter = kFALSE;
1770 const Int_t kCompression = 1;
1772 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1774 AliError("the program alimdc was not found");
1775 if (fStopOnError) return kFALSE;
1780 AliInfo(Form("converting DATE file %s to root file %s",
1781 dateFileName, rootFileName));
1783 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1784 const char* tagDBFS = "/tmp/mdc1/tags";
1786 // User defined file system locations
1787 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1788 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1789 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1790 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1791 if (gSystem->Getenv("ALIMDC_TAGDB"))
1792 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1794 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1795 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1796 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1798 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1799 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1800 gSystem->Exec(Form("mkdir %s",tagDBFS));
1802 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1803 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1804 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1806 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1807 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1808 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1810 return (result == 0);
1814 //_____________________________________________________________________________
1815 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1817 // delete existing run loaders, open a new one and load gAlice
1819 delete AliRunLoader::Instance();
1820 AliRunLoader* runLoader =
1821 AliRunLoader::Open(fGAliceFileName.Data(),
1822 AliConfig::GetDefaultEventFolderName(), mode);
1824 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1827 runLoader->LoadgAlice();
1828 runLoader->LoadHeader();
1829 gAlice = runLoader->GetAliRun();
1831 AliError(Form("no gAlice object found in file %s",
1832 fGAliceFileName.Data()));
1838 //_____________________________________________________________________________
1839 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1841 // get or calculate the number of signal events per background event
1843 if (!fBkgrdFileNames) return 1;
1844 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1845 if (nBkgrdFiles == 0) return 1;
1847 // get the number of signal events
1849 AliRunLoader* runLoader =
1850 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1851 if (!runLoader) return 1;
1853 nEvents = runLoader->GetNumberOfEvents();
1858 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1859 // get the number of background events
1860 const char* fileName = ((TObjString*)
1861 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1862 AliRunLoader* runLoader =
1863 AliRunLoader::Open(fileName, "BKGRD");
1864 if (!runLoader) continue;
1865 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1868 // get or calculate the number of signal per background events
1869 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1870 if (nSignalPerBkgrd <= 0) {
1871 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1872 } else if (result && (result != nSignalPerBkgrd)) {
1873 AliInfo(Form("the number of signal events per background event "
1874 "will be changed from %d to %d for stream %d",
1875 nSignalPerBkgrd, result, iBkgrdFile+1));
1876 nSignalPerBkgrd = result;
1879 if (!result) result = nSignalPerBkgrd;
1880 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1881 AliWarning(Form("not enough background events (%d) for %d signal events "
1882 "using %d signal per background events for stream %d",
1883 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1890 //_____________________________________________________________________________
1891 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1893 // check whether detName is contained in detectors
1894 // if yes, it is removed from detectors
1896 // check if all detectors are selected
1897 if ((detectors.CompareTo("ALL") == 0) ||
1898 detectors.BeginsWith("ALL ") ||
1899 detectors.EndsWith(" ALL") ||
1900 detectors.Contains(" ALL ")) {
1905 // search for the given detector
1906 Bool_t result = kFALSE;
1907 if ((detectors.CompareTo(detName) == 0) ||
1908 detectors.BeginsWith(detName+" ") ||
1909 detectors.EndsWith(" "+detName) ||
1910 detectors.Contains(" "+detName+" ")) {
1911 detectors.ReplaceAll(detName, "");
1915 // clean up the detectors string
1916 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1917 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1918 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1923 //_____________________________________________________________________________
1924 Int_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName, Int_t N, Int_t nSkip)
1927 // Steering routine to convert raw data in directory rawDirectory/ to fake SDigits.
1928 // These can be used for embedding of MC tracks into RAW data using the standard
1929 // merging procedure.
1931 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1934 AliError("no gAlice object. Restart aliroot and try again.");
1937 if (gAlice->Modules()->GetEntries() > 0) {
1938 AliError("gAlice was already run. Restart aliroot and try again.");
1942 AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1946 gROOT->LoadMacro(fConfigFileName.Data());
1947 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1949 if(AliCDBManager::Instance()->GetRun() >= 0) {
1950 SetRunNumber(AliCDBManager::Instance()->GetRun());
1952 AliWarning("Run number not initialized!!");
1955 AliRunLoader::Instance()->CdGAFile();
1957 AliPDG::AddParticlesToPdgDataBase();
1959 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1961 gAlice->GetMCApp()->Init();
1963 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1964 gAlice->InitLoaders();
1965 AliRunLoader::Instance()->MakeTree("E");
1966 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1967 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1968 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1971 // Save stuff at the beginning of the file to avoid file corruption
1972 AliRunLoader::Instance()->CdGAFile();
1977 //AliCDBManager* man = AliCDBManager::Instance();
1978 //man->SetRun(0); // Should this come from rawdata header ?
1982 // Get the runloader
1983 AliRunLoader* runLoader = AliRunLoader::Instance();
1985 // Open esd file if available
1988 AliESDEvent* esd = 0;
1989 if (esdFileName && (strlen(esdFileName)>0)) {
1990 esdFile = TFile::Open(esdFileName);
1992 esd = new AliESDEvent();
1993 esdFile->GetObject("esdTree", treeESD);
1995 esd->ReadFromTree(treeESD);
1997 AliInfo(Form("Asking to skip first %d ESDs events",nSkip));
2006 // Create the RawReader
2007 TString fileName(rawDirectory);
2008 AliRawReader* rawReader = AliRawReader::Create(fileName.Data());
2009 if (!rawReader) return (kFALSE);
2011 // if (!fEquipIdMap.IsNull() && fRawReader)
2012 // fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
2014 // Get list of detectors
2015 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
2018 AliHeader* header = runLoader->GetHeader();
2022 if (!(rawReader->NextEvent())) break;
2023 runLoader->SetEventNumber(nev);
2024 runLoader->GetHeader()->Reset(rawReader->GetRunNumber(),
2026 runLoader->GetEvent(nev);
2027 AliInfo(Form("We are at event %d",nev));
2030 TString detStr = fMakeSDigits;
2031 for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
2032 AliModule* det = (AliModule*) detArray->At(iDet);
2033 if (!det || !det->IsActive()) continue;
2034 if (IsSelected(det->GetName(), detStr)) {
2035 AliInfo(Form("Calling Raw2SDigits for %s", det->GetName()));
2036 det->Raw2SDigits(rawReader);
2043 // If ESD information available obtain reconstructed vertex and store in header.
2045 AliInfo(Form("Selected event %d correspond to event %d in raw and to %d in esd",nev,rawReader->GetEventIndex(),nSkip+rawReader->GetEventIndex()));
2046 treeESD->GetEvent(nSkip+rawReader->GetEventIndex());
2047 const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
2048 Double_t position[3];
2049 esdVertex->GetXYZ(position);
2050 AliGenEventHeader* mcHeader = new AliGenEventHeader("ESD");
2053 for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
2054 mcHeader->SetPrimaryVertex(mcV);
2055 header->Reset(0,nev);
2056 header->SetGenEventHeader(mcHeader);
2057 AliInfo(Form("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]));
2061 runLoader->TreeE()->Fill();
2062 AliInfo(Form("Finished event %d",nev));
2071 runLoader->CdGAFile();
2072 runLoader->WriteHeader("OVERWRITE");
2073 runLoader->WriteRunLoader();
2078 //_____________________________________________________________________________
2079 void AliSimulation::FinishRun()
2082 // Called at the end of the run.
2087 AliDebug(1, "Finish Lego");
2088 AliRunLoader::Instance()->CdGAFile();
2092 // Clean detector information
2093 TIter next(gAlice->Modules());
2094 AliModule *detector;
2095 while((detector = dynamic_cast<AliModule*>(next()))) {
2096 AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
2097 detector->FinishRun();
2100 AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
2101 AliRunLoader::Instance()->WriteHeader("OVERWRITE");
2103 // Write AliRun info and all detectors parameters
2104 AliRunLoader::Instance()->CdGAFile();
2105 gAlice->Write(0,TObject::kOverwrite);//write AliRun
2106 AliRunLoader::Instance()->Write(0,TObject::kOverwrite);//write RunLoader itself
2108 if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();
2109 AliRunLoader::Instance()->Synchronize();
2112 //_____________________________________________________________________________
2113 Int_t AliSimulation::GetDetIndex(const char* detector)
2115 // return the detector index corresponding to detector
2117 for (index = 0; index < fgkNDetectors ; index++) {
2118 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2124 //_____________________________________________________________________________
2125 Bool_t AliSimulation::CreateHLT()
2127 // Init the HLT simulation.
2128 // The function loads the library and creates the instance of AliHLTSimulation.
2129 // the main reason for the decoupled creation is to set the transient OCDB
2130 // objects before the OCDB is locked
2132 // load the library dynamically
2133 gSystem->Load(ALIHLTSIMULATION_LIBRARY);
2135 // check for the library version
2136 AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
2138 AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
2141 if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
2142 AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
2145 // print compile info
2146 typedef void (*CompileInfo)( const char*& date, const char*& time);
2147 CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
2149 const char* date="";
2150 const char* time="";
2151 (*fctInfo)(date, time);
2152 if (!date) date="unknown";
2153 if (!time) time="unknown";
2154 AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
2156 AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
2159 // create instance of the HLT simulation
2160 AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
2161 if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
2162 AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
2166 TString specObjects;
2167 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
2168 if (specObjects.Length()>0) specObjects+=" ";
2169 specObjects+=fSpecCDBUri[i]->GetName();
2172 AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
2173 if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
2174 AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
2180 //_____________________________________________________________________________
2181 Bool_t AliSimulation::RunHLT()
2183 // Run the HLT simulation
2184 // HLT simulation is implemented in HLT/sim/AliHLTSimulation
2185 // Disabled if fRunHLT is empty, default vaule is "default".
2186 // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
2187 // The default simulation depends on the HLT component libraries and their
2188 // corresponding agents which define components and chains to run. See
2189 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
2190 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
2192 // The libraries to be loaded can be specified as an option.
2194 // AliSimulation sim;
2195 // sim.SetRunHLT("libAliHLTSample.so");
2197 // will only load <tt>libAliHLTSample.so</tt>
2199 // Other available options:
2200 // \li loglevel=<i>level</i> <br>
2201 // logging level for this processing
2203 // disable redirection of log messages to AliLog class
2204 // \li config=<i>macro</i>
2205 // configuration macro
2206 // \li chains=<i>configuration</i>
2207 // comma separated list of configurations to be run during simulation
2208 // \li rawfile=<i>file</i>
2209 // source for the RawReader to be created, the default is <i>./</i> if
2210 // raw data is simulated
2214 if (!fpHLT && !CreateHLT()) {
2217 AliHLTSimulation* pHLT=fpHLT;
2219 AliRunLoader* pRunLoader = LoadRun("READ");
2220 if (!pRunLoader) return kFALSE;
2222 // initialize CDB storage, run number, set CDB lock
2223 // thats for the case of running HLT simulation without all the other steps
2224 // multiple calls are handled by the function, so we can just call
2226 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
2229 // init the HLT simulation
2231 if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2232 TString detStr = fWriteRawData;
2233 if (!IsSelected("HLT", detStr)) {
2234 options+=" writerawfiles=";
2236 options+=" writerawfiles=HLT";
2239 if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2240 // as a matter of fact, HLT will run reconstruction and needs the RawReader
2241 // in order to get detector data. By default, RawReaderFile is used to read
2242 // the already simulated ddl files. Date and Root files from the raw data
2243 // are generated after the HLT simulation.
2244 options+=" rawfile=./";
2247 AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2248 if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2249 AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2251 // run the HLT simulation
2252 AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2253 if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2254 AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2258 // delete the instance
2259 AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2260 if (fctDelete==NULL || fctDelete(pHLT)<0) {
2261 AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2265 return iResult>=0?kTRUE:kFALSE;
2268 //_____________________________________________________________________________
2269 Bool_t AliSimulation::RunQA()
2271 // run the QA on summable hits, digits or digits
2273 //if(!gAlice) return kFALSE;
2274 AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2276 TString detectorsw("") ;
2278 AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2279 detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ;
2280 if ( detectorsw.IsNull() )
2285 //_____________________________________________________________________________
2286 Bool_t AliSimulation::SetRunQA(TString detAndAction)
2288 // Allows to run QA for a selected set of detectors
2289 // and a selected set of tasks among HITS, SDIGITS and DIGITS
2290 // all selected detectors run the same selected tasks
2292 if (!detAndAction.Contains(":")) {
2293 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2297 Int_t colon = detAndAction.Index(":") ;
2298 fQADetectors = detAndAction(0, colon) ;
2299 if (fQADetectors.Contains("ALL") ){
2300 TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2301 Int_t minus = fQADetectors.Last('-') ;
2302 TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2303 TString toRemove("") ;
2304 while (minus >= 0) {
2305 toRemove = fQADetectors(minus+1, fQADetectors.Length()) ;
2306 toRemove = toRemove.Strip() ;
2307 toKeep.ReplaceAll(toRemove, "") ;
2308 fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ;
2309 minus = fQADetectors.Last('-') ;
2311 fQADetectors = toKeep ;
2313 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2314 if (fQATasks.Contains("ALL") ) {
2315 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
2317 fQATasks.ToUpper() ;
2319 if ( fQATasks.Contains("HIT") )
2320 tempo = Form("%d ", AliQAv1::kHITS) ;
2321 if ( fQATasks.Contains("SDIGIT") )
2322 tempo += Form("%d ", AliQAv1::kSDIGITS) ;
2323 if ( fQATasks.Contains("DIGIT") )
2324 tempo += Form("%d ", AliQAv1::kDIGITS) ;
2326 if (fQATasks.IsNull()) {
2327 AliInfo("No QA requested\n") ;
2332 TString tempo(fQATasks) ;
2333 tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS)) ;
2334 tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;
2335 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;
2336 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2338 AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ;
2339 AliQAManager::QAManager()->SetTasks(fQATasks) ;
2340 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++)
2341 AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2346 //_____________________________________________________________________________
2347 void AliSimulation::ProcessEnvironmentVars()
2349 // Extract run number and random generator seed from env variables
2351 AliInfo("Processing environment variables");
2353 // Random Number seed
2355 // first check that seed is not already set
2357 if (gSystem->Getenv("CONFIG_SEED")) {
2358 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2361 if (gSystem->Getenv("CONFIG_SEED")) {
2362 AliInfo(Form("Seed for random number generation already set (%d)"
2363 ": CONFIG_SEED variable ignored!", fSeed));
2367 AliInfo(Form("Seed for random number generation = %d ", fSeed));
2371 // first check that run number is not already set
2373 if (gSystem->Getenv("DC_RUN")) {
2374 fRun = atoi(gSystem->Getenv("DC_RUN"));
2377 if (gSystem->Getenv("DC_RUN")) {
2378 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2382 AliInfo(Form("Run number = %d", fRun));
2385 //---------------------------------------------------------------------
2386 void AliSimulation::WriteGRPEntry()
2388 // Get the necessary information from galice (generator, trigger etc) and
2389 // write a GRP entry corresponding to the settings in the Config.C used
2390 // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2393 AliInfo("Writing global run parameters entry into the OCDB");
2395 AliGRPObject* grpObj = new AliGRPObject();
2397 grpObj->SetRunType("PHYSICS");
2398 grpObj->SetTimeStart(fTimeStart);
2399 grpObj->SetTimeEnd(fTimeEnd);
2400 grpObj->SetBeamEnergyIsSqrtSHalfGeV(); // new format of GRP: store sqrt(s)/2 in GeV
2402 const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2408 gen->GetProjectile(projectile,a,z);
2410 gen->GetTarget(target,a,z);
2411 TString beamType = projectile + "-" + target;
2412 beamType.ReplaceAll(" ","");
2413 if (!beamType.CompareTo("-")) {
2414 grpObj->SetBeamType("UNKNOWN");
2415 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
2418 grpObj->SetBeamType(beamType);
2420 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 * a / z);
2422 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 );
2424 // Heavy ion run, the event specie is set to kHighMult
2425 fEventSpecie = AliRecoParam::kHighMult;
2426 if ((strcmp(beamType,"p-p") == 0) ||
2427 (strcmp(beamType,"p-") == 0) ||
2428 (strcmp(beamType,"-p") == 0) ||
2429 (strcmp(beamType,"P-P") == 0) ||
2430 (strcmp(beamType,"P-") == 0) ||
2431 (strcmp(beamType,"-P") == 0)) {
2432 // Proton run, the event specie is set to kLowMult
2433 fEventSpecie = AliRecoParam::kLowMult;
2437 AliWarning("Unknown beam type and energy! Setting energy to 0");
2438 grpObj->SetBeamEnergy(0);
2439 grpObj->SetBeamType("UNKNOWN");
2442 UInt_t detectorPattern = 0;
2444 TObjArray *detArray = gAlice->Detectors();
2445 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2446 if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2447 AliDebug(1, Form("Detector #%d found: %s", iDet, AliDAQ::OfflineModuleName(iDet)));
2448 detectorPattern |= (1 << iDet);
2453 if (!fTriggerConfig.IsNull())
2454 detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2457 if (!fRunHLT.IsNull())
2458 detectorPattern |= (1 << AliDAQ::kHLTId);
2460 grpObj->SetNumberOfDetectors((Char_t)nDets);
2461 grpObj->SetDetectorMask((Int_t)detectorPattern);
2462 grpObj->SetLHCPeriod("LHC08c");
2463 grpObj->SetLHCState("STABLE_BEAMS");
2465 AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2466 Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2468 Float_t factorSol = field ? field->GetFactorSol() : 0;
2469 Float_t currentSol = TMath::Abs(factorSol)>1E-6 ?
2470 TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
2472 Float_t factorDip = field ? field->GetFactorDip() : 0;
2473 Float_t currentDip = 6000.*TMath::Abs(factorDip);
2475 grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
2476 grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);
2477 grpObj->SetL3Polarity(factorSol>0 ? 0:1);
2478 grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
2479 if (field) grpObj->SetUniformBMap(field->IsUniform()); // for special MC with k5kGUniform map
2480 grpObj->SetPolarityConventionLHC(); // LHC convention +/+ current -> -/- field main components
2482 grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2484 //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2486 // Now store the entry in OCDB
2487 AliCDBManager* man = AliCDBManager::Instance();
2489 man->SetLock(0, fKey);
2491 AliCDBStorage* sto = man->GetStorage(fGRPWriteLocation.Data());
2494 AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2495 AliCDBMetaData *metadata= new AliCDBMetaData();
2497 metadata->SetResponsible("alice-off@cern.ch");
2498 metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2500 sto->Put(grpObj,id,metadata);
2501 man->SetLock(1, fKey);
2504 //_____________________________________________________________________________
2505 time_t AliSimulation::GenerateTimeStamp() const
2507 // Generate event time-stamp according to
2508 // SOR/EOR time from GRP
2509 if (fUseTimeStampFromCDB)
2510 return fTimeStart + gRandom->Integer(fTimeEnd-fTimeStart);