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 64623 2013-10-21 13:38:58Z 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",
160 "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE","AD",
163 //_____________________________________________________________________________
164 AliSimulation::AliSimulation(const char* configFileName,
165 const char* name, const char* title) :
168 fRunGeneratorOnly(kFALSE),
169 fRunGeneration(kTRUE),
170 fRunSimulation(kTRUE),
171 fLoadAlignFromCDB(kTRUE),
172 fLoadAlObjsListOfDets("ALL"),
176 fMakeDigitsFromHits(""),
178 fRawDataFileName(""),
179 fDeleteIntermediateFiles(kFALSE),
180 fWriteSelRawData(kFALSE),
181 fStopOnError(kFALSE),
182 fUseMonitoring(kFALSE),
184 fConfigFileName(configFileName),
185 fGAliceFileName("galice.root"),
187 fBkgrdFileNames(NULL),
188 fAlignObjArray(NULL),
189 fUseBkgrdVertex(kTRUE),
190 fRegionOfInterest(kFALSE),
196 fInitCDBCalled(kFALSE),
197 fInitRunNumberCalled(kFALSE),
198 fSetRunNumberFromDataCalled(kFALSE),
199 fEmbeddingFlag(kFALSE),
202 fUseVertexFromCDB(0),
203 fUseMagFieldFromGRP(0),
204 fGRPWriteLocation(Form("local://%s", gSystem->pwd())),
205 fUseTimeStampFromCDB(0),
211 fEventSpecie(AliRecoParam::kDefault),
212 fWriteQAExpertData(kTRUE),
216 fWriteGRPEntry(kTRUE)
218 // create simulation object with default parameters
220 SetGAliceFile("galice.root");
223 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
224 qam->SetActiveDetectors(fQADetectors) ;
225 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
226 qam->SetTasks(fQATasks) ;
229 //_____________________________________________________________________________
230 AliSimulation::~AliSimulation()
234 fEventsPerFile.Delete();
235 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
236 // delete fAlignObjArray; fAlignObjArray=0;
238 if (fBkgrdFileNames) {
239 fBkgrdFileNames->Delete();
240 delete fBkgrdFileNames;
243 fSpecCDBUri.Delete();
244 if (fgInstance==this) fgInstance = 0;
246 AliQAManager::QAManager()->ShowQA() ;
247 AliQAManager::Destroy() ;
248 AliCodeTimer::Instance()->Print();
252 //_____________________________________________________________________________
253 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
255 // set the number of events for one run
260 //_____________________________________________________________________________
261 void AliSimulation::InitQA()
263 // activate a default CDB storage
264 // First check if we have any CDB storage set, because it is used
265 // to retrieve the calibration and alignment constants
267 if (fInitCDBCalled) return;
268 fInitCDBCalled = kTRUE;
270 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
271 qam->SetActiveDetectors(fQADetectors) ;
272 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
273 qam->SetTasks(fQATasks) ;
274 if (fWriteQAExpertData)
275 qam->SetWriteExpert() ;
277 if (qam->IsDefaultStorageSet()) {
278 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
279 AliWarning("Default QA reference storage has been already set !");
280 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fQARefUri.Data()));
281 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
282 fQARefUri = qam->GetDefaultStorage()->GetURI();
284 if (fQARefUri.Length() > 0) {
285 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
286 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
287 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
289 fQARefUri="local://$ALICE_ROOT/QARef";
290 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
291 AliWarning("Default QA reference storage not yet set !!!!");
292 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
293 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
295 qam->SetDefaultStorage(fQARefUri);
299 //_____________________________________________________________________________
300 void AliSimulation::InitCDB()
302 // activate a default CDB storage
303 // First check if we have any CDB storage set, because it is used
304 // to retrieve the calibration and alignment constants
306 if (fInitCDBCalled) return;
307 fInitCDBCalled = kTRUE;
309 AliCDBManager* man = AliCDBManager::Instance();
310 if (man->IsDefaultStorageSet())
312 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
313 AliWarning("Default CDB storage has been already set !");
314 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
315 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
316 fCDBUri = man->GetDefaultStorage()->GetURI();
319 if (fCDBUri.Length() > 0)
321 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
322 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
323 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
325 fCDBUri="local://$ALICE_ROOT/OCDB";
326 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
327 AliWarning("Default CDB storage not yet set !!!!");
328 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
329 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
332 man->SetDefaultStorage(fCDBUri);
335 // Now activate the detector specific CDB storage locations
336 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
337 TObject* obj = fSpecCDBUri[i];
339 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
340 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
341 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
342 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
347 //_____________________________________________________________________________
348 void AliSimulation::InitRunNumber(){
349 // check run number. If not set, set it to 0 !!!!
351 if (fInitRunNumberCalled) return;
352 fInitRunNumberCalled = kTRUE;
354 AliCDBManager* man = AliCDBManager::Instance();
355 if (man->GetRun() >= 0)
357 AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
358 "Use external variable DC_RUN or AliSimulation::SetRun()!"));
362 AliDebug(2, Form("Setting CDB run number to: %d",fRun));
365 AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
374 //_____________________________________________________________________________
375 void AliSimulation::SetCDBLock() {
376 // Set CDB lock: from now on it is forbidden to reset the run number
377 // or the default storage or to activate any further storage!
379 ULong64_t key = AliCDBManager::Instance()->SetLock(1);
383 //_____________________________________________________________________________
384 void AliSimulation::SetDefaultStorage(const char* uri) {
385 // Store the desired default CDB storage location
386 // Activate it later within the Run() method
392 //_____________________________________________________________________________
393 void AliSimulation::SetQARefDefaultStorage(const char* uri) {
394 // Store the desired default CDB storage location
395 // Activate it later within the Run() method
398 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
401 //_____________________________________________________________________________
402 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
403 // Store a detector-specific CDB storage location
404 // Activate it later within the Run() method
406 AliCDBPath aPath(calibType);
407 if(!aPath.IsValid()){
408 AliError(Form("Not a valid path: %s", calibType));
412 TObject* obj = fSpecCDBUri.FindObject(calibType);
413 if (obj) fSpecCDBUri.Remove(obj);
414 fSpecCDBUri.Add(new TNamed(calibType, uri));
418 //_____________________________________________________________________________
419 void AliSimulation::SetRunNumber(Int_t run)
422 // Activate it later within the Run() method
427 //_____________________________________________________________________________
428 void AliSimulation::SetSeed(Int_t seed)
431 // Activate it later within the Run() method
436 //_____________________________________________________________________________
437 Bool_t AliSimulation::SetRunNumberFromData()
439 // Set the CDB manager run number
440 // The run number is retrieved from gAlice
442 if (fSetRunNumberFromDataCalled) return kTRUE;
443 fSetRunNumberFromDataCalled = kTRUE;
445 AliCDBManager* man = AliCDBManager::Instance();
446 Int_t runData = -1, runCDB = -1;
448 AliRunLoader* runLoader = LoadRun("READ");
449 if (!runLoader) return kFALSE;
451 runData = runLoader->GetHeader()->GetRun();
455 runCDB = man->GetRun();
457 if (runCDB != runData) {
458 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
459 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
460 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
461 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
466 man->SetRun(runData);
469 if(man->GetRun() < 0) {
470 AliError("Run number not properly initalized!");
479 //_____________________________________________________________________________
480 void AliSimulation::SetConfigFile(const char* fileName)
482 // set the name of the config file
484 fConfigFileName = fileName;
487 //_____________________________________________________________________________
488 void AliSimulation::SetGAliceFile(const char* fileName)
490 // set the name of the galice file
491 // the path is converted to an absolute one if it is relative
493 fGAliceFileName = fileName;
494 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
495 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
497 fGAliceFileName = absFileName;
498 delete[] absFileName;
501 AliDebug(2, Form("galice file name set to %s", fileName));
504 //_____________________________________________________________________________
505 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
508 // set the number of events per file for the given detector and data type
509 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
511 TNamed* obj = new TNamed(detector, type);
512 obj->SetUniqueID(nEvents);
513 fEventsPerFile.Add(obj);
516 //_____________________________________________________________________________
517 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
519 // Read the alignment objects from CDB.
520 // Each detector is supposed to have the
521 // alignment objects in DET/Align/Data CDB path.
522 // All the detector objects are then collected,
523 // sorted by geometry level (starting from ALIC) and
524 // then applied to the TGeo geometry.
525 // Finally an overlaps check is performed.
527 if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
528 AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
532 // initialize CDB storage, run number, set CDB lock
534 // if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
537 Bool_t delRunLoader = kFALSE;
539 runLoader = LoadRun("READ");
540 if (!runLoader) return kFALSE;
541 delRunLoader = kTRUE;
544 // Export ideal geometry
545 if(!IsGeometryFromFile()) AliGeomManager::GetGeometry()->Export("geometry.root");
547 // Load alignment data from CDB and apply to geometry through AliGeomManager
548 if(fLoadAlignFromCDB){
550 TString detStr = fLoadAlObjsListOfDets;
551 TString loadAlObjsListOfDets = "";
553 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
554 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
555 AliModule* det = (AliModule*) detArray->At(iDet);
556 if (!det || !det->IsActive()) continue;
557 if (IsSelected(det->GetName(), detStr)) {
558 //add det to list of dets to be aligned from CDB
559 loadAlObjsListOfDets += det->GetName();
560 loadAlObjsListOfDets += " ";
562 } // end loop over detectors
563 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
564 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
566 // Check if the array with alignment objects was
567 // provided by the user. If yes, apply the objects
568 // to the present TGeo geometry
569 if (fAlignObjArray) {
570 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
571 AliError("The misalignment of one or more volumes failed!"
572 "Compare the list of simulated detectors and the list of detector alignment data!");
573 if (delRunLoader) delete runLoader;
579 // Update the internal geometry of modules (ITS needs it)
580 TString detStr = fLoadAlObjsListOfDets;
581 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
582 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
584 AliModule* det = (AliModule*) detArray->At(iDet);
585 if (!det || !det->IsActive()) continue;
586 if (IsSelected(det->GetName(), detStr)) {
587 det->UpdateInternalGeometry();
589 } // end loop over detectors
592 if (delRunLoader) delete runLoader;
597 //_____________________________________________________________________________
598 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
600 // add a file with background events for merging
602 TObjString* fileNameStr = new TObjString(fileName);
603 fileNameStr->SetUniqueID(nSignalPerBkgrd);
604 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
605 fBkgrdFileNames->Add(fileNameStr);
608 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
610 // add a file with background events for embeddin
611 MergeWith(fileName, nSignalPerBkgrd);
612 fEmbeddingFlag = kTRUE;
615 //_____________________________________________________________________________
616 Bool_t AliSimulation::Run(Int_t nEvents)
618 // run the generation, simulation and digitization
621 AliCodeTimerAuto("",0)
622 AliSysInfo::AddStamp("Start_Run");
624 // Load run number and seed from environmental vars
625 ProcessEnvironmentVars();
626 AliSysInfo::AddStamp("ProcessEnvironmentVars");
628 gRandom->SetSeed(fSeed);
630 if (nEvents > 0) fNEvents = nEvents;
632 // Run generator-only code on demand
633 if (fRunGeneratorOnly)
635 if(!RunGeneratorOnly())
637 if (fStopOnError) return kFALSE;
643 // create and setup the HLT instance
644 if (!fRunHLT.IsNull() && !CreateHLT()) {
645 if (fStopOnError) return kFALSE;
650 // generation and simulation -> hits
651 if (fRunGeneration) {
652 if (!RunSimulation()) if (fStopOnError) return kFALSE;
654 AliSysInfo::AddStamp("RunSimulation");
656 // initialize CDB storage from external environment
657 // (either CDB manager or AliSimulation setters),
658 // if not already done in RunSimulation()
660 AliSysInfo::AddStamp("InitCDB");
662 // Set run number in CDBManager from data
663 // From this point on the run number must be always loaded from data!
664 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
666 // Set CDB lock: from now on it is forbidden to reset the run number
667 // or the default storage or to activate any further storage!
670 // If RunSimulation was not called, load the geometry and misalign it
671 if (!AliGeomManager::GetGeometry()) {
672 // Initialize the geometry manager
673 AliGeomManager::LoadGeometry("geometry.root");
674 AliSysInfo::AddStamp("GetGeometry");
675 // // Check that the consistency of symbolic names for the activated subdetectors
676 // // in the geometry loaded by AliGeomManager
677 // AliRunLoader* runLoader = LoadRun("READ");
678 // if (!runLoader) return kFALSE;
680 // TString detsToBeChecked = "";
681 // TObjArray* detArray = runLoader->GetAliRun()->Detectors();
682 // for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
683 // AliModule* det = (AliModule*) detArray->At(iDet);
684 // if (!det || !det->IsActive()) continue;
685 // detsToBeChecked += det->GetName();
686 // detsToBeChecked += " ";
687 // } // end loop over detectors
688 // if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
689 if(!AliGeomManager::CheckSymNamesLUT("ALL"))
690 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
692 if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
694 if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
696 AliSysInfo::AddStamp("MissalignGeometry");
699 // hits -> summable digits
700 AliSysInfo::AddStamp("Start_sdigitization");
701 if (!fMakeSDigits.IsNull()) {
702 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
705 AliSysInfo::AddStamp("Stop_sdigitization");
707 AliSysInfo::AddStamp("Start_digitization");
708 // summable digits -> digits
709 if (!fMakeDigits.IsNull()) {
710 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
711 if (fStopOnError) return kFALSE;
714 AliSysInfo::AddStamp("Stop_digitization");
719 if (!fMakeDigitsFromHits.IsNull()) {
720 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
721 AliWarning(Form("Merging and direct creation of digits from hits "
722 "was selected for some detectors. "
723 "No merging will be done for the following detectors: %s",
724 fMakeDigitsFromHits.Data()));
726 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
727 if (fStopOnError) return kFALSE;
731 AliSysInfo::AddStamp("Hits2Digits");
735 if (!fTriggerConfig.IsNull() && !RunTrigger(fTriggerConfig,fMakeDigits)) {
736 if (fStopOnError) return kFALSE;
739 AliSysInfo::AddStamp("RunTrigger");
742 // digits -> raw data
743 if (!fWriteRawData.IsNull()) {
744 if (!WriteRawData(fWriteRawData, fRawDataFileName,
745 fDeleteIntermediateFiles,fWriteSelRawData)) {
746 if (fStopOnError) return kFALSE;
750 AliSysInfo::AddStamp("WriteRaw");
752 // run HLT simulation on simulated digit data if raw data is not
753 // simulated, otherwise its called as part of WriteRawData
754 if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
756 if (fStopOnError) return kFALSE;
760 AliSysInfo::AddStamp("RunHLT");
764 Bool_t rv = RunQA() ;
770 AliSysInfo::AddStamp("RunQA");
772 TString snapshotFileOut("");
773 if(TString(gSystem->Getenv("OCDB_SNAPSHOT_CREATE")) == TString("kTRUE")){
774 AliInfo(" ******** Creating the snapshot! *********");
775 TString snapshotFile(gSystem->Getenv("OCDB_SNAPSHOT_FILENAME"));
776 if(!(snapshotFile.IsNull() || snapshotFile.IsWhitespace()))
777 snapshotFileOut = snapshotFile;
779 snapshotFileOut="OCDB.root";
780 AliCDBManager::Instance()->DumpToSnapshotFile(snapshotFileOut.Data(),kFALSE);
783 // Cleanup of CDB manager: cache and active storages!
784 AliCDBManager::Instance()->ClearCache();
789 //_______________________________________________________________________
790 Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
791 Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
792 Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
795 // Generates lego plots of:
796 // - radiation length map phi vs theta
797 // - radiation length map phi vs eta
798 // - interaction length map
799 // - g/cm2 length map
801 // ntheta bins in theta, eta
802 // themin minimum angle in theta (degrees)
803 // themax maximum angle in theta (degrees)
805 // phimin minimum angle in phi (degrees)
806 // phimax maximum angle in phi (degrees)
807 // rmin minimum radius
808 // rmax maximum radius
811 // The number of events generated = ntheta*nphi
812 // run input parameters in macro setup (default="Config.C")
814 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
817 <img src="picts/AliRunLego1.gif">
822 <img src="picts/AliRunLego2.gif">
827 <img src="picts/AliRunLego3.gif">
832 // run the generation and simulation
834 AliCodeTimerAuto("",0)
836 // initialize CDB storage and run number from external environment
837 // (either CDB manager or AliSimulation setters)
843 AliError("no gAlice object. Restart aliroot and try again.");
846 if (gAlice->Modules()->GetEntries() > 0) {
847 AliError("gAlice was already run. Restart aliroot and try again.");
851 AliInfo(Form("initializing gAlice with config file %s",
852 fConfigFileName.Data()));
855 if (nev == -1) nev = nc1 * nc2;
857 // check if initialisation has been done
858 // If runloader has been initialized, set the number of events per file to nc1 * nc2
861 if (!gener) gener = new AliLegoGenerator();
863 // Configure Generator
865 gener->SetRadiusRange(rmin, rmax);
866 gener->SetZMax(zmax);
867 gener->SetCoor1Range(nc1, c1min, c1max);
868 gener->SetCoor2Range(nc2, c2min, c2max);
872 fLego = new AliLego("lego",gener);
874 //__________________________________________________________________________
878 gROOT->LoadMacro(setup);
879 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
881 if(AliCDBManager::Instance()->GetRun() >= 0) {
882 SetRunNumber(AliCDBManager::Instance()->GetRun());
884 AliWarning("Run number not initialized!!");
887 AliRunLoader::Instance()->CdGAFile();
889 AliPDG::AddParticlesToPdgDataBase();
891 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
893 gAlice->GetMCApp()->Init();
896 //Must be here because some MCs (G4) adds detectors here and not in Config.C
897 gAlice->InitLoaders();
898 AliRunLoader::Instance()->MakeTree("E");
901 // Save stuff at the beginning of the file to avoid file corruption
902 AliRunLoader::Instance()->CdGAFile();
905 //Save current generator
906 AliGenerator *gen=gAlice->GetMCApp()->Generator();
907 gAlice->GetMCApp()->ResetGenerator(gener);
908 //Prepare MC for Lego Run
914 AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
915 gMC->ProcessRun(nev);
917 // End of this run, close files
919 // Restore current generator
920 gAlice->GetMCApp()->ResetGenerator(gen);
921 // Delete Lego Object
927 //_____________________________________________________________________________
928 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
932 AliCodeTimerAuto("",0)
934 // initialize CDB storage from external environment
935 // (either CDB manager or AliSimulation setters),
936 // if not already done in RunSimulation()
939 // Set run number in CDBManager from data
940 // From this point on the run number must be always loaded from data!
941 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
943 // Set CDB lock: from now on it is forbidden to reset the run number
944 // or the default storage or to activate any further storage!
947 AliRunLoader* runLoader = LoadRun("READ");
948 if (!runLoader) return kFALSE;
949 TString trconfiguration = config;
951 if (trconfiguration.IsNull()) {
952 if(!fTriggerConfig.IsNull()) {
953 trconfiguration = fTriggerConfig;
956 AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
959 runLoader->MakeTree( "GG" );
960 AliCentralTrigger* aCTP = runLoader->GetTrigger();
961 // Load Configuration
962 if (!aCTP->LoadConfiguration( trconfiguration ))
966 if( !aCTP->RunTrigger( runLoader , detectors ) ) {
978 //_____________________________________________________________________________
979 Bool_t AliSimulation::WriteTriggerRawData()
981 // Writes the CTP (trigger) DDL raw data
982 // Details of the format are given in the
983 // trigger TDR - pages 134 and 135.
984 AliCTPRawData writer;
990 //_____________________________________________________________________________
991 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
993 // run the generation and simulation
995 AliCodeTimerAuto("",0)
997 // initialize CDB storage and run number from external environment
998 // (either CDB manager or AliSimulation setters)
999 AliSysInfo::AddStamp("RunSimulation_Begin");
1001 AliSysInfo::AddStamp("RunSimulation_InitCDB");
1005 AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
1008 AliError("no gAlice object. Restart aliroot and try again.");
1011 if (gAlice->Modules()->GetEntries() > 0) {
1012 AliError("gAlice was already run. Restart aliroot and try again.");
1016 // Setup monitoring if requested
1017 gAlice->GetMCApp()->SetUseMonitoring(fUseMonitoring);
1019 AliInfo(Form("initializing gAlice with config file %s",
1020 fConfigFileName.Data()));
1023 // Initialize ALICE Simulation run
1028 // If requested set the mag. field from the GRP entry.
1029 // After this the field is loccked and cannot be changed by Config.C
1030 if (fUseMagFieldFromGRP) {
1032 grpM.ReadGRPEntry();
1034 AliInfo("Field is locked now. It cannot be changed in Config.C");
1038 TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
1039 gROOT->LoadMacro(fConfigFileName.Data());
1040 Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
1041 if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
1042 AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
1044 AliSysInfo::AddStamp("RunSimulation_Config");
1047 // If requested obtain the vertex position and vertex sigma_z from the CDB
1048 // This overwrites the settings from the Config.C
1049 if (fUseVertexFromCDB) {
1050 Double_t vtxPos[3] = {0., 0., 0.};
1051 Double_t vtxSig[3] = {0., 0., 0.};
1052 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1054 AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1056 if(vertex->GetXRes()>2.8) { // > pipe radius --> it's a dummy object, don't use it
1057 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1058 if (entry) vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1062 vertex->GetXYZ(vtxPos);
1063 vertex->GetSigmaXYZ(vtxSig);
1064 AliInfo("Overwriting Config.C vertex settings !");
1065 AliInfo(Form("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
1066 vtxPos[0], vtxPos[1], vtxPos[2], vtxSig[2]));
1068 AliGenerator *gen = gAlice->GetMCApp()->Generator();
1069 gen->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]); // vertex position
1070 gen->SetSigmaZ(vtxSig[2]);
1075 // If requested we take the SOR and EOR time-stamps from the GRP and use them
1076 // in order to generate the event time-stamps
1077 if (fUseTimeStampFromCDB) {
1079 grpM.ReadGRPEntry();
1080 const AliGRPObject *grpObj = grpM.GetGRPData();
1081 if (!grpObj || (grpObj->GetTimeEnd() <= grpObj->GetTimeStart())) {
1082 AliError("Missing GRP or bad SOR/EOR time-stamps! Switching off the time-stamp generation from GRP!");
1083 fTimeStart = fTimeEnd = 0;
1084 fUseTimeStampFromCDB = kFALSE;
1087 fTimeStart = grpObj->GetTimeStart();
1088 fTimeEnd = grpObj->GetTimeEnd();
1092 if(AliCDBManager::Instance()->GetRun() >= 0) {
1093 AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
1094 AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
1096 AliWarning("Run number not initialized!!");
1099 AliRunLoader::Instance()->CdGAFile();
1102 AliPDG::AddParticlesToPdgDataBase();
1104 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1105 AliSysInfo::AddStamp("RunSimulation_GetField");
1106 gAlice->GetMCApp()->Init();
1107 AliSysInfo::AddStamp("RunSimulation_InitMCApp");
1109 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1110 gAlice->InitLoaders();
1111 AliRunLoader::Instance()->MakeTree("E");
1112 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1113 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1114 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1116 // Save stuff at the beginning of the file to avoid file corruption
1117 AliRunLoader::Instance()->CdGAFile();
1119 gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1120 AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1121 //___________________________________________________________________________________________
1123 AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
1125 // Set run number in CDBManager
1126 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
1128 AliRunLoader* runLoader = AliRunLoader::Instance();
1130 AliError(Form("gAlice has no run loader object. "
1131 "Check your config file: %s", fConfigFileName.Data()));
1134 SetGAliceFile(runLoader->GetFileName());
1136 // Misalign geometry
1137 #if ROOT_VERSION_CODE < 331527
1138 AliGeomManager::SetGeometry(gGeoManager);
1140 // Check that the consistency of symbolic names for the activated subdetectors
1141 // in the geometry loaded by AliGeomManager
1142 TString detsToBeChecked = "";
1143 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1144 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1145 AliModule* det = (AliModule*) detArray->At(iDet);
1146 if (!det || !det->IsActive()) continue;
1147 detsToBeChecked += det->GetName();
1148 detsToBeChecked += " ";
1149 } // end loop over detectors
1150 if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
1151 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
1152 MisalignGeometry(runLoader);
1153 AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
1156 // AliRunLoader* runLoader = AliRunLoader::Instance();
1157 // if (!runLoader) {
1158 // AliError(Form("gAlice has no run loader object. "
1159 // "Check your config file: %s", fConfigFileName.Data()));
1162 // SetGAliceFile(runLoader->GetFileName());
1164 if (!gAlice->GetMCApp()->Generator()) {
1165 AliError(Form("gAlice has no generator object. "
1166 "Check your config file: %s", fConfigFileName.Data()));
1170 // Write GRP entry corresponding to the setting found in Cofig.C
1173 AliSysInfo::AddStamp("RunSimulation_WriteGRP");
1175 if (nEvents <= 0) nEvents = fNEvents;
1177 // get vertex from background file in case of merging
1178 if (fUseBkgrdVertex &&
1179 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
1180 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
1181 const char* fileName = ((TObjString*)
1182 (fBkgrdFileNames->At(0)))->GetName();
1183 AliInfo(Form("The vertex will be taken from the background "
1184 "file %s with nSignalPerBackground = %d",
1185 fileName, signalPerBkgrd));
1186 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
1187 gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
1190 if (!fRunSimulation) {
1191 gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
1194 // set the number of events per file for given detectors and data types
1195 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
1196 if (!fEventsPerFile[i]) continue;
1197 const char* detName = fEventsPerFile[i]->GetName();
1198 const char* typeName = fEventsPerFile[i]->GetTitle();
1199 TString loaderName(detName);
1200 loaderName += "Loader";
1201 AliLoader* loader = runLoader->GetLoader(loaderName);
1203 AliError(Form("RunSimulation no loader for %s found\n Number of events per file not set for %s %s",
1204 detName, typeName, detName));
1207 AliDataLoader* dataLoader =
1208 loader->GetDataLoader(typeName);
1210 AliError(Form("no data loader for %s found\n"
1211 "Number of events per file not set for %s %s",
1212 typeName, detName, typeName));
1215 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
1216 AliDebug(1, Form("number of events per file set to %d for %s %s",
1217 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
1220 AliInfo("running gAlice");
1221 AliSysInfo::AddStamp("Start_ProcessRun");
1223 // Create the Root Tree with one branch per detector
1224 //Hits moved to begin event -> now we are crating separate tree for each event
1225 gMC->ProcessRun(nEvents);
1227 // End of this run, close files
1228 if(nEvents>0) FinishRun();
1230 AliSysInfo::AddStamp("Stop_ProcessRun");
1236 //_____________________________________________________________________________
1237 Bool_t AliSimulation::RunGeneratorOnly()
1240 TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
1241 gROOT->LoadMacro(fConfigFileName.Data());
1242 Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
1243 if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
1244 AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
1247 // Setup the runloader and generator, check if everything is OK
1248 AliRunLoader* runLoader = AliRunLoader::Instance();
1249 AliGenerator* generator = gAlice->GetMCApp()->Generator();
1251 AliError(Form("gAlice has no run loader object. "
1252 "Check your config file: %s", fConfigFileName.Data()));
1256 AliError(Form("gAlice has no generator object. "
1257 "Check your config file: %s", fConfigFileName.Data()));
1261 runLoader->LoadKinematics("RECREATE");
1262 runLoader->MakeTree("E");
1264 // Create stack and header
1265 runLoader->MakeStack();
1266 AliStack* stack = runLoader->Stack();
1267 AliHeader* header = runLoader->GetHeader();
1269 // Intialize generator
1271 generator->SetStack(stack);
1273 // Run main generator loop
1275 for (Int_t iev=0; iev<fNEvents; iev++)
1278 header->Reset(0,iev);
1279 runLoader->SetEventNumber(iev);
1281 runLoader->MakeTree("K");
1284 generator->Generate();
1287 header->SetNprimary(stack->GetNprimary());
1288 header->SetNtrack(stack->GetNtrack());
1289 stack->FinishEvent();
1290 header->SetStack(stack);
1291 runLoader->TreeE()->Fill();
1292 runLoader->WriteKinematics("OVERWRITE");
1296 generator->FinishRun();
1298 runLoader->WriteHeader("OVERWRITE");
1305 //_____________________________________________________________________________
1306 Bool_t AliSimulation::RunSDigitization(const char* detectors)
1308 // run the digitization and produce summable digits
1309 static Int_t eventNr=0;
1310 AliCodeTimerAuto("",0) ;
1312 // initialize CDB storage, run number, set CDB lock
1314 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1317 AliRunLoader* runLoader = LoadRun();
1318 if (!runLoader) return kFALSE;
1320 TString detStr = detectors;
1321 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1322 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1323 AliModule* det = (AliModule*) detArray->At(iDet);
1324 if (!det || !det->IsActive()) continue;
1325 if (IsSelected(det->GetName(), detStr)) {
1326 AliInfo(Form("creating summable digits for %s", det->GetName()));
1327 AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
1328 det->Hits2SDigits();
1329 AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
1330 AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
1334 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1335 AliError(Form("the following detectors were not found: %s",
1337 if (fStopOnError) return kFALSE;
1346 //_____________________________________________________________________________
1347 Bool_t AliSimulation::RunDigitization(const char* detectors,
1348 const char* excludeDetectors)
1350 // run the digitization and produce digits from sdigits
1352 AliCodeTimerAuto("",0)
1354 // initialize CDB storage, run number, set CDB lock
1356 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1359 delete AliRunLoader::Instance();
1364 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1365 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1366 AliDigitizationInput digInp(nStreams, signalPerBkgrd);
1367 // digInp.SetEmbeddingFlag(fEmbeddingFlag);
1368 digInp.SetRegionOfInterest(fRegionOfInterest);
1369 digInp.SetInputStream(0, fGAliceFileName.Data());
1370 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1371 const char* fileName = ((TObjString*)(fBkgrdFileNames->At(iStream-1)))->GetName();
1372 digInp.SetInputStream(iStream, fileName);
1375 detArr.SetOwner(kTRUE);
1376 TString detStr = detectors;
1377 TString detExcl = excludeDetectors;
1378 if (!static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice()) {
1379 AliError("Error occured while getting gAlice from Input 0");
1382 AliRunLoader* runLoader = AliRunLoader::GetRunLoader(digInp.GetInputStream(0)->GetFolderName());
1383 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1384 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1385 AliModule* det = (AliModule*) detArray->At(iDet);
1386 if (!det || !det->IsActive()) continue;
1387 if (!IsSelected(det->GetName(), detStr) || IsSelected(det->GetName(), detExcl)) continue;
1388 AliDigitizer* digitizer = det->CreateDigitizer(&digInp);
1389 if (!digitizer || !digitizer->Init()) {
1390 AliError(Form("no digitizer for %s", det->GetName()));
1391 if (fStopOnError) return kFALSE;
1394 detArr.AddLast(digitizer);
1395 AliInfo(Form("Created digitizer from SDigits -> Digits for %s", det->GetName()));
1398 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1399 AliError(Form("the following detectors were not found: %s", detStr.Data()));
1400 if (fStopOnError) return kFALSE;
1403 Int_t ndigs = detArr.GetEntriesFast();
1404 Int_t eventsCreated = 0;
1405 AliRunLoader* outRl = digInp.GetOutRunLoader();
1406 while ((eventsCreated++ < fNEvents) || (fNEvents < 0)) {
1407 if (!digInp.ConnectInputTrees()) break;
1408 digInp.InitEvent(); //this must be after call of Connect Input tress.
1409 if (outRl) outRl->SetEventNumber(eventsCreated-1);
1410 static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice(); // use gAlice of the first input stream
1411 for (int id=0;id<ndigs;id++) ((AliDigitizer*)detArr[id])->Digitize("");
1412 digInp.FinishEvent();
1414 digInp.FinishGlobal();
1419 //_____________________________________________________________________________
1420 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1422 // run the digitization and produce digits from hits
1424 AliCodeTimerAuto("",0)
1426 // initialize CDB storage, run number, set CDB lock
1428 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1431 AliRunLoader* runLoader = LoadRun("READ");
1432 if (!runLoader) return kFALSE;
1434 TString detStr = detectors;
1435 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1436 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1437 AliModule* det = (AliModule*) detArray->At(iDet);
1438 if (!det || !det->IsActive()) continue;
1439 if (IsSelected(det->GetName(), detStr)) {
1440 AliInfo(Form("creating digits from hits for %s", det->GetName()));
1445 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1446 AliError(Form("the following detectors were not found: %s",
1448 if (fStopOnError) return kFALSE;
1454 //_____________________________________________________________________________
1455 Bool_t AliSimulation::WriteRawData(const char* detectors,
1456 const char* fileName,
1457 Bool_t deleteIntermediateFiles,
1460 // convert the digits to raw data
1461 // First DDL raw data files for the given detectors are created.
1462 // If a file name is given, the DDL files are then converted to a DATE file.
1463 // If deleteIntermediateFiles is true, the DDL raw files are deleted
1465 // If the file name has the extension ".root", the DATE file is converted
1467 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1468 // 'selrawdata' flag can be used to enable writing of detectors raw data
1469 // accoring to the trigger cluster.
1471 AliCodeTimerAuto("",0)
1472 AliSysInfo::AddStamp("WriteRawData_Start");
1474 TString detStr = detectors;
1475 if (!WriteRawFiles(detStr.Data())) {
1476 if (fStopOnError) return kFALSE;
1478 AliSysInfo::AddStamp("WriteRawFiles");
1480 // run HLT simulation on simulated DDL raw files
1481 // and produce HLT ddl raw files to be included in date/root file
1482 // bugfix 2009-06-26: the decision whether to write HLT raw data
1483 // is taken in RunHLT. Here HLT always needs to be run in order to
1484 // create HLT digits, unless its switched off. This is due to the
1485 // special placement of the HLT between the generation of DDL files
1486 // and conversion to DATE/Root file.
1487 detStr.ReplaceAll("HLT", "");
1488 if (!fRunHLT.IsNull()) {
1490 if (fStopOnError) return kFALSE;
1493 AliSysInfo::AddStamp("WriteRawData_RunHLT");
1495 TString dateFileName(fileName);
1496 if (!dateFileName.IsNull()) {
1497 Bool_t rootOutput = dateFileName.EndsWith(".root");
1498 if (rootOutput) dateFileName += ".date";
1499 TString selDateFileName;
1501 selDateFileName = "selected.";
1502 selDateFileName+= dateFileName;
1504 if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1505 if (fStopOnError) return kFALSE;
1507 AliSysInfo::AddStamp("ConvertRawFilesToDate");
1508 if (deleteIntermediateFiles) {
1509 AliRunLoader* runLoader = LoadRun("READ");
1510 if (runLoader) for (Int_t iEvent = 0;
1511 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1513 snprintf(command, 256, "rm -r raw%d", iEvent);
1514 gSystem->Exec(command);
1520 if (!ConvertDateToRoot(dateFileName, fileName)) {
1521 if (fStopOnError) return kFALSE;
1523 AliSysInfo::AddStamp("ConvertDateToRoot");
1524 if (deleteIntermediateFiles) {
1525 gSystem->Unlink(dateFileName);
1528 TString selFileName = "selected.";
1529 selFileName += fileName;
1530 if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1531 if (fStopOnError) return kFALSE;
1533 if (deleteIntermediateFiles) {
1534 gSystem->Unlink(selDateFileName);
1543 //_____________________________________________________________________________
1544 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1546 // convert the digits to raw data DDL files
1548 AliCodeTimerAuto("",0)
1550 AliRunLoader* runLoader = LoadRun("READ");
1551 if (!runLoader) return kFALSE;
1553 // write raw data to DDL files
1554 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1555 AliInfo(Form("processing event %d", iEvent));
1556 runLoader->GetEvent(iEvent);
1557 TString baseDir = gSystem->WorkingDirectory();
1559 snprintf(dirName, 256, "raw%d", iEvent);
1560 gSystem->MakeDirectory(dirName);
1561 if (!gSystem->ChangeDirectory(dirName)) {
1562 AliError(Form("couldn't change to directory %s", dirName));
1563 if (fStopOnError) return kFALSE; else continue;
1566 ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1569 TString detStr = detectors;
1570 if (IsSelected("HLT", detStr)) {
1571 // Do nothing. "HLT" will be removed from detStr and HLT raw
1572 // data files are generated in RunHLT.
1575 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1576 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1577 AliModule* det = (AliModule*) detArray->At(iDet);
1578 if (!det || !det->IsActive()) continue;
1579 if (IsSelected(det->GetName(), detStr)) {
1580 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1585 if (!WriteTriggerRawData())
1586 if (fStopOnError) return kFALSE;
1588 gSystem->ChangeDirectory(baseDir);
1589 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1590 AliError(Form("the following detectors were not found: %s",
1592 if (fStopOnError) return kFALSE;
1601 //_____________________________________________________________________________
1602 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1603 const char* selDateFileName)
1605 // convert raw data DDL files to a DATE file with the program "dateStream"
1606 // The second argument is not empty when the user decides to write
1607 // the detectors raw data according to the trigger cluster.
1609 AliCodeTimerAuto("",0)
1611 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1613 AliError("the program dateStream was not found");
1614 if (fStopOnError) return kFALSE;
1619 AliRunLoader* runLoader = LoadRun("READ");
1620 if (!runLoader) return kFALSE;
1622 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1623 Bool_t selrawdata = kFALSE;
1624 if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1627 // Note the option -s. It is used in order to avoid
1628 // the generation of SOR/EOR events.
1629 snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d",
1630 dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1631 FILE* pipe = gSystem->OpenPipe(command, "w");
1634 AliError(Form("Cannot execute command: %s",command));
1638 Int_t selEvents = 0;
1639 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1641 UInt_t detectorPattern = 0;
1642 runLoader->GetEvent(iEvent);
1643 if (!runLoader->LoadTrigger()) {
1644 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1645 detectorPattern = aCTP->GetClusterMask();
1646 // Check if the event was triggered by CTP
1648 if (aCTP->GetClassMask()) selEvents++;
1652 AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1654 AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1655 selrawdata = kFALSE;
1659 fprintf(pipe, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
1663 // loop over detectors and DDLs
1664 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1665 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1667 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1668 Int_t ldcID = Int_t(ldc + 0.0001);
1669 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1671 char rawFileName[256];
1672 snprintf(rawFileName, 256, "raw%d/%s",
1673 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1675 // check existence and size of raw data file
1676 FILE* file = fopen(rawFileName, "rb");
1677 if (!file) continue;
1678 fseek(file, 0, SEEK_END);
1679 unsigned long size = ftell(file);
1681 if (!size) continue;
1683 if (ldcID != prevLDC) {
1684 fprintf(pipe, " LDC Id %d\n", ldcID);
1687 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1692 Int_t result = gSystem->ClosePipe(pipe);
1694 if (!(selrawdata && selEvents > 0)) {
1696 return (result == 0);
1699 AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1701 snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d",
1702 selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1703 FILE* pipe2 = gSystem->OpenPipe(command, "w");
1705 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1707 // Get the trigger decision and cluster
1708 UInt_t detectorPattern = 0;
1710 runLoader->GetEvent(iEvent);
1711 if (!runLoader->LoadTrigger()) {
1712 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1713 if (aCTP->GetClassMask() == 0) continue;
1714 detectorPattern = aCTP->GetClusterMask();
1715 detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1716 AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1719 fprintf(pipe2, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
1723 // loop over detectors and DDLs
1724 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1725 // Write only raw data from detectors that
1726 // are contained in the trigger cluster(s)
1727 if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1729 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1731 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1732 Int_t ldcID = Int_t(ldc + 0.0001);
1733 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1735 char rawFileName[256];
1736 snprintf(rawFileName, 256, "raw%d/%s",
1737 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1739 // check existence and size of raw data file
1740 FILE* file = fopen(rawFileName, "rb");
1741 if (!file) continue;
1742 fseek(file, 0, SEEK_END);
1743 unsigned long size = ftell(file);
1745 if (!size) continue;
1747 if (ldcID != prevLDC) {
1748 fprintf(pipe2, " LDC Id %d\n", ldcID);
1751 fprintf(pipe2, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1756 Int_t result2 = gSystem->ClosePipe(pipe2);
1759 return ((result == 0) && (result2 == 0));
1762 //_____________________________________________________________________________
1763 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1764 const char* rootFileName)
1766 // convert a DATE file to a root file with the program "alimdc"
1769 const Int_t kDBSize = 2000000000;
1770 const Int_t kTagDBSize = 1000000000;
1771 const Bool_t kFilter = kFALSE;
1772 const Int_t kCompression = 1;
1774 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1776 AliError("the program alimdc was not found");
1777 if (fStopOnError) return kFALSE;
1782 AliInfo(Form("converting DATE file %s to root file %s",
1783 dateFileName, rootFileName));
1785 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1786 const char* tagDBFS = "/tmp/mdc1/tags";
1788 // User defined file system locations
1789 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1790 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1791 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1792 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1793 if (gSystem->Getenv("ALIMDC_TAGDB"))
1794 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1796 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1797 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1798 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1800 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1801 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1802 gSystem->Exec(Form("mkdir %s",tagDBFS));
1804 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1805 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1806 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1808 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1809 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1810 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1812 return (result == 0);
1816 //_____________________________________________________________________________
1817 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1819 // delete existing run loaders, open a new one and load gAlice
1821 delete AliRunLoader::Instance();
1822 AliRunLoader* runLoader =
1823 AliRunLoader::Open(fGAliceFileName.Data(),
1824 AliConfig::GetDefaultEventFolderName(), mode);
1826 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1829 runLoader->LoadgAlice();
1830 runLoader->LoadHeader();
1831 gAlice = runLoader->GetAliRun();
1833 AliError(Form("no gAlice object found in file %s",
1834 fGAliceFileName.Data()));
1840 //_____________________________________________________________________________
1841 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1843 // get or calculate the number of signal events per background event
1845 if (!fBkgrdFileNames) return 1;
1846 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1847 if (nBkgrdFiles == 0) return 1;
1849 // get the number of signal events
1851 AliRunLoader* runLoader =
1852 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1853 if (!runLoader) return 1;
1855 nEvents = runLoader->GetNumberOfEvents();
1860 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1861 // get the number of background events
1862 const char* fileName = ((TObjString*)
1863 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1864 AliRunLoader* runLoader =
1865 AliRunLoader::Open(fileName, "BKGRD");
1866 if (!runLoader) continue;
1867 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1870 // get or calculate the number of signal per background events
1871 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1872 if (nSignalPerBkgrd <= 0) {
1873 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1874 } else if (result && (result != nSignalPerBkgrd)) {
1875 AliInfo(Form("the number of signal events per background event "
1876 "will be changed from %d to %d for stream %d",
1877 nSignalPerBkgrd, result, iBkgrdFile+1));
1878 nSignalPerBkgrd = result;
1881 if (!result) result = nSignalPerBkgrd;
1882 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1883 AliWarning(Form("not enough background events (%d) for %d signal events "
1884 "using %d signal per background events for stream %d",
1885 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1892 //_____________________________________________________________________________
1893 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1895 // check whether detName is contained in detectors
1896 // if yes, it is removed from detectors
1898 // check if all detectors are selected
1899 if ((detectors.CompareTo("ALL") == 0) ||
1900 detectors.BeginsWith("ALL ") ||
1901 detectors.EndsWith(" ALL") ||
1902 detectors.Contains(" ALL ")) {
1907 // search for the given detector
1908 Bool_t result = kFALSE;
1909 if ((detectors.CompareTo(detName) == 0) ||
1910 detectors.BeginsWith(detName+" ") ||
1911 detectors.EndsWith(" "+detName) ||
1912 detectors.Contains(" "+detName+" ")) {
1913 detectors.ReplaceAll(detName, "");
1917 // clean up the detectors string
1918 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1919 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1920 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1925 //_____________________________________________________________________________
1926 Int_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName, Int_t N, Int_t nSkip)
1929 // Steering routine to convert raw data in directory rawDirectory/ to fake SDigits.
1930 // These can be used for embedding of MC tracks into RAW data using the standard
1931 // merging procedure.
1933 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1936 AliError("no gAlice object. Restart aliroot and try again.");
1939 if (gAlice->Modules()->GetEntries() > 0) {
1940 AliError("gAlice was already run. Restart aliroot and try again.");
1944 AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1948 gROOT->LoadMacro(fConfigFileName.Data());
1949 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1951 if(AliCDBManager::Instance()->GetRun() >= 0) {
1952 SetRunNumber(AliCDBManager::Instance()->GetRun());
1954 AliWarning("Run number not initialized!!");
1957 AliRunLoader::Instance()->CdGAFile();
1959 AliPDG::AddParticlesToPdgDataBase();
1961 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1963 gAlice->GetMCApp()->Init();
1965 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1966 gAlice->InitLoaders();
1967 AliRunLoader::Instance()->MakeTree("E");
1968 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1969 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1970 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1973 // Save stuff at the beginning of the file to avoid file corruption
1974 AliRunLoader::Instance()->CdGAFile();
1979 //AliCDBManager* man = AliCDBManager::Instance();
1980 //man->SetRun(0); // Should this come from rawdata header ?
1984 // Get the runloader
1985 AliRunLoader* runLoader = AliRunLoader::Instance();
1987 // Open esd file if available
1990 AliESDEvent* esd = 0;
1991 if (esdFileName && (strlen(esdFileName)>0)) {
1992 esdFile = TFile::Open(esdFileName);
1994 esd = new AliESDEvent();
1995 esdFile->GetObject("esdTree", treeESD);
1997 esd->ReadFromTree(treeESD);
1999 AliInfo(Form("Asking to skip first %d ESDs events",nSkip));
2008 // Create the RawReader
2009 TString fileName(rawDirectory);
2010 AliRawReader* rawReader = AliRawReader::Create(fileName.Data());
2011 if (!rawReader) return (kFALSE);
2013 // if (!fEquipIdMap.IsNull() && fRawReader)
2014 // fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
2016 // Get list of detectors
2017 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
2020 AliHeader* header = runLoader->GetHeader();
2024 if (!(rawReader->NextEvent())) break;
2025 runLoader->SetEventNumber(nev);
2026 runLoader->GetHeader()->Reset(rawReader->GetRunNumber(),
2028 runLoader->GetEvent(nev);
2029 AliInfo(Form("We are at event %d",nev));
2032 TString detStr = fMakeSDigits;
2033 for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
2034 AliModule* det = (AliModule*) detArray->At(iDet);
2035 if (!det || !det->IsActive()) continue;
2036 if (IsSelected(det->GetName(), detStr)) {
2037 AliInfo(Form("Calling Raw2SDigits for %s", det->GetName()));
2038 det->Raw2SDigits(rawReader);
2045 // If ESD information available obtain reconstructed vertex and store in header.
2047 AliInfo(Form("Selected event %d correspond to event %d in raw and to %d in esd",nev,rawReader->GetEventIndex(),nSkip+rawReader->GetEventIndex()));
2048 treeESD->GetEvent(nSkip+rawReader->GetEventIndex());
2049 const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
2050 Double_t position[3];
2051 esdVertex->GetXYZ(position);
2052 AliGenEventHeader* mcHeader = new AliGenEventHeader("ESD");
2055 for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
2056 mcHeader->SetPrimaryVertex(mcV);
2057 header->Reset(0,nev);
2058 header->SetGenEventHeader(mcHeader);
2059 AliInfo(Form("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]));
2063 runLoader->TreeE()->Fill();
2064 AliInfo(Form("Finished event %d",nev));
2073 runLoader->CdGAFile();
2074 runLoader->WriteHeader("OVERWRITE");
2075 runLoader->WriteRunLoader();
2080 //_____________________________________________________________________________
2081 void AliSimulation::FinishRun()
2084 // Called at the end of the run.
2089 AliDebug(1, "Finish Lego");
2090 AliRunLoader::Instance()->CdGAFile();
2094 // Clean detector information
2095 TIter next(gAlice->Modules());
2096 AliModule *detector;
2097 while((detector = dynamic_cast<AliModule*>(next()))) {
2098 AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
2099 detector->FinishRun();
2102 AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
2103 AliRunLoader::Instance()->WriteHeader("OVERWRITE");
2105 // Write AliRun info and all detectors parameters
2106 AliRunLoader::Instance()->CdGAFile();
2107 gAlice->Write(0,TObject::kOverwrite);//write AliRun
2108 AliRunLoader::Instance()->Write(0,TObject::kOverwrite);//write RunLoader itself
2112 if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();
2113 AliRunLoader::Instance()->Synchronize();
2116 //_____________________________________________________________________________
2117 Int_t AliSimulation::GetDetIndex(const char* detector)
2119 // return the detector index corresponding to detector
2121 for (index = 0; index < fgkNDetectors ; index++) {
2122 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2128 //_____________________________________________________________________________
2129 Bool_t AliSimulation::CreateHLT()
2131 // Init the HLT simulation.
2132 // The function loads the library and creates the instance of AliHLTSimulation.
2133 // the main reason for the decoupled creation is to set the transient OCDB
2134 // objects before the OCDB is locked
2136 // load the library dynamically
2137 gSystem->Load(ALIHLTSIMULATION_LIBRARY);
2139 // check for the library version
2140 AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
2142 AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
2145 if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
2146 AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
2149 // print compile info
2150 typedef void (*CompileInfo)( const char*& date, const char*& time);
2151 CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
2153 const char* date="";
2154 const char* time="";
2155 (*fctInfo)(date, time);
2156 if (!date) date="unknown";
2157 if (!time) time="unknown";
2158 AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
2160 AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
2163 // create instance of the HLT simulation
2164 AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
2165 if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
2166 AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
2170 TString specObjects;
2171 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
2172 if (specObjects.Length()>0) specObjects+=" ";
2173 specObjects+=fSpecCDBUri[i]->GetName();
2176 AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
2177 if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
2178 AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
2184 //_____________________________________________________________________________
2185 Bool_t AliSimulation::RunHLT()
2187 // Run the HLT simulation
2188 // HLT simulation is implemented in HLT/sim/AliHLTSimulation
2189 // Disabled if fRunHLT is empty, default vaule is "default".
2190 // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
2191 // The default simulation depends on the HLT component libraries and their
2192 // corresponding agents which define components and chains to run. See
2193 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
2194 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
2196 // The libraries to be loaded can be specified as an option.
2198 // AliSimulation sim;
2199 // sim.SetRunHLT("libAliHLTSample.so");
2201 // will only load <tt>libAliHLTSample.so</tt>
2203 // Other available options:
2204 // \li loglevel=<i>level</i> <br>
2205 // logging level for this processing
2207 // disable redirection of log messages to AliLog class
2208 // \li config=<i>macro</i>
2209 // configuration macro
2210 // \li chains=<i>configuration</i>
2211 // comma separated list of configurations to be run during simulation
2212 // \li rawfile=<i>file</i>
2213 // source for the RawReader to be created, the default is <i>./</i> if
2214 // raw data is simulated
2218 if (!fpHLT && !CreateHLT()) {
2221 AliHLTSimulation* pHLT=fpHLT;
2223 AliRunLoader* pRunLoader = LoadRun("READ");
2224 if (!pRunLoader) return kFALSE;
2226 // initialize CDB storage, run number, set CDB lock
2227 // thats for the case of running HLT simulation without all the other steps
2228 // multiple calls are handled by the function, so we can just call
2230 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
2233 // init the HLT simulation
2235 if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2236 TString detStr = fWriteRawData;
2237 if (!IsSelected("HLT", detStr)) {
2238 options+=" writerawfiles=";
2240 options+=" writerawfiles=HLT";
2243 if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2244 // as a matter of fact, HLT will run reconstruction and needs the RawReader
2245 // in order to get detector data. By default, RawReaderFile is used to read
2246 // the already simulated ddl files. Date and Root files from the raw data
2247 // are generated after the HLT simulation.
2248 options+=" rawfile=./";
2251 AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2252 if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2253 AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2255 // run the HLT simulation
2256 AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2257 if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2258 AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2262 // delete the instance
2263 AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2264 if (fctDelete==NULL || fctDelete(pHLT)<0) {
2265 AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2269 return iResult>=0?kTRUE:kFALSE;
2272 //_____________________________________________________________________________
2273 Bool_t AliSimulation::RunQA()
2275 // run the QA on summable hits, digits or digits
2277 //if(!gAlice) return kFALSE;
2278 AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2280 TString detectorsw("") ;
2282 AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2283 detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ;
2284 if ( detectorsw.IsNull() )
2289 //_____________________________________________________________________________
2290 Bool_t AliSimulation::SetRunQA(TString detAndAction)
2292 // Allows to run QA for a selected set of detectors
2293 // and a selected set of tasks among HITS, SDIGITS and DIGITS
2294 // all selected detectors run the same selected tasks
2296 if (!detAndAction.Contains(":")) {
2297 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2301 Int_t colon = detAndAction.Index(":") ;
2302 fQADetectors = detAndAction(0, colon) ;
2303 if (fQADetectors.Contains("ALL") ){
2304 TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2305 Int_t minus = fQADetectors.Last('-') ;
2306 TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2307 TString toRemove("") ;
2308 while (minus >= 0) {
2309 toRemove = fQADetectors(minus+1, fQADetectors.Length()) ;
2310 toRemove = toRemove.Strip() ;
2311 toKeep.ReplaceAll(toRemove, "") ;
2312 fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ;
2313 minus = fQADetectors.Last('-') ;
2315 fQADetectors = toKeep ;
2317 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2318 if (fQATasks.Contains("ALL") ) {
2319 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
2321 fQATasks.ToUpper() ;
2323 if ( fQATasks.Contains("HIT") )
2324 tempo = Form("%d ", AliQAv1::kHITS) ;
2325 if ( fQATasks.Contains("SDIGIT") )
2326 tempo += Form("%d ", AliQAv1::kSDIGITS) ;
2327 if ( fQATasks.Contains("DIGIT") )
2328 tempo += Form("%d ", AliQAv1::kDIGITS) ;
2330 if (fQATasks.IsNull()) {
2331 AliInfo("No QA requested\n") ;
2336 TString tempo(fQATasks) ;
2337 tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS)) ;
2338 tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;
2339 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;
2340 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2342 AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ;
2343 AliQAManager::QAManager()->SetTasks(fQATasks) ;
2344 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++)
2345 AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2350 //_____________________________________________________________________________
2351 void AliSimulation::ProcessEnvironmentVars()
2353 // Extract run number and random generator seed from env variables
2355 AliInfo("Processing environment variables");
2357 // Random Number seed
2359 // first check that seed is not already set
2361 if (gSystem->Getenv("CONFIG_SEED")) {
2362 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2365 if (gSystem->Getenv("CONFIG_SEED")) {
2366 AliInfo(Form("Seed for random number generation already set (%d)"
2367 ": CONFIG_SEED variable ignored!", fSeed));
2371 AliInfo(Form("Seed for random number generation = %d ", fSeed));
2375 // first check that run number is not already set
2377 if (gSystem->Getenv("DC_RUN")) {
2378 fRun = atoi(gSystem->Getenv("DC_RUN"));
2381 if (gSystem->Getenv("DC_RUN")) {
2382 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2386 AliInfo(Form("Run number = %d", fRun));
2389 //---------------------------------------------------------------------
2390 void AliSimulation::WriteGRPEntry()
2392 // Get the necessary information from galice (generator, trigger etc) and
2393 // write a GRP entry corresponding to the settings in the Config.C used
2394 // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2397 AliInfo("Writing global run parameters entry into the OCDB");
2399 AliGRPObject* grpObj = new AliGRPObject();
2401 grpObj->SetRunType("PHYSICS");
2402 grpObj->SetTimeStart(fTimeStart);
2403 grpObj->SetTimeEnd(fTimeEnd);
2404 grpObj->SetBeamEnergyIsSqrtSHalfGeV(); // new format of GRP: store sqrt(s)/2 in GeV
2406 const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2412 gen->GetProjectile(projectile,a,z);
2414 gen->GetTarget(target,a,z);
2415 TString beamType = projectile + "-" + target;
2416 beamType.ReplaceAll(" ","");
2417 if (!beamType.CompareTo("-")) {
2418 grpObj->SetBeamType("UNKNOWN");
2419 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
2422 grpObj->SetBeamType(beamType);
2424 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 * a / z);
2426 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 );
2428 // Heavy ion run, the event specie is set to kHighMult
2429 fEventSpecie = AliRecoParam::kHighMult;
2430 if ((strcmp(beamType,"p-p") == 0) ||
2431 (strcmp(beamType,"p-") == 0) ||
2432 (strcmp(beamType,"-p") == 0) ||
2433 (strcmp(beamType,"P-P") == 0) ||
2434 (strcmp(beamType,"P-") == 0) ||
2435 (strcmp(beamType,"-P") == 0)) {
2436 // Proton run, the event specie is set to kLowMult
2437 fEventSpecie = AliRecoParam::kLowMult;
2441 AliWarning("Unknown beam type and energy! Setting energy to 0");
2442 grpObj->SetBeamEnergy(0);
2443 grpObj->SetBeamType("UNKNOWN");
2446 UInt_t detectorPattern = 0;
2448 TObjArray *detArray = gAlice->Detectors();
2449 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2450 if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2451 AliDebug(1, Form("Detector #%d found: %s", iDet, AliDAQ::OfflineModuleName(iDet)));
2452 detectorPattern |= (1 << iDet);
2457 if (!fTriggerConfig.IsNull())
2458 detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2461 if (!fRunHLT.IsNull())
2462 detectorPattern |= (1 << AliDAQ::kHLTId);
2464 grpObj->SetNumberOfDetectors((Char_t)nDets);
2465 grpObj->SetDetectorMask((Int_t)detectorPattern);
2466 grpObj->SetLHCPeriod("LHC08c");
2467 grpObj->SetLHCState("STABLE_BEAMS");
2469 AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2470 Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2472 Float_t factorSol = field ? field->GetFactorSol() : 0;
2473 Float_t currentSol = TMath::Abs(factorSol)>1E-6 ?
2474 TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
2476 Float_t factorDip = field ? field->GetFactorDip() : 0;
2477 Float_t currentDip = 6000.*TMath::Abs(factorDip);
2479 grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
2480 grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);
2481 grpObj->SetL3Polarity(factorSol>0 ? 0:1);
2482 grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
2483 if (field) grpObj->SetUniformBMap(field->IsUniform()); // for special MC with k5kGUniform map
2484 grpObj->SetPolarityConventionLHC(); // LHC convention +/+ current -> -/- field main components
2486 grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2488 //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2490 // Now store the entry in OCDB
2491 AliCDBManager* man = AliCDBManager::Instance();
2493 man->SetLock(0, fKey);
2495 AliCDBStorage* sto = man->GetStorage(fGRPWriteLocation.Data());
2498 AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2499 AliCDBMetaData *metadata= new AliCDBMetaData();
2501 metadata->SetResponsible("alice-off@cern.ch");
2502 metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2504 sto->Put(grpObj,id,metadata);
2505 man->SetLock(1, fKey);
2508 //_____________________________________________________________________________
2509 time_t AliSimulation::GenerateTimeStamp() const
2511 // Generate event time-stamp according to
2512 // SOR/EOR time from GRP
2513 if (fUseTimeStampFromCDB)
2514 return fTimeStart + gRandom->Integer(fTimeEnd-fTimeStart);
2519 //_____________________________________________________________________________
2520 void AliSimulation::StoreUsedCDBMaps() const
2522 // write in galice.root maps with used CDB paths
2524 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
2525 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
2527 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
2528 cdbMapCopy->SetOwner(1);
2529 // cdbMapCopy->SetName("cdbMap");
2530 TIter iter(cdbMap->GetTable());
2533 while((pair = dynamic_cast<TPair*> (iter.Next()))){
2534 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
2535 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
2536 if (keyStr && valStr)
2537 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
2540 TList *cdbListCopy = new TList();
2541 cdbListCopy->SetOwner(1);
2542 // cdbListCopy->SetName("cdbList");
2544 TIter iter2(cdbList);
2547 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
2548 cdbListCopy->Add(new TObjString(id->ToString().Data()));
2551 AliRunLoader::Instance()->CdGAFile();
2552 gDirectory->WriteObject(cdbMapCopy,"cdbMap","kSingleKey");
2553 gDirectory->WriteObject(cdbListCopy,"cdbList","kSingleKey");
2555 AliInfo(Form("Stored used OCDB entries as TMap %s and TList %s in %s",
2556 cdbMapCopy->GetName(),
2557 cdbListCopy->GetName(),
2558 fGAliceFileName.Data()));