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 // The trigger configuration is set by the method SetTriggerConfig(X) //
107 // X can take three kinds of values //
109 // - The exact string "none" - case insensitive. In this case, not trigger //
110 // information is generated from the digits. //
111 // - The empty string or "ocdb" - case insensitive. In this case the //
112 // trigger configuration is read from OCDB //
113 // - Some string - say "p-p" - in which case the configuration is read from //
114 // fixed files in $ALICE_ROOT/GRP/CTP/ - say $ALICE_ROOT/GRP/CTP/p-p.cfg //
116 // Default is to read from OCDB. //
118 ///////////////////////////////////////////////////////////////////////////////
121 #include <TGeoGlobalMagField.h>
122 #include <TGeoManager.h>
123 #include <TObjString.h>
126 #include <TVirtualMC.h>
127 #include <TVirtualMCApplication.h>
129 #include <TInterpreter.h>
131 #include "AliAlignObj.h"
132 #include "AliCDBEntry.h"
133 #include "AliCDBManager.h"
134 #include "AliGRPManager.h"
135 #include "AliCDBStorage.h"
136 #include "AliCTPRawData.h"
137 #include "AliCentralTrigger.h"
138 #include "AliCentralTrigger.h"
139 #include "AliCodeTimer.h"
141 #include "AliDigitizer.h"
142 #include "AliESDEvent.h"
143 #include "AliGRPObject.h"
144 #include "AliGenEventHeader.h"
145 #include "AliGenerator.h"
146 #include "AliGeomManager.h"
147 #include "AliHLTSimulation.h"
148 #include "AliHeader.h"
150 #include "AliLegoGenerator.h"
154 #include "AliModule.h"
156 #include "AliRawReaderDate.h"
157 #include "AliRawReaderFile.h"
158 #include "AliRawReaderRoot.h"
160 #include "AliDigitizationInput.h"
161 #include "AliRunLoader.h"
162 #include "AliStack.h"
163 #include "AliSimulation.h"
164 #include "AliSysInfo.h"
165 #include "AliVertexGenFile.h"
168 ClassImp(AliSimulation)
170 AliSimulation *AliSimulation::fgInstance = 0;
171 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD",
172 "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE","AD",
175 //_____________________________________________________________________________
176 AliSimulation::AliSimulation(const char* configFileName,
177 const char* name, const char* title) :
180 fRunGeneratorOnly(kFALSE),
181 fRunGeneration(kTRUE),
182 fRunSimulation(kTRUE),
183 fLoadAlignFromCDB(kTRUE),
184 fLoadAlObjsListOfDets("ALL"),
188 fMakeDigitsFromHits(""),
190 fRawDataFileName(""),
191 fDeleteIntermediateFiles(kFALSE),
192 fWriteSelRawData(kFALSE),
193 fStopOnError(kFALSE),
194 fUseMonitoring(kFALSE),
196 fConfigFileName(configFileName),
197 fGAliceFileName("galice.root"),
199 fBkgrdFileNames(NULL),
200 fAlignObjArray(NULL),
201 fUseBkgrdVertex(kTRUE),
202 fRegionOfInterest(kFALSE),
208 fInitCDBCalled(kFALSE),
209 fInitRunNumberCalled(kFALSE),
210 fSetRunNumberFromDataCalled(kFALSE),
211 fEmbeddingFlag(kFALSE),
214 fUseVertexFromCDB(0),
215 fUseMagFieldFromGRP(0),
216 fGRPWriteLocation(Form("local://%s", gSystem->pwd())),
217 fUseTimeStampFromCDB(0),
223 fEventSpecie(AliRecoParam::kDefault),
224 fWriteQAExpertData(kTRUE),
228 fWriteGRPEntry(kTRUE)
230 // create simulation object with default parameters
232 SetGAliceFile("galice.root");
235 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
236 qam->SetActiveDetectors(fQADetectors) ;
237 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
238 qam->SetTasks(fQATasks) ;
241 //_____________________________________________________________________________
242 AliSimulation::~AliSimulation()
246 fEventsPerFile.Delete();
247 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
248 // delete fAlignObjArray; fAlignObjArray=0;
250 if (fBkgrdFileNames) {
251 fBkgrdFileNames->Delete();
252 delete fBkgrdFileNames;
255 fSpecCDBUri.Delete();
256 if (fgInstance==this) fgInstance = 0;
258 AliQAManager::QAManager()->ShowQA() ;
259 AliQAManager::Destroy() ;
260 AliCodeTimer::Instance()->Print();
264 //_____________________________________________________________________________
265 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
267 // set the number of events for one run
272 //_____________________________________________________________________________
273 void AliSimulation::InitQA()
275 // activate a default CDB storage
276 // First check if we have any CDB storage set, because it is used
277 // to retrieve the calibration and alignment constants
279 if (fInitCDBCalled) return;
280 fInitCDBCalled = kTRUE;
282 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
283 qam->SetActiveDetectors(fQADetectors) ;
284 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
285 qam->SetTasks(fQATasks) ;
286 if (fWriteQAExpertData)
287 qam->SetWriteExpert() ;
289 if (qam->IsDefaultStorageSet()) {
290 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
291 AliWarning("Default QA reference storage has been already set !");
292 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fQARefUri.Data()));
293 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
294 fQARefUri = qam->GetDefaultStorage()->GetURI();
296 if (fQARefUri.Length() > 0) {
297 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
298 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
299 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
301 fQARefUri="local://$ALICE_ROOT/QARef";
302 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
303 AliWarning("Default QA reference storage not yet set !!!!");
304 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
305 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
307 qam->SetDefaultStorage(fQARefUri);
311 //_____________________________________________________________________________
312 void AliSimulation::InitCDB()
314 // activate a default CDB storage
315 // First check if we have any CDB storage set, because it is used
316 // to retrieve the calibration and alignment constants
318 if (fInitCDBCalled) return;
319 fInitCDBCalled = kTRUE;
321 AliCDBManager* man = AliCDBManager::Instance();
322 if (man->IsDefaultStorageSet())
324 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
325 AliWarning("Default CDB storage has been already set !");
326 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
327 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
328 fCDBUri = man->GetDefaultStorage()->GetURI();
331 if (fCDBUri.Length() > 0)
333 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
334 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
335 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
337 fCDBUri="local://$ALICE_ROOT/OCDB";
338 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
339 AliWarning("Default CDB storage not yet set !!!!");
340 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
341 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
344 man->SetDefaultStorage(fCDBUri);
347 // Now activate the detector specific CDB storage locations
348 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
349 TObject* obj = fSpecCDBUri[i];
351 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
352 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
353 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
354 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
359 //_____________________________________________________________________________
360 void AliSimulation::InitRunNumber(){
361 // check run number. If not set, set it to 0 !!!!
363 if (fInitRunNumberCalled) return;
364 fInitRunNumberCalled = kTRUE;
366 AliCDBManager* man = AliCDBManager::Instance();
367 if (man->GetRun() >= 0)
369 AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
370 "Use external variable DC_RUN or AliSimulation::SetRun()!"));
374 AliDebug(2, Form("Setting CDB run number to: %d",fRun));
377 AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
386 //_____________________________________________________________________________
387 void AliSimulation::SetCDBLock() {
388 // Set CDB lock: from now on it is forbidden to reset the run number
389 // or the default storage or to activate any further storage!
391 ULong64_t key = AliCDBManager::Instance()->SetLock(1);
395 //_____________________________________________________________________________
396 void AliSimulation::SetDefaultStorage(const char* uri) {
397 // Store the desired default CDB storage location
398 // Activate it later within the Run() method
404 //_____________________________________________________________________________
405 void AliSimulation::SetQARefDefaultStorage(const char* uri) {
406 // Store the desired default CDB storage location
407 // Activate it later within the Run() method
410 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
413 //_____________________________________________________________________________
414 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
415 // Store a detector-specific CDB storage location
416 // Activate it later within the Run() method
418 AliCDBPath aPath(calibType);
419 if(!aPath.IsValid()){
420 AliError(Form("Not a valid path: %s", calibType));
424 TObject* obj = fSpecCDBUri.FindObject(calibType);
425 if (obj) fSpecCDBUri.Remove(obj);
426 fSpecCDBUri.Add(new TNamed(calibType, uri));
430 //_____________________________________________________________________________
431 void AliSimulation::SetRunNumber(Int_t run)
434 // Activate it later within the Run() method
439 //_____________________________________________________________________________
440 void AliSimulation::SetSeed(Int_t seed)
443 // Activate it later within the Run() method
448 //_____________________________________________________________________________
449 Bool_t AliSimulation::SetRunNumberFromData()
451 // Set the CDB manager run number
452 // The run number is retrieved from gAlice
454 if (fSetRunNumberFromDataCalled) return kTRUE;
455 fSetRunNumberFromDataCalled = kTRUE;
457 AliCDBManager* man = AliCDBManager::Instance();
458 Int_t runData = -1, runCDB = -1;
460 AliRunLoader* runLoader = LoadRun("READ");
461 if (!runLoader) return kFALSE;
463 runData = runLoader->GetHeader()->GetRun();
467 runCDB = man->GetRun();
469 if (runCDB != runData) {
470 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
471 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
472 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
473 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
478 man->SetRun(runData);
481 if(man->GetRun() < 0) {
482 AliError("Run number not properly initalized!");
491 //_____________________________________________________________________________
492 void AliSimulation::SetConfigFile(const char* fileName)
494 // set the name of the config file
496 fConfigFileName = fileName;
499 //_____________________________________________________________________________
500 void AliSimulation::SetGAliceFile(const char* fileName)
502 // set the name of the galice file
503 // the path is converted to an absolute one if it is relative
505 fGAliceFileName = fileName;
506 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
507 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
509 fGAliceFileName = absFileName;
510 delete[] absFileName;
513 AliDebug(2, Form("galice file name set to %s", fileName));
516 //_____________________________________________________________________________
517 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
520 // set the number of events per file for the given detector and data type
521 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
523 TNamed* obj = new TNamed(detector, type);
524 obj->SetUniqueID(nEvents);
525 fEventsPerFile.Add(obj);
528 //_____________________________________________________________________________
529 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
531 // Read the alignment objects from CDB.
532 // Each detector is supposed to have the
533 // alignment objects in DET/Align/Data CDB path.
534 // All the detector objects are then collected,
535 // sorted by geometry level (starting from ALIC) and
536 // then applied to the TGeo geometry.
537 // Finally an overlaps check is performed.
539 if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
540 AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
544 // initialize CDB storage, run number, set CDB lock
546 // if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
549 Bool_t delRunLoader = kFALSE;
551 runLoader = LoadRun("READ");
552 if (!runLoader) return kFALSE;
553 delRunLoader = kTRUE;
556 // Export ideal geometry
557 if(!IsGeometryFromFile()) AliGeomManager::GetGeometry()->Export("geometry.root");
559 // Load alignment data from CDB and apply to geometry through AliGeomManager
560 if(fLoadAlignFromCDB){
562 TString detStr = fLoadAlObjsListOfDets;
563 TString loadAlObjsListOfDets = "";
565 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
566 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
567 AliModule* det = (AliModule*) detArray->At(iDet);
568 if (!det || !det->IsActive()) continue;
569 if (IsSelected(det->GetName(), detStr)) {
570 //add det to list of dets to be aligned from CDB
571 loadAlObjsListOfDets += det->GetName();
572 loadAlObjsListOfDets += " ";
574 } // end loop over detectors
575 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
576 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
578 // Check if the array with alignment objects was
579 // provided by the user. If yes, apply the objects
580 // to the present TGeo geometry
581 if (fAlignObjArray) {
582 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
583 AliError("The misalignment of one or more volumes failed!"
584 "Compare the list of simulated detectors and the list of detector alignment data!");
585 if (delRunLoader) delete runLoader;
591 // Update the internal geometry of modules (ITS needs it)
592 TString detStr = fLoadAlObjsListOfDets;
593 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
594 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
596 AliModule* det = (AliModule*) detArray->At(iDet);
597 if (!det || !det->IsActive()) continue;
598 if (IsSelected(det->GetName(), detStr)) {
599 det->UpdateInternalGeometry();
601 } // end loop over detectors
604 if (delRunLoader) delete runLoader;
609 //_____________________________________________________________________________
610 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
612 // add a file with background events for merging
614 TObjString* fileNameStr = new TObjString(fileName);
615 fileNameStr->SetUniqueID(nSignalPerBkgrd);
616 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
617 fBkgrdFileNames->Add(fileNameStr);
620 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
622 // add a file with background events for embeddin
623 MergeWith(fileName, nSignalPerBkgrd);
624 fEmbeddingFlag = kTRUE;
627 //_____________________________________________________________________________
628 Bool_t AliSimulation::Run(Int_t nEvents)
630 // run the generation, simulation and digitization
633 AliCodeTimerAuto("",0)
634 AliSysInfo::AddStamp("Start_Run");
636 // Load run number and seed from environmental vars
637 ProcessEnvironmentVars();
638 AliSysInfo::AddStamp("ProcessEnvironmentVars");
640 gRandom->SetSeed(fSeed);
642 if (nEvents > 0) fNEvents = nEvents;
644 // Run generator-only code on demand
645 if (fRunGeneratorOnly)
647 if(!RunGeneratorOnly())
649 if (fStopOnError) return kFALSE;
655 // create and setup the HLT instance
656 if (!fRunHLT.IsNull() && !CreateHLT()) {
657 if (fStopOnError) return kFALSE;
662 // generation and simulation -> hits
663 if (fRunGeneration) {
664 if (!RunSimulation()) if (fStopOnError) return kFALSE;
666 AliSysInfo::AddStamp("RunSimulation");
668 // initialize CDB storage from external environment
669 // (either CDB manager or AliSimulation setters),
670 // if not already done in RunSimulation()
672 AliSysInfo::AddStamp("InitCDB");
674 // Set run number in CDBManager from data
675 // From this point on the run number must be always loaded from data!
676 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
678 // Set CDB lock: from now on it is forbidden to reset the run number
679 // or the default storage or to activate any further storage!
682 // If RunSimulation was not called, load the geometry and misalign it
683 if (!AliGeomManager::GetGeometry()) {
684 // Initialize the geometry manager
685 AliGeomManager::LoadGeometry("geometry.root");
686 AliSysInfo::AddStamp("GetGeometry");
687 // // Check that the consistency of symbolic names for the activated subdetectors
688 // // in the geometry loaded by AliGeomManager
689 // AliRunLoader* runLoader = LoadRun("READ");
690 // if (!runLoader) return kFALSE;
692 // TString detsToBeChecked = "";
693 // TObjArray* detArray = runLoader->GetAliRun()->Detectors();
694 // for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
695 // AliModule* det = (AliModule*) detArray->At(iDet);
696 // if (!det || !det->IsActive()) continue;
697 // detsToBeChecked += det->GetName();
698 // detsToBeChecked += " ";
699 // } // end loop over detectors
700 // if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
701 if(!AliGeomManager::CheckSymNamesLUT("ALL"))
702 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
704 if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
706 if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
708 AliSysInfo::AddStamp("MissalignGeometry");
711 // hits -> summable digits
712 AliSysInfo::AddStamp("Start_sdigitization");
713 if (!fMakeSDigits.IsNull()) {
714 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
717 AliSysInfo::AddStamp("Stop_sdigitization");
719 AliSysInfo::AddStamp("Start_digitization");
720 // summable digits -> digits
721 if (!fMakeDigits.IsNull()) {
722 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
723 if (fStopOnError) return kFALSE;
726 AliSysInfo::AddStamp("Stop_digitization");
731 if (!fMakeDigitsFromHits.IsNull()) {
732 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
733 AliWarning(Form("Merging and direct creation of digits from hits "
734 "was selected for some detectors. "
735 "No merging will be done for the following detectors: %s",
736 fMakeDigitsFromHits.Data()));
738 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
739 if (fStopOnError) return kFALSE;
743 AliSysInfo::AddStamp("Hits2Digits");
746 // digits -> trigger. Set trigger configuration to "none" - any
747 // case - to not generate the trigger information. Set the trigger
748 // configuration to some string X to read from file at
749 // $ALICE_ROOT/GRP/CTP/X. Set the trigger configuration to the
750 // empty string or "ocdb" - any case - to read from OCDB.
751 if (!fTriggerConfig.EqualTo("none",TString::kIgnoreCase) &&
752 !RunTrigger(fTriggerConfig,fMakeDigits)) {
753 if (fStopOnError) return kFALSE;
756 AliSysInfo::AddStamp("RunTrigger");
759 // digits -> raw data
760 if (!fWriteRawData.IsNull()) {
761 if (!WriteRawData(fWriteRawData, fRawDataFileName,
762 fDeleteIntermediateFiles,fWriteSelRawData)) {
763 if (fStopOnError) return kFALSE;
767 AliSysInfo::AddStamp("WriteRaw");
769 // run HLT simulation on simulated digit data if raw data is not
770 // simulated, otherwise its called as part of WriteRawData
771 if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
773 if (fStopOnError) return kFALSE;
777 AliSysInfo::AddStamp("RunHLT");
781 Bool_t rv = RunQA() ;
787 AliSysInfo::AddStamp("RunQA");
791 TString snapshotFileOut("");
792 if(TString(gSystem->Getenv("OCDB_SNAPSHOT_CREATE")) == TString("kTRUE")){
793 AliInfo(" ******** Creating the snapshot! *********");
794 TString snapshotFile(gSystem->Getenv("OCDB_SNAPSHOT_FILENAME"));
795 if(!(snapshotFile.IsNull() || snapshotFile.IsWhitespace()))
796 snapshotFileOut = snapshotFile;
798 snapshotFileOut="OCDB.root";
799 AliCDBManager::Instance()->DumpToSnapshotFile(snapshotFileOut.Data(),kFALSE);
802 // Cleanup of CDB manager: cache and active storages!
803 AliCDBManager::Instance()->ClearCache();
808 //_______________________________________________________________________
809 Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
810 Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
811 Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
814 // Generates lego plots of:
815 // - radiation length map phi vs theta
816 // - radiation length map phi vs eta
817 // - interaction length map
818 // - g/cm2 length map
820 // ntheta bins in theta, eta
821 // themin minimum angle in theta (degrees)
822 // themax maximum angle in theta (degrees)
824 // phimin minimum angle in phi (degrees)
825 // phimax maximum angle in phi (degrees)
826 // rmin minimum radius
827 // rmax maximum radius
830 // The number of events generated = ntheta*nphi
831 // run input parameters in macro setup (default="Config.C")
833 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
836 <img src="picts/AliRunLego1.gif">
841 <img src="picts/AliRunLego2.gif">
846 <img src="picts/AliRunLego3.gif">
851 // run the generation and simulation
853 AliCodeTimerAuto("",0)
855 // initialize CDB storage and run number from external environment
856 // (either CDB manager or AliSimulation setters)
862 AliError("no gAlice object. Restart aliroot and try again.");
865 if (gAlice->Modules()->GetEntries() > 0) {
866 AliError("gAlice was already run. Restart aliroot and try again.");
870 AliInfo(Form("initializing gAlice with config file %s",
871 fConfigFileName.Data()));
874 if (nev == -1) nev = nc1 * nc2;
876 // check if initialisation has been done
877 // If runloader has been initialized, set the number of events per file to nc1 * nc2
880 if (!gener) gener = new AliLegoGenerator();
882 // Configure Generator
884 gener->SetRadiusRange(rmin, rmax);
885 gener->SetZMax(zmax);
886 gener->SetCoor1Range(nc1, c1min, c1max);
887 gener->SetCoor2Range(nc2, c2min, c2max);
891 fLego = new AliLego("lego",gener);
893 //__________________________________________________________________________
897 gROOT->LoadMacro(setup);
898 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
900 if(AliCDBManager::Instance()->GetRun() >= 0) {
901 SetRunNumber(AliCDBManager::Instance()->GetRun());
903 AliWarning("Run number not initialized!!");
906 AliRunLoader::Instance()->CdGAFile();
908 AliPDG::AddParticlesToPdgDataBase();
910 TVirtualMC::GetMC()->SetMagField(TGeoGlobalMagField::Instance()->GetField());
912 gAlice->GetMCApp()->Init();
915 //Must be here because some MCs (G4) adds detectors here and not in Config.C
916 gAlice->InitLoaders();
917 AliRunLoader::Instance()->MakeTree("E");
920 // Save stuff at the beginning of the file to avoid file corruption
921 AliRunLoader::Instance()->CdGAFile();
924 //Save current generator
925 AliGenerator *gen=gAlice->GetMCApp()->Generator();
926 gAlice->GetMCApp()->ResetGenerator(gener);
927 //Prepare MC for Lego Run
928 TVirtualMC::GetMC()->InitLego();
933 AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
934 TVirtualMC::GetMC()->ProcessRun(nev);
936 // End of this run, close files
938 // Restore current generator
939 gAlice->GetMCApp()->ResetGenerator(gen);
940 // Delete Lego Object
946 //_____________________________________________________________________________
947 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
951 AliCodeTimerAuto("",0)
953 // initialize CDB storage from external environment
954 // (either CDB manager or AliSimulation setters),
955 // if not already done in RunSimulation()
958 // Set run number in CDBManager from data
959 // From this point on the run number must be always loaded from data!
960 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
962 // Set CDB lock: from now on it is forbidden to reset the run number
963 // or the default storage or to activate any further storage!
966 AliRunLoader* runLoader = LoadRun("READ");
967 if (!runLoader) return kFALSE;
968 TString trconfiguration = config;
970 if (trconfiguration.IsNull()) {
971 if(!fTriggerConfig.IsNull()) {
972 trconfiguration = fTriggerConfig;
975 AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
978 runLoader->MakeTree( "GG" );
979 AliCentralTrigger* aCTP = runLoader->GetTrigger();
980 // Load Configuration
981 if (!aCTP->LoadConfiguration( trconfiguration ))
985 if( !aCTP->RunTrigger( runLoader , detectors ) ) {
997 //_____________________________________________________________________________
998 Bool_t AliSimulation::WriteTriggerRawData()
1000 // Writes the CTP (trigger) DDL raw data
1001 // Details of the format are given in the
1002 // trigger TDR - pages 134 and 135.
1003 AliCTPRawData writer;
1005 writer.RawDataRun2();
1010 //_____________________________________________________________________________
1011 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
1013 // run the generation and simulation
1015 AliCodeTimerAuto("",0)
1017 // initialize CDB storage and run number from external environment
1018 // (either CDB manager or AliSimulation setters)
1019 AliSysInfo::AddStamp("RunSimulation_Begin");
1021 AliSysInfo::AddStamp("RunSimulation_InitCDB");
1025 AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
1028 AliError("no gAlice object. Restart aliroot and try again.");
1031 if (gAlice->Modules()->GetEntries() > 0) {
1032 AliError("gAlice was already run. Restart aliroot and try again.");
1036 // Setup monitoring if requested
1037 gAlice->GetMCApp()->SetUseMonitoring(fUseMonitoring);
1039 AliInfo(Form("initializing gAlice with config file %s",
1040 fConfigFileName.Data()));
1043 // Initialize ALICE Simulation run
1048 // If requested set the mag. field from the GRP entry.
1049 // After this the field is loccked and cannot be changed by Config.C
1050 if (fUseMagFieldFromGRP) {
1052 grpM.ReadGRPEntry();
1054 AliInfo("Field is locked now. It cannot be changed in Config.C");
1058 TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
1059 gROOT->LoadMacro(fConfigFileName.Data());
1060 Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
1061 if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
1062 AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
1064 AliSysInfo::AddStamp("RunSimulation_Config");
1067 // If requested obtain the vertex position and vertex sigma_z from the CDB
1068 // This overwrites the settings from the Config.C
1069 if (fUseVertexFromCDB) {
1070 Double_t vtxPos[3] = {0., 0., 0.};
1071 Double_t vtxSig[3] = {0., 0., 0.};
1072 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1074 AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1076 if(vertex->GetXRes()>2.8) { // > pipe radius --> it's a dummy object, don't use it
1077 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1078 if (entry) vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
1082 vertex->GetXYZ(vtxPos);
1083 vertex->GetSigmaXYZ(vtxSig);
1084 AliInfo("Overwriting Config.C vertex settings !");
1085 AliInfo(Form("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
1086 vtxPos[0], vtxPos[1], vtxPos[2], vtxSig[2]));
1088 AliGenerator *gen = gAlice->GetMCApp()->Generator();
1089 gen->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]); // vertex position
1090 gen->SetSigmaZ(vtxSig[2]);
1095 // If requested we take the SOR and EOR time-stamps from the GRP and use them
1096 // in order to generate the event time-stamps
1097 if (fUseTimeStampFromCDB) {
1099 grpM.ReadGRPEntry();
1100 const AliGRPObject *grpObj = grpM.GetGRPData();
1101 if (!grpObj || (grpObj->GetTimeEnd() <= grpObj->GetTimeStart())) {
1102 AliError("Missing GRP or bad SOR/EOR time-stamps! Switching off the time-stamp generation from GRP!");
1103 fTimeStart = fTimeEnd = 0;
1104 fUseTimeStampFromCDB = kFALSE;
1107 fTimeStart = grpObj->GetTimeStart();
1108 fTimeEnd = grpObj->GetTimeEnd();
1112 if(AliCDBManager::Instance()->GetRun() >= 0) {
1113 AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
1114 AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
1116 AliWarning("Run number not initialized!!");
1119 AliRunLoader::Instance()->CdGAFile();
1122 AliPDG::AddParticlesToPdgDataBase();
1124 TVirtualMC::GetMC()->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1125 AliSysInfo::AddStamp("RunSimulation_GetField");
1126 gAlice->GetMCApp()->Init();
1127 AliSysInfo::AddStamp("RunSimulation_InitMCApp");
1129 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1130 gAlice->InitLoaders();
1131 AliRunLoader::Instance()->MakeTree("E");
1132 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1133 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1134 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1136 // Save stuff at the beginning of the file to avoid file corruption
1137 AliRunLoader::Instance()->CdGAFile();
1139 gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1140 AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1141 //___________________________________________________________________________________________
1143 AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
1145 // Set run number in CDBManager
1146 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
1148 AliRunLoader* runLoader = AliRunLoader::Instance();
1150 AliError(Form("gAlice has no run loader object. "
1151 "Check your config file: %s", fConfigFileName.Data()));
1154 SetGAliceFile(runLoader->GetFileName());
1156 // Misalign geometry
1157 #if ROOT_VERSION_CODE < 331527
1158 AliGeomManager::SetGeometry(gGeoManager);
1160 // Check that the consistency of symbolic names for the activated subdetectors
1161 // in the geometry loaded by AliGeomManager
1162 TString detsToBeChecked = "";
1163 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1164 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1165 AliModule* det = (AliModule*) detArray->At(iDet);
1166 if (!det || !det->IsActive()) continue;
1167 detsToBeChecked += det->GetName();
1168 detsToBeChecked += " ";
1169 } // end loop over detectors
1170 if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
1171 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
1172 MisalignGeometry(runLoader);
1173 AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
1176 // AliRunLoader* runLoader = AliRunLoader::Instance();
1177 // if (!runLoader) {
1178 // AliError(Form("gAlice has no run loader object. "
1179 // "Check your config file: %s", fConfigFileName.Data()));
1182 // SetGAliceFile(runLoader->GetFileName());
1184 if (!gAlice->GetMCApp()->Generator()) {
1185 AliError(Form("gAlice has no generator object. "
1186 "Check your config file: %s", fConfigFileName.Data()));
1190 // Write GRP entry corresponding to the setting found in Cofig.C
1193 AliSysInfo::AddStamp("RunSimulation_WriteGRP");
1195 if (nEvents <= 0) nEvents = fNEvents;
1197 // get vertex from background file in case of merging
1198 if (fUseBkgrdVertex &&
1199 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
1200 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
1201 const char* fileName = ((TObjString*)
1202 (fBkgrdFileNames->At(0)))->GetName();
1203 AliInfo(Form("The vertex will be taken from the background "
1204 "file %s with nSignalPerBackground = %d",
1205 fileName, signalPerBkgrd));
1206 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
1207 gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
1210 if (!fRunSimulation) {
1211 gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
1214 // set the number of events per file for given detectors and data types
1215 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
1216 if (!fEventsPerFile[i]) continue;
1217 const char* detName = fEventsPerFile[i]->GetName();
1218 const char* typeName = fEventsPerFile[i]->GetTitle();
1219 TString loaderName(detName);
1220 loaderName += "Loader";
1221 AliLoader* loader = runLoader->GetLoader(loaderName);
1223 AliError(Form("RunSimulation no loader for %s found\n Number of events per file not set for %s %s",
1224 detName, typeName, detName));
1227 AliDataLoader* dataLoader =
1228 loader->GetDataLoader(typeName);
1230 AliError(Form("no data loader for %s found\n"
1231 "Number of events per file not set for %s %s",
1232 typeName, detName, typeName));
1235 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
1236 AliDebug(1, Form("number of events per file set to %d for %s %s",
1237 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
1240 AliInfo("running gAlice");
1241 AliSysInfo::AddStamp("Start_ProcessRun");
1243 // Create the Root Tree with one branch per detector
1244 //Hits moved to begin event -> now we are crating separate tree for each event
1245 TVirtualMC::GetMC()->ProcessRun(nEvents);
1247 // End of this run, close files
1248 if(nEvents>0) FinishRun();
1250 AliSysInfo::AddStamp("Stop_ProcessRun");
1256 //_____________________________________________________________________________
1257 Bool_t AliSimulation::RunGeneratorOnly()
1260 TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
1261 gROOT->LoadMacro(fConfigFileName.Data());
1262 Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
1263 if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
1264 AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
1267 // Setup the runloader and generator, check if everything is OK
1268 AliRunLoader* runLoader = AliRunLoader::Instance();
1269 AliGenerator* generator = gAlice->GetMCApp()->Generator();
1271 AliError(Form("gAlice has no run loader object. "
1272 "Check your config file: %s", fConfigFileName.Data()));
1276 AliError(Form("gAlice has no generator object. "
1277 "Check your config file: %s", fConfigFileName.Data()));
1281 runLoader->LoadKinematics("RECREATE");
1282 runLoader->MakeTree("E");
1284 // Create stack and header
1285 runLoader->MakeStack();
1286 AliStack* stack = runLoader->Stack();
1287 AliHeader* header = runLoader->GetHeader();
1289 // Intialize generator
1291 generator->SetStack(stack);
1293 // Run main generator loop
1295 for (Int_t iev=0; iev<fNEvents; iev++)
1298 header->Reset(0,iev);
1299 runLoader->SetEventNumber(iev);
1301 runLoader->MakeTree("K");
1304 generator->Generate();
1307 header->SetNprimary(stack->GetNprimary());
1308 header->SetNtrack(stack->GetNtrack());
1309 stack->FinishEvent();
1310 header->SetStack(stack);
1311 runLoader->TreeE()->Fill();
1312 runLoader->WriteKinematics("OVERWRITE");
1316 generator->FinishRun();
1318 runLoader->WriteHeader("OVERWRITE");
1325 //_____________________________________________________________________________
1326 Bool_t AliSimulation::RunSDigitization(const char* detectors)
1328 // run the digitization and produce summable digits
1329 static Int_t eventNr=0;
1330 AliCodeTimerAuto("",0) ;
1332 // initialize CDB storage, run number, set CDB lock
1334 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1337 AliRunLoader* runLoader = LoadRun();
1338 if (!runLoader) return kFALSE;
1340 TString detStr = detectors;
1341 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1342 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1343 AliModule* det = (AliModule*) detArray->At(iDet);
1344 if (!det || !det->IsActive()) continue;
1345 if (IsSelected(det->GetName(), detStr)) {
1346 AliInfo(Form("creating summable digits for %s", det->GetName()));
1347 AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
1348 det->Hits2SDigits();
1349 AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
1350 AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
1354 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1355 AliError(Form("the following detectors were not found: %s",
1357 if (fStopOnError) return kFALSE;
1366 //_____________________________________________________________________________
1367 Bool_t AliSimulation::RunDigitization(const char* detectors,
1368 const char* excludeDetectors)
1370 // run the digitization and produce digits from sdigits
1372 AliCodeTimerAuto("",0)
1374 // initialize CDB storage, run number, set CDB lock
1376 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1379 delete AliRunLoader::Instance();
1384 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1385 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1386 AliDigitizationInput digInp(nStreams, signalPerBkgrd);
1387 // digInp.SetEmbeddingFlag(fEmbeddingFlag);
1388 digInp.SetRegionOfInterest(fRegionOfInterest);
1389 digInp.SetInputStream(0, fGAliceFileName.Data());
1390 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1391 const char* fileName = ((TObjString*)(fBkgrdFileNames->At(iStream-1)))->GetName();
1392 digInp.SetInputStream(iStream, fileName);
1395 detArr.SetOwner(kTRUE);
1396 TString detStr = detectors;
1397 TString detExcl = excludeDetectors;
1398 if (!static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice()) {
1399 AliError("Error occured while getting gAlice from Input 0");
1402 AliRunLoader* runLoader = AliRunLoader::GetRunLoader(digInp.GetInputStream(0)->GetFolderName());
1403 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1404 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1405 AliModule* det = (AliModule*) detArray->At(iDet);
1406 if (!det || !det->IsActive()) continue;
1407 if (!IsSelected(det->GetName(), detStr) || IsSelected(det->GetName(), detExcl)) continue;
1408 AliDigitizer* digitizer = det->CreateDigitizer(&digInp);
1409 if (!digitizer || !digitizer->Init()) {
1410 AliError(Form("no digitizer for %s", det->GetName()));
1411 if (fStopOnError) return kFALSE;
1414 detArr.AddLast(digitizer);
1415 AliInfo(Form("Created digitizer from SDigits -> Digits for %s", det->GetName()));
1418 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1419 AliError(Form("the following detectors were not found: %s", detStr.Data()));
1420 if (fStopOnError) return kFALSE;
1423 Int_t ndigs = detArr.GetEntriesFast();
1424 Int_t eventsCreated = 0;
1425 AliRunLoader* outRl = digInp.GetOutRunLoader();
1426 while ((eventsCreated++ < fNEvents) || (fNEvents < 0)) {
1427 if (!digInp.ConnectInputTrees()) break;
1428 digInp.InitEvent(); //this must be after call of Connect Input tress.
1429 if (outRl) outRl->SetEventNumber(eventsCreated-1);
1430 static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice(); // use gAlice of the first input stream
1431 for (int id=0;id<ndigs;id++) ((AliDigitizer*)detArr[id])->Digitize("");
1432 digInp.FinishEvent();
1434 digInp.FinishGlobal();
1439 //_____________________________________________________________________________
1440 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1442 // run the digitization and produce digits from hits
1444 AliCodeTimerAuto("",0)
1446 // initialize CDB storage, run number, set CDB lock
1448 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1451 AliRunLoader* runLoader = LoadRun("READ");
1452 if (!runLoader) return kFALSE;
1454 TString detStr = detectors;
1455 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1456 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1457 AliModule* det = (AliModule*) detArray->At(iDet);
1458 if (!det || !det->IsActive()) continue;
1459 if (IsSelected(det->GetName(), detStr)) {
1460 AliInfo(Form("creating digits from hits for %s", det->GetName()));
1465 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1466 AliError(Form("the following detectors were not found: %s",
1468 if (fStopOnError) return kFALSE;
1474 //_____________________________________________________________________________
1475 Bool_t AliSimulation::WriteRawData(const char* detectors,
1476 const char* fileName,
1477 Bool_t deleteIntermediateFiles,
1480 // convert the digits to raw data
1481 // First DDL raw data files for the given detectors are created.
1482 // If a file name is given, the DDL files are then converted to a DATE file.
1483 // If deleteIntermediateFiles is true, the DDL raw files are deleted
1485 // If the file name has the extension ".root", the DATE file is converted
1487 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1488 // 'selrawdata' flag can be used to enable writing of detectors raw data
1489 // accoring to the trigger cluster.
1491 AliCodeTimerAuto("",0)
1492 AliSysInfo::AddStamp("WriteRawData_Start");
1494 TString detStr = detectors;
1495 if (!WriteRawFiles(detStr.Data())) {
1496 if (fStopOnError) return kFALSE;
1498 AliSysInfo::AddStamp("WriteRawFiles");
1500 // run HLT simulation on simulated DDL raw files
1501 // and produce HLT ddl raw files to be included in date/root file
1502 // bugfix 2009-06-26: the decision whether to write HLT raw data
1503 // is taken in RunHLT. Here HLT always needs to be run in order to
1504 // create HLT digits, unless its switched off. This is due to the
1505 // special placement of the HLT between the generation of DDL files
1506 // and conversion to DATE/Root file.
1507 detStr.ReplaceAll("HLT", "");
1508 if (!fRunHLT.IsNull()) {
1510 if (fStopOnError) return kFALSE;
1513 AliSysInfo::AddStamp("WriteRawData_RunHLT");
1515 TString dateFileName(fileName);
1516 if (!dateFileName.IsNull()) {
1517 Bool_t rootOutput = dateFileName.EndsWith(".root");
1518 if (rootOutput) dateFileName += ".date";
1519 TString selDateFileName;
1521 selDateFileName = "selected.";
1522 selDateFileName+= dateFileName;
1524 if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1525 if (fStopOnError) return kFALSE;
1527 AliSysInfo::AddStamp("ConvertRawFilesToDate");
1528 if (deleteIntermediateFiles) {
1529 AliRunLoader* runLoader = LoadRun("READ");
1530 if (runLoader) for (Int_t iEvent = 0;
1531 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1533 snprintf(command, 256, "rm -r raw%d", iEvent);
1534 gSystem->Exec(command);
1540 if (!ConvertDateToRoot(dateFileName, fileName)) {
1541 if (fStopOnError) return kFALSE;
1543 AliSysInfo::AddStamp("ConvertDateToRoot");
1544 if (deleteIntermediateFiles) {
1545 gSystem->Unlink(dateFileName);
1548 TString selFileName = "selected.";
1549 selFileName += fileName;
1550 if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1551 if (fStopOnError) return kFALSE;
1553 if (deleteIntermediateFiles) {
1554 gSystem->Unlink(selDateFileName);
1563 //_____________________________________________________________________________
1564 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1566 // convert the digits to raw data DDL files
1568 AliCodeTimerAuto("",0)
1570 AliRunLoader* runLoader = LoadRun("READ");
1571 if (!runLoader) return kFALSE;
1573 // write raw data to DDL files
1574 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1575 AliInfo(Form("processing event %d", iEvent));
1576 runLoader->GetEvent(iEvent);
1577 TString baseDir = gSystem->WorkingDirectory();
1579 snprintf(dirName, 256, "raw%d", iEvent);
1580 gSystem->MakeDirectory(dirName);
1581 if (!gSystem->ChangeDirectory(dirName)) {
1582 AliError(Form("couldn't change to directory %s", dirName));
1583 if (fStopOnError) return kFALSE; else continue;
1586 ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1589 TString detStr = detectors;
1590 if (IsSelected("HLT", detStr)) {
1591 // Do nothing. "HLT" will be removed from detStr and HLT raw
1592 // data files are generated in RunHLT.
1595 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1596 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1597 AliModule* det = (AliModule*) detArray->At(iDet);
1598 if (!det || !det->IsActive()) continue;
1599 if (IsSelected(det->GetName(), detStr)) {
1600 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1605 if (!WriteTriggerRawData())
1606 if (fStopOnError) return kFALSE;
1608 gSystem->ChangeDirectory(baseDir);
1609 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1610 AliError(Form("the following detectors were not found: %s",
1612 if (fStopOnError) return kFALSE;
1621 //_____________________________________________________________________________
1622 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1623 const char* selDateFileName)
1625 // convert raw data DDL files to a DATE file with the program "dateStream"
1626 // The second argument is not empty when the user decides to write
1627 // the detectors raw data according to the trigger cluster.
1629 AliCodeTimerAuto("",0)
1631 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1633 AliError("the program dateStream was not found");
1634 if (fStopOnError) return kFALSE;
1639 AliRunLoader* runLoader = LoadRun("READ");
1640 if (!runLoader) return kFALSE;
1642 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1643 Bool_t selrawdata = kFALSE;
1644 if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1647 // Note the option -s. It is used in order to avoid
1648 // the generation of SOR/EOR events.
1649 snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d",
1650 dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1651 FILE* pipe = gSystem->OpenPipe(command, "w");
1654 AliError(Form("Cannot execute command: %s",command));
1658 Int_t selEvents = 0;
1659 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1661 UInt_t detectorPattern = 0;
1662 runLoader->GetEvent(iEvent);
1663 if (!runLoader->LoadTrigger()) {
1664 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1665 detectorPattern = aCTP->GetClusterMask();
1666 // Check if the event was triggered by CTP
1668 if (aCTP->GetClassMask()) selEvents++;
1672 AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1674 AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1675 selrawdata = kFALSE;
1679 fprintf(pipe, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
1683 // loop over detectors and DDLs
1684 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1685 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1687 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1688 Int_t ldcID = Int_t(ldc + 0.0001);
1689 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1691 char rawFileName[256];
1692 snprintf(rawFileName, 256, "raw%d/%s",
1693 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1695 // check existence and size of raw data file
1696 FILE* file = fopen(rawFileName, "rb");
1697 if (!file) continue;
1698 fseek(file, 0, SEEK_END);
1699 unsigned long size = ftell(file);
1701 if (!size) continue;
1703 if (ldcID != prevLDC) {
1704 fprintf(pipe, " LDC Id %d\n", ldcID);
1707 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1712 Int_t result = gSystem->ClosePipe(pipe);
1714 if (!(selrawdata && selEvents > 0)) {
1716 return (result == 0);
1719 AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1721 snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d",
1722 selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1723 FILE* pipe2 = gSystem->OpenPipe(command, "w");
1725 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1727 // Get the trigger decision and cluster
1728 UInt_t detectorPattern = 0;
1730 runLoader->GetEvent(iEvent);
1731 if (!runLoader->LoadTrigger()) {
1732 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1733 if (aCTP->GetClassMask() == 0) continue;
1734 detectorPattern = aCTP->GetClusterMask();
1735 detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1736 AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1739 fprintf(pipe2, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
1743 // loop over detectors and DDLs
1744 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1745 // Write only raw data from detectors that
1746 // are contained in the trigger cluster(s)
1747 if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1749 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1751 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1752 Int_t ldcID = Int_t(ldc + 0.0001);
1753 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1755 char rawFileName[256];
1756 snprintf(rawFileName, 256, "raw%d/%s",
1757 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1759 // check existence and size of raw data file
1760 FILE* file = fopen(rawFileName, "rb");
1761 if (!file) continue;
1762 fseek(file, 0, SEEK_END);
1763 unsigned long size = ftell(file);
1765 if (!size) continue;
1767 if (ldcID != prevLDC) {
1768 fprintf(pipe2, " LDC Id %d\n", ldcID);
1771 fprintf(pipe2, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1776 Int_t result2 = gSystem->ClosePipe(pipe2);
1779 return ((result == 0) && (result2 == 0));
1782 //_____________________________________________________________________________
1783 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1784 const char* rootFileName)
1786 // convert a DATE file to a root file with the program "alimdc"
1789 const Int_t kDBSize = 2000000000;
1790 const Int_t kTagDBSize = 1000000000;
1791 const Bool_t kFilter = kFALSE;
1792 const Int_t kCompression = 1;
1794 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1796 AliError("the program alimdc was not found");
1797 if (fStopOnError) return kFALSE;
1802 AliInfo(Form("converting DATE file %s to root file %s",
1803 dateFileName, rootFileName));
1805 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1806 const char* tagDBFS = "/tmp/mdc1/tags";
1808 // User defined file system locations
1809 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1810 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1811 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1812 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1813 if (gSystem->Getenv("ALIMDC_TAGDB"))
1814 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1816 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1817 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1818 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1820 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1821 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1822 gSystem->Exec(Form("mkdir %s",tagDBFS));
1824 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1825 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1826 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1828 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1829 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1830 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1832 return (result == 0);
1836 //_____________________________________________________________________________
1837 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1839 // delete existing run loaders, open a new one and load gAlice
1841 delete AliRunLoader::Instance();
1842 AliRunLoader* runLoader =
1843 AliRunLoader::Open(fGAliceFileName.Data(),
1844 AliConfig::GetDefaultEventFolderName(), mode);
1846 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1849 runLoader->LoadgAlice();
1850 runLoader->LoadHeader();
1851 gAlice = runLoader->GetAliRun();
1853 AliError(Form("no gAlice object found in file %s",
1854 fGAliceFileName.Data()));
1860 //_____________________________________________________________________________
1861 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1863 // get or calculate the number of signal events per background event
1865 if (!fBkgrdFileNames) return 1;
1866 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1867 if (nBkgrdFiles == 0) return 1;
1869 // get the number of signal events
1871 AliRunLoader* runLoader =
1872 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1873 if (!runLoader) return 1;
1875 nEvents = runLoader->GetNumberOfEvents();
1880 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1881 // get the number of background events
1882 const char* fileName = ((TObjString*)
1883 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1884 AliRunLoader* runLoader =
1885 AliRunLoader::Open(fileName, "BKGRD");
1886 if (!runLoader) continue;
1887 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1890 // get or calculate the number of signal per background events
1891 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1892 if (nSignalPerBkgrd <= 0) {
1893 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1894 } else if (result && (result != nSignalPerBkgrd)) {
1895 AliInfo(Form("the number of signal events per background event "
1896 "will be changed from %d to %d for stream %d",
1897 nSignalPerBkgrd, result, iBkgrdFile+1));
1898 nSignalPerBkgrd = result;
1901 if (!result) result = nSignalPerBkgrd;
1902 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1903 AliWarning(Form("not enough background events (%d) for %d signal events "
1904 "using %d signal per background events for stream %d",
1905 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1912 //_____________________________________________________________________________
1913 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1915 // check whether detName is contained in detectors
1916 // if yes, it is removed from detectors
1918 // check if all detectors are selected
1919 if ((detectors.CompareTo("ALL") == 0) ||
1920 detectors.BeginsWith("ALL ") ||
1921 detectors.EndsWith(" ALL") ||
1922 detectors.Contains(" ALL ")) {
1927 // search for the given detector
1928 Bool_t result = kFALSE;
1929 if ((detectors.CompareTo(detName) == 0) ||
1930 detectors.BeginsWith(detName+" ") ||
1931 detectors.EndsWith(" "+detName) ||
1932 detectors.Contains(" "+detName+" ")) {
1933 detectors.ReplaceAll(detName, "");
1937 // clean up the detectors string
1938 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1939 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1940 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1945 //_____________________________________________________________________________
1946 Int_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName, Int_t N, Int_t nSkip)
1949 // Steering routine to convert raw data in directory rawDirectory/ to fake SDigits.
1950 // These can be used for embedding of MC tracks into RAW data using the standard
1951 // merging procedure.
1953 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1956 AliError("no gAlice object. Restart aliroot and try again.");
1959 if (gAlice->Modules()->GetEntries() > 0) {
1960 AliError("gAlice was already run. Restart aliroot and try again.");
1964 AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1968 gROOT->LoadMacro(fConfigFileName.Data());
1969 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1971 if(AliCDBManager::Instance()->GetRun() >= 0) {
1972 SetRunNumber(AliCDBManager::Instance()->GetRun());
1974 AliWarning("Run number not initialized!!");
1977 AliRunLoader::Instance()->CdGAFile();
1979 AliPDG::AddParticlesToPdgDataBase();
1981 TVirtualMC::GetMC()->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1983 gAlice->GetMCApp()->Init();
1985 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1986 gAlice->InitLoaders();
1987 AliRunLoader::Instance()->MakeTree("E");
1988 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1989 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1990 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1993 // Save stuff at the beginning of the file to avoid file corruption
1994 AliRunLoader::Instance()->CdGAFile();
1999 //AliCDBManager* man = AliCDBManager::Instance();
2000 //man->SetRun(0); // Should this come from rawdata header ?
2004 // Get the runloader
2005 AliRunLoader* runLoader = AliRunLoader::Instance();
2007 // Open esd file if available
2010 AliESDEvent* esd = 0;
2011 if (esdFileName && (strlen(esdFileName)>0)) {
2012 esdFile = TFile::Open(esdFileName);
2014 esd = new AliESDEvent();
2015 esdFile->GetObject("esdTree", treeESD);
2017 esd->ReadFromTree(treeESD);
2019 AliInfo(Form("Asking to skip first %d ESDs events",nSkip));
2028 // Create the RawReader
2029 TString fileName(rawDirectory);
2030 AliRawReader* rawReader = AliRawReader::Create(fileName.Data());
2031 if (!rawReader) return (kFALSE);
2033 // if (!fEquipIdMap.IsNull() && fRawReader)
2034 // fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
2036 // Get list of detectors
2037 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
2040 AliHeader* header = runLoader->GetHeader();
2044 if (!(rawReader->NextEvent())) break;
2045 runLoader->SetEventNumber(nev);
2046 runLoader->GetHeader()->Reset(rawReader->GetRunNumber(),
2048 runLoader->GetEvent(nev);
2049 AliInfo(Form("We are at event %d",nev));
2052 TString detStr = fMakeSDigits;
2053 for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
2054 AliModule* det = (AliModule*) detArray->At(iDet);
2055 if (!det || !det->IsActive()) continue;
2056 if (IsSelected(det->GetName(), detStr)) {
2057 AliInfo(Form("Calling Raw2SDigits for %s", det->GetName()));
2058 det->Raw2SDigits(rawReader);
2065 // If ESD information available obtain reconstructed vertex and store in header.
2067 AliInfo(Form("Selected event %d correspond to event %d in raw and to %d in esd",nev,rawReader->GetEventIndex(),nSkip+rawReader->GetEventIndex()));
2068 treeESD->GetEvent(nSkip+rawReader->GetEventIndex());
2069 const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
2070 Double_t position[3];
2071 esdVertex->GetXYZ(position);
2072 AliGenEventHeader* mcHeader = new AliGenEventHeader("ESD");
2075 for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
2076 mcHeader->SetPrimaryVertex(mcV);
2077 header->Reset(0,nev);
2078 header->SetGenEventHeader(mcHeader);
2079 AliInfo(Form("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]));
2083 runLoader->TreeE()->Fill();
2084 AliInfo(Form("Finished event %d",nev));
2093 runLoader->CdGAFile();
2094 runLoader->WriteHeader("OVERWRITE");
2095 runLoader->WriteRunLoader();
2100 //_____________________________________________________________________________
2101 void AliSimulation::FinishRun()
2104 // Called at the end of the run.
2109 AliDebug(1, "Finish Lego");
2110 AliRunLoader::Instance()->CdGAFile();
2114 // Clean detector information
2115 TIter next(gAlice->Modules());
2116 AliModule *detector;
2117 while((detector = dynamic_cast<AliModule*>(next()))) {
2118 AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
2119 detector->FinishRun();
2122 AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
2123 AliRunLoader::Instance()->WriteHeader("OVERWRITE");
2125 // Write AliRun info and all detectors parameters
2126 AliRunLoader::Instance()->CdGAFile();
2127 gAlice->Write(0,TObject::kOverwrite);//write AliRun
2128 AliRunLoader::Instance()->Write(0,TObject::kOverwrite);//write RunLoader itself
2130 if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();
2131 AliRunLoader::Instance()->Synchronize();
2134 //_____________________________________________________________________________
2135 Int_t AliSimulation::GetDetIndex(const char* detector)
2137 // return the detector index corresponding to detector
2139 for (index = 0; index < fgkNDetectors ; index++) {
2140 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2146 //_____________________________________________________________________________
2147 Bool_t AliSimulation::CreateHLT()
2149 // Init the HLT simulation.
2150 // The function loads the library and creates the instance of AliHLTSimulation.
2151 // the main reason for the decoupled creation is to set the transient OCDB
2152 // objects before the OCDB is locked
2154 // load the library dynamically
2155 gSystem->Load(ALIHLTSIMULATION_LIBRARY);
2157 // check for the library version
2158 AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
2160 AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
2163 if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
2164 AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
2167 // print compile info
2168 typedef void (*CompileInfo)( const char*& date, const char*& time);
2169 CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
2171 const char* date="";
2172 const char* time="";
2173 (*fctInfo)(date, time);
2174 if (!date) date="unknown";
2175 if (!time) time="unknown";
2176 AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
2178 AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
2181 // create instance of the HLT simulation
2182 AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
2183 if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
2184 AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
2188 TString specObjects;
2189 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
2190 if (specObjects.Length()>0) specObjects+=" ";
2191 specObjects+=fSpecCDBUri[i]->GetName();
2194 AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
2195 if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
2196 AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
2202 //_____________________________________________________________________________
2203 Bool_t AliSimulation::RunHLT()
2205 // Run the HLT simulation
2206 // HLT simulation is implemented in HLT/sim/AliHLTSimulation
2207 // Disabled if fRunHLT is empty, default vaule is "default".
2208 // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
2209 // The default simulation depends on the HLT component libraries and their
2210 // corresponding agents which define components and chains to run. See
2211 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
2212 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
2214 // The libraries to be loaded can be specified as an option.
2216 // AliSimulation sim;
2217 // sim.SetRunHLT("libAliHLTSample.so");
2219 // will only load <tt>libAliHLTSample.so</tt>
2221 // Other available options:
2222 // \li loglevel=<i>level</i> <br>
2223 // logging level for this processing
2225 // disable redirection of log messages to AliLog class
2226 // \li config=<i>macro</i>
2227 // configuration macro
2228 // \li chains=<i>configuration</i>
2229 // comma separated list of configurations to be run during simulation
2230 // \li rawfile=<i>file</i>
2231 // source for the RawReader to be created, the default is <i>./</i> if
2232 // raw data is simulated
2236 if (!fpHLT && !CreateHLT()) {
2239 AliHLTSimulation* pHLT=fpHLT;
2241 AliRunLoader* pRunLoader = LoadRun("READ");
2242 if (!pRunLoader) return kFALSE;
2244 // initialize CDB storage, run number, set CDB lock
2245 // thats for the case of running HLT simulation without all the other steps
2246 // multiple calls are handled by the function, so we can just call
2248 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
2251 // init the HLT simulation
2253 if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2254 TString detStr = fWriteRawData;
2255 if (!IsSelected("HLT", detStr)) {
2256 options+=" writerawfiles=";
2258 options+=" writerawfiles=HLT";
2261 if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2262 // as a matter of fact, HLT will run reconstruction and needs the RawReader
2263 // in order to get detector data. By default, RawReaderFile is used to read
2264 // the already simulated ddl files. Date and Root files from the raw data
2265 // are generated after the HLT simulation.
2266 options+=" rawfile=./";
2269 AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2270 if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2271 AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2273 // run the HLT simulation
2274 AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2275 if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2276 AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2280 // delete the instance
2281 AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2282 if (fctDelete==NULL || fctDelete(pHLT)<0) {
2283 AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2287 return iResult>=0?kTRUE:kFALSE;
2290 //_____________________________________________________________________________
2291 Bool_t AliSimulation::RunQA()
2293 // run the QA on summable hits, digits or digits
2295 //if(!gAlice) return kFALSE;
2296 AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2298 TString detectorsw("") ;
2300 AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2301 detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ;
2302 if ( detectorsw.IsNull() )
2307 //_____________________________________________________________________________
2308 Bool_t AliSimulation::SetRunQA(TString detAndAction)
2310 // Allows to run QA for a selected set of detectors
2311 // and a selected set of tasks among HITS, SDIGITS and DIGITS
2312 // all selected detectors run the same selected tasks
2314 if (!detAndAction.Contains(":")) {
2315 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2319 Int_t colon = detAndAction.Index(":") ;
2320 fQADetectors = detAndAction(0, colon) ;
2321 if (fQADetectors.Contains("ALL") ){
2322 TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2323 Int_t minus = fQADetectors.Last('-') ;
2324 TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2325 TString toRemove("") ;
2326 while (minus >= 0) {
2327 toRemove = fQADetectors(minus+1, fQADetectors.Length()) ;
2328 toRemove = toRemove.Strip() ;
2329 toKeep.ReplaceAll(toRemove, "") ;
2330 fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ;
2331 minus = fQADetectors.Last('-') ;
2333 fQADetectors = toKeep ;
2335 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2336 if (fQATasks.Contains("ALL") ) {
2337 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
2339 fQATasks.ToUpper() ;
2341 if ( fQATasks.Contains("HIT") )
2342 tempo = Form("%d ", AliQAv1::kHITS) ;
2343 if ( fQATasks.Contains("SDIGIT") )
2344 tempo += Form("%d ", AliQAv1::kSDIGITS) ;
2345 if ( fQATasks.Contains("DIGIT") )
2346 tempo += Form("%d ", AliQAv1::kDIGITS) ;
2348 if (fQATasks.IsNull()) {
2349 AliInfo("No QA requested\n") ;
2354 TString tempo(fQATasks) ;
2355 tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS)) ;
2356 tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;
2357 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;
2358 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2360 AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ;
2361 AliQAManager::QAManager()->SetTasks(fQATasks) ;
2362 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++)
2363 AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2368 //_____________________________________________________________________________
2369 void AliSimulation::ProcessEnvironmentVars()
2371 // Extract run number and random generator seed from env variables
2373 AliInfo("Processing environment variables");
2375 // Random Number seed
2377 // first check that seed is not already set
2379 if (gSystem->Getenv("CONFIG_SEED")) {
2380 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2383 if (gSystem->Getenv("CONFIG_SEED")) {
2384 AliInfo(Form("Seed for random number generation already set (%d)"
2385 ": CONFIG_SEED variable ignored!", fSeed));
2389 AliInfo(Form("Seed for random number generation = %d ", fSeed));
2393 // first check that run number is not already set
2395 if (gSystem->Getenv("DC_RUN")) {
2396 fRun = atoi(gSystem->Getenv("DC_RUN"));
2399 if (gSystem->Getenv("DC_RUN")) {
2400 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2404 AliInfo(Form("Run number = %d", fRun));
2407 //---------------------------------------------------------------------
2408 void AliSimulation::WriteGRPEntry()
2410 // Get the necessary information from galice (generator, trigger etc) and
2411 // write a GRP entry corresponding to the settings in the Config.C used
2412 // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2415 AliInfo("Writing global run parameters entry into the OCDB");
2417 AliGRPObject* grpObj = new AliGRPObject();
2419 grpObj->SetRunType("PHYSICS");
2420 grpObj->SetTimeStart(fTimeStart);
2421 grpObj->SetTimeEnd(fTimeEnd);
2422 grpObj->SetBeamEnergyIsSqrtSHalfGeV(); // new format of GRP: store sqrt(s)/2 in GeV
2424 const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2430 gen->GetProjectile(projectile,a,z);
2432 gen->GetTarget(target,a,z);
2433 TString beamType = projectile + "-" + target;
2434 beamType.ReplaceAll(" ","");
2435 if (!beamType.CompareTo("-")) {
2436 grpObj->SetBeamType("UNKNOWN");
2437 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
2440 grpObj->SetBeamType(beamType);
2442 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 * a / z);
2444 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 );
2446 // Heavy ion run, the event specie is set to kHighMult
2447 fEventSpecie = AliRecoParam::kHighMult;
2448 if ((strcmp(beamType,"p-p") == 0) ||
2449 (strcmp(beamType,"p-") == 0) ||
2450 (strcmp(beamType,"-p") == 0) ||
2451 (strcmp(beamType,"P-P") == 0) ||
2452 (strcmp(beamType,"P-") == 0) ||
2453 (strcmp(beamType,"-P") == 0)) {
2454 // Proton run, the event specie is set to kLowMult
2455 fEventSpecie = AliRecoParam::kLowMult;
2459 AliWarning("Unknown beam type and energy! Setting energy to 0");
2460 grpObj->SetBeamEnergy(0);
2461 grpObj->SetBeamType("UNKNOWN");
2464 UInt_t detectorPattern = 0;
2466 TObjArray *detArray = gAlice->Detectors();
2467 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2468 if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2469 AliDebug(1, Form("Detector #%d found: %s", iDet, AliDAQ::OfflineModuleName(iDet)));
2470 detectorPattern |= (1 << iDet);
2475 if (!fTriggerConfig.IsNull())
2476 detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2479 if (!fRunHLT.IsNull())
2480 detectorPattern |= (1 << AliDAQ::kHLTId);
2482 grpObj->SetNumberOfDetectors((Char_t)nDets);
2483 grpObj->SetDetectorMask((Int_t)detectorPattern);
2484 grpObj->SetLHCPeriod("LHC08c");
2485 grpObj->SetLHCState("STABLE_BEAMS");
2487 AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2488 Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2490 Float_t factorSol = field ? field->GetFactorSol() : 0;
2491 Float_t currentSol = TMath::Abs(factorSol)>1E-6 ?
2492 TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
2494 Float_t factorDip = field ? field->GetFactorDip() : 0;
2495 Float_t currentDip = 6000.*TMath::Abs(factorDip);
2497 grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
2498 grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);
2499 grpObj->SetL3Polarity(factorSol>0 ? 0:1);
2500 grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
2501 if (field) grpObj->SetUniformBMap(field->IsUniform()); // for special MC with k5kGUniform map
2502 grpObj->SetPolarityConventionLHC(); // LHC convention +/+ current -> -/- field main components
2504 grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2506 //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2508 // Now store the entry in OCDB
2509 AliCDBManager* man = AliCDBManager::Instance();
2511 man->SetLock(0, fKey);
2513 AliCDBStorage* sto = man->GetStorage(fGRPWriteLocation.Data());
2516 AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2517 AliCDBMetaData *metadata= new AliCDBMetaData();
2519 metadata->SetResponsible("alice-off@cern.ch");
2520 metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2522 sto->Put(grpObj,id,metadata);
2523 man->SetLock(1, fKey);
2526 //_____________________________________________________________________________
2527 time_t AliSimulation::GenerateTimeStamp() const
2529 // Generate event time-stamp according to
2530 // SOR/EOR time from GRP
2531 if (fUseTimeStampFromCDB)
2532 return fTimeStart + gRandom->Integer(fTimeEnd-fTimeStart);
2537 //_____________________________________________________________________________
2538 void AliSimulation::StoreUsedCDBMaps() const
2540 // write in galice.root maps with used CDB paths
2543 AliRunLoader* runLoader = LoadRun();
2545 AliError("Failed to open gAlice.root in write mode");
2549 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
2550 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
2552 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
2553 cdbMapCopy->SetOwner(1);
2554 // cdbMapCopy->SetName("cdbMap");
2555 TIter iter(cdbMap->GetTable());
2558 while((pair = dynamic_cast<TPair*> (iter.Next()))){
2559 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
2560 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
2561 if (keyStr && valStr)
2562 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
2565 TList *cdbListCopy = new TList();
2566 cdbListCopy->SetOwner(1);
2567 // cdbListCopy->SetName("cdbList");
2569 TIter iter2(cdbList);
2572 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
2573 cdbListCopy->Add(new TObjString(id->ToString().Data()));
2576 AliRunLoader::Instance()->CdGAFile();
2577 gDirectory->WriteObject(cdbMapCopy,"cdbMap","kSingleKey");
2578 gDirectory->WriteObject(cdbListCopy,"cdbList","kSingleKey");
2581 AliInfo(Form("Stored used OCDB entries as TMap %s and TList %s in %s",
2582 cdbMapCopy->GetName(),
2583 cdbListCopy->GetName(),
2584 fGAliceFileName.Data()));