1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running generation, simulation and digitization //
22 // Hits, sdigits and digits are created for all detectors by typing: //
24 // AliSimulation sim; //
27 // The Run method returns kTRUE in case of successful execution. //
28 // The number of events can be given as argument to the Run method or it //
31 // sim.SetNumberOfEvents(n); //
33 // The name of the configuration file can be passed as argument to the //
34 // AliSimulation constructor or can be specified by //
36 // sim.SetConfigFile("..."); //
38 // The generation of particles and the simulation of detector hits can be //
39 // switched on or off by //
41 // sim.SetRunGeneration(kTRUE); // generation of primary particles //
42 // sim.SetRunSimulation(kFALSE); // but no tracking //
44 // For which detectors sdigits and digits will be created, can be steered //
47 // sim.SetMakeSDigits("ALL"); // make sdigits for all detectors //
48 // sim.SetMakeDigits("ITS TPC"); // make digits only for ITS and TPC //
50 // The argument is a (case sensitive) string with the names of the //
51 // detectors separated by a space. An empty string ("") can be used to //
52 // disable the creation of sdigits or digits. The special string "ALL" //
53 // selects all available detectors. This is the default. //
55 // The creation of digits from hits instead of from sdigits can be selected //
58 // sim.SetMakeDigitsFromHits("TRD"); //
60 // The argument is again a string with the selected detectors. Be aware that //
61 // this feature is not available for all detectors and that merging is not //
62 // possible, when digits are created directly from hits. //
64 // Background events can be merged by calling //
66 // sim.MergeWith("background/galice.root", 2); //
68 // The first argument is the file name of the background galice file. The //
69 // second argument is the number of signal events per background event. //
70 // By default this number is calculated from the number of available //
71 // background events. MergeWith can be called several times to merge more //
72 // than two event streams. It is assumed that the sdigits were already //
73 // produced for the background events. //
75 // The output of raw data can be switched on by calling //
77 // sim.SetWriteRawData("MUON"); // write raw data for MUON //
79 // The default output format of the raw data are DDL files. They are //
80 // converted to a DATE file, if a file name is given as second argument. //
81 // For this conversion the program "dateStream" is required. If the file //
82 // name has the extension ".root", the DATE file is converted to a root //
83 // file. The program "alimdc" is used for this purpose. For the conversion //
84 // to DATE and root format the two conversion programs have to be installed. //
85 // Only the raw data in the final format is kept if the third argument is //
88 // The methods RunSimulation, RunSDigitization, RunDigitization, //
89 // RunHitsDigitization and WriteRawData can be used to run only parts of //
90 // the full simulation chain. The creation of raw data DDL files and their //
91 // conversion to the DATE or root format can be run directly by calling //
92 // the methods WriteRawFiles, ConvertRawFilesToDate and ConvertDateToRoot. //
94 // The default number of events per file, which is usually set in the //
95 // config file, can be changed for individual detectors and data types //
98 // sim.SetEventsPerFile("PHOS", "Reconstructed Points", 3); //
100 // The first argument is the detector, the second one the data type and the //
101 // last one the number of events per file. Valid data types are "Hits", //
102 // "Summable Digits", "Digits", "Reconstructed Points" and "Tracks". //
103 // The number of events per file has to be set before the simulation of //
104 // hits. Otherwise it has no effect. //
106 ///////////////////////////////////////////////////////////////////////////////
110 #include <TGeoGlobalMagField.h>
111 #include <TGeoManager.h>
112 #include <TObjString.h>
115 #include <TVirtualMC.h>
116 #include <TVirtualMCApplication.h>
119 #include "AliAlignObj.h"
120 #include "AliCDBEntry.h"
121 #include "AliCDBManager.h"
122 #include "AliCDBStorage.h"
123 #include "AliCTPRawData.h"
124 #include "AliCentralTrigger.h"
125 #include "AliCentralTrigger.h"
126 #include "AliCodeTimer.h"
128 #include "AliDigitizer.h"
129 #include "AliESDEvent.h"
130 #include "AliGRPObject.h"
131 #include "AliGenEventHeader.h"
132 #include "AliGenerator.h"
133 #include "AliGeomManager.h"
134 #include "AliHLTSimulation.h"
135 #include "AliHeader.h"
137 #include "AliLegoGenerator.h"
141 #include "AliModule.h"
143 #include "AliRawReaderDate.h"
144 #include "AliRawReaderFile.h"
145 #include "AliRawReaderRoot.h"
147 #include "AliRunDigitizer.h"
148 #include "AliRunLoader.h"
149 #include "AliSimulation.h"
150 #include "AliSysInfo.h"
151 #include "AliVertexGenFile.h"
153 ClassImp(AliSimulation)
155 AliSimulation *AliSimulation::fgInstance = 0;
156 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
158 //_____________________________________________________________________________
159 AliSimulation::AliSimulation(const char* configFileName,
160 const char* name, const char* title) :
163 fRunGeneration(kTRUE),
164 fRunSimulation(kTRUE),
165 fLoadAlignFromCDB(kTRUE),
166 fLoadAlObjsListOfDets("ALL"),
170 fMakeDigitsFromHits(""),
172 fRawDataFileName(""),
173 fDeleteIntermediateFiles(kFALSE),
174 fWriteSelRawData(kFALSE),
175 fStopOnError(kFALSE),
177 fConfigFileName(configFileName),
178 fGAliceFileName("galice.root"),
180 fBkgrdFileNames(NULL),
181 fAlignObjArray(NULL),
182 fUseBkgrdVertex(kTRUE),
183 fRegionOfInterest(kFALSE),
189 fInitCDBCalled(kFALSE),
190 fInitRunNumberCalled(kFALSE),
191 fSetRunNumberFromDataCalled(kFALSE),
192 fEmbeddingFlag(kFALSE),
197 fEventSpecie(AliRecoParam::kDefault),
198 fWriteQAExpertData(kTRUE),
202 fWriteGRPEntry(kTRUE)
204 // create simulation object with default parameters
206 SetGAliceFile("galice.root");
209 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
210 qam->SetActiveDetectors(fQADetectors) ;
211 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
212 qam->SetTasks(fQATasks) ;
215 //_____________________________________________________________________________
216 AliSimulation::~AliSimulation()
220 fEventsPerFile.Delete();
221 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
222 // delete fAlignObjArray; fAlignObjArray=0;
224 if (fBkgrdFileNames) {
225 fBkgrdFileNames->Delete();
226 delete fBkgrdFileNames;
229 fSpecCDBUri.Delete();
230 if (fgInstance==this) fgInstance = 0;
232 AliQAManager::QAManager()->ShowQA() ;
233 AliQAManager::Destroy() ;
234 AliCodeTimer::Instance()->Print();
238 //_____________________________________________________________________________
239 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
241 // set the number of events for one run
246 //_____________________________________________________________________________
247 void AliSimulation::InitQA()
249 // activate a default CDB storage
250 // First check if we have any CDB storage set, because it is used
251 // to retrieve the calibration and alignment constants
253 if (fInitCDBCalled) return;
254 fInitCDBCalled = kTRUE;
256 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
257 qam->SetActiveDetectors(fQADetectors) ;
258 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
259 qam->SetTasks(fQATasks) ;
260 if (fWriteQAExpertData)
261 qam->SetWriteExpert() ;
263 if (qam->IsDefaultStorageSet()) {
264 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
265 AliWarning("Default QA reference storage has been already set !");
266 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fQARefUri.Data()));
267 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
268 fQARefUri = qam->GetDefaultStorage()->GetURI();
270 if (fQARefUri.Length() > 0) {
271 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
272 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
273 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
275 fQARefUri="local://$ALICE_ROOT/QARef";
276 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
277 AliWarning("Default QA reference storage not yet set !!!!");
278 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
279 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
281 qam->SetDefaultStorage(fQARefUri);
285 //_____________________________________________________________________________
286 void AliSimulation::InitCDB()
288 // activate a default CDB storage
289 // First check if we have any CDB storage set, because it is used
290 // to retrieve the calibration and alignment constants
292 if (fInitCDBCalled) return;
293 fInitCDBCalled = kTRUE;
295 AliCDBManager* man = AliCDBManager::Instance();
296 if (man->IsDefaultStorageSet())
298 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
299 AliWarning("Default CDB storage has been already set !");
300 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
301 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
302 fCDBUri = man->GetDefaultStorage()->GetURI();
305 if (fCDBUri.Length() > 0)
307 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
308 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
309 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
311 fCDBUri="local://$ALICE_ROOT/OCDB";
312 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
313 AliWarning("Default CDB storage not yet set !!!!");
314 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
315 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
318 man->SetDefaultStorage(fCDBUri);
321 // Now activate the detector specific CDB storage locations
322 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
323 TObject* obj = fSpecCDBUri[i];
325 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
326 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
327 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
328 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
333 //_____________________________________________________________________________
334 void AliSimulation::InitRunNumber(){
335 // check run number. If not set, set it to 0 !!!!
337 if (fInitRunNumberCalled) return;
338 fInitRunNumberCalled = kTRUE;
340 AliCDBManager* man = AliCDBManager::Instance();
341 if (man->GetRun() >= 0)
343 AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
344 "Use external variable DC_RUN or AliSimulation::SetRun()!"));
348 AliDebug(2, Form("Setting CDB run number to: %d",fRun));
351 AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
360 //_____________________________________________________________________________
361 void AliSimulation::SetCDBLock() {
362 // Set CDB lock: from now on it is forbidden to reset the run number
363 // or the default storage or to activate any further storage!
365 AliCDBManager::Instance()->SetLock(1);
368 //_____________________________________________________________________________
369 void AliSimulation::SetDefaultStorage(const char* uri) {
370 // Store the desired default CDB storage location
371 // Activate it later within the Run() method
377 //_____________________________________________________________________________
378 void AliSimulation::SetQARefDefaultStorage(const char* uri) {
379 // Store the desired default CDB storage location
380 // Activate it later within the Run() method
383 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
386 //_____________________________________________________________________________
387 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
388 // Store a detector-specific CDB storage location
389 // Activate it later within the Run() method
391 AliCDBPath aPath(calibType);
392 if(!aPath.IsValid()){
393 AliError(Form("Not a valid path: %s", calibType));
397 TObject* obj = fSpecCDBUri.FindObject(calibType);
398 if (obj) fSpecCDBUri.Remove(obj);
399 fSpecCDBUri.Add(new TNamed(calibType, uri));
403 //_____________________________________________________________________________
404 void AliSimulation::SetRunNumber(Int_t run)
407 // Activate it later within the Run() method
412 //_____________________________________________________________________________
413 void AliSimulation::SetSeed(Int_t seed)
416 // Activate it later within the Run() method
421 //_____________________________________________________________________________
422 Bool_t AliSimulation::SetRunNumberFromData()
424 // Set the CDB manager run number
425 // The run number is retrieved from gAlice
427 if (fSetRunNumberFromDataCalled) return kTRUE;
428 fSetRunNumberFromDataCalled = kTRUE;
430 AliCDBManager* man = AliCDBManager::Instance();
431 Int_t runData = -1, runCDB = -1;
433 AliRunLoader* runLoader = LoadRun("READ");
434 if (!runLoader) return kFALSE;
436 runData = runLoader->GetHeader()->GetRun();
440 runCDB = man->GetRun();
442 if (runCDB != runData) {
443 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
444 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
445 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
446 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
451 man->SetRun(runData);
454 if(man->GetRun() < 0) {
455 AliError("Run number not properly initalized!");
464 //_____________________________________________________________________________
465 void AliSimulation::SetConfigFile(const char* fileName)
467 // set the name of the config file
469 fConfigFileName = fileName;
472 //_____________________________________________________________________________
473 void AliSimulation::SetGAliceFile(const char* fileName)
475 // set the name of the galice file
476 // the path is converted to an absolute one if it is relative
478 fGAliceFileName = fileName;
479 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
480 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
482 fGAliceFileName = absFileName;
483 delete[] absFileName;
486 AliDebug(2, Form("galice file name set to %s", fileName));
489 //_____________________________________________________________________________
490 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
493 // set the number of events per file for the given detector and data type
494 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
496 TNamed* obj = new TNamed(detector, type);
497 obj->SetUniqueID(nEvents);
498 fEventsPerFile.Add(obj);
501 //_____________________________________________________________________________
502 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
504 // Read the alignment objects from CDB.
505 // Each detector is supposed to have the
506 // alignment objects in DET/Align/Data CDB path.
507 // All the detector objects are then collected,
508 // sorted by geometry level (starting from ALIC) and
509 // then applied to the TGeo geometry.
510 // Finally an overlaps check is performed.
512 if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
513 AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
517 // initialize CDB storage, run number, set CDB lock
519 // if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
522 Bool_t delRunLoader = kFALSE;
524 runLoader = LoadRun("READ");
525 if (!runLoader) return kFALSE;
526 delRunLoader = kTRUE;
529 // Export ideal geometry
530 if(!IsGeometryFromFile()) AliGeomManager::GetGeometry()->Export("geometry.root");
532 // Load alignment data from CDB and apply to geometry through AliGeomManager
533 if(fLoadAlignFromCDB){
535 TString detStr = fLoadAlObjsListOfDets;
536 TString loadAlObjsListOfDets = "";
538 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
539 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
540 AliModule* det = (AliModule*) detArray->At(iDet);
541 if (!det || !det->IsActive()) continue;
542 if (IsSelected(det->GetName(), detStr)) {
543 //add det to list of dets to be aligned from CDB
544 loadAlObjsListOfDets += det->GetName();
545 loadAlObjsListOfDets += " ";
547 } // end loop over detectors
548 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
549 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
551 // Check if the array with alignment objects was
552 // provided by the user. If yes, apply the objects
553 // to the present TGeo geometry
554 if (fAlignObjArray) {
555 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
556 AliError("The misalignment of one or more volumes failed!"
557 "Compare the list of simulated detectors and the list of detector alignment data!");
558 if (delRunLoader) delete runLoader;
564 // Update the internal geometry of modules (ITS needs it)
565 TString detStr = fLoadAlObjsListOfDets;
566 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
567 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
569 AliModule* det = (AliModule*) detArray->At(iDet);
570 if (!det || !det->IsActive()) continue;
571 if (IsSelected(det->GetName(), detStr)) {
572 det->UpdateInternalGeometry();
574 } // end loop over detectors
577 if (delRunLoader) delete runLoader;
582 //_____________________________________________________________________________
583 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
585 // add a file with background events for merging
587 TObjString* fileNameStr = new TObjString(fileName);
588 fileNameStr->SetUniqueID(nSignalPerBkgrd);
589 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
590 fBkgrdFileNames->Add(fileNameStr);
593 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
595 // add a file with background events for embeddin
596 MergeWith(fileName, nSignalPerBkgrd);
597 fEmbeddingFlag = kTRUE;
600 //_____________________________________________________________________________
601 Bool_t AliSimulation::Run(Int_t nEvents)
603 // run the generation, simulation and digitization
606 AliCodeTimerAuto("",0)
607 AliSysInfo::AddStamp("Start_Run");
609 // Load run number and seed from environmental vars
610 ProcessEnvironmentVars();
611 AliSysInfo::AddStamp("ProcessEnvironmentVars");
613 gRandom->SetSeed(fSeed);
615 if (nEvents > 0) fNEvents = nEvents;
617 // create and setup the HLT instance
618 if (!fRunHLT.IsNull() && !CreateHLT()) {
619 if (fStopOnError) return kFALSE;
624 // generation and simulation -> hits
625 if (fRunGeneration) {
626 if (!RunSimulation()) if (fStopOnError) return kFALSE;
628 AliSysInfo::AddStamp("RunSimulation");
630 // initialize CDB storage from external environment
631 // (either CDB manager or AliSimulation setters),
632 // if not already done in RunSimulation()
634 AliSysInfo::AddStamp("InitCDB");
636 // Set run number in CDBManager from data
637 // From this point on the run number must be always loaded from data!
638 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
640 // Set CDB lock: from now on it is forbidden to reset the run number
641 // or the default storage or to activate any further storage!
644 // If RunSimulation was not called, load the geometry and misalign it
645 if (!AliGeomManager::GetGeometry()) {
646 // Initialize the geometry manager
647 AliGeomManager::LoadGeometry("geometry.root");
648 AliSysInfo::AddStamp("GetGeometry");
651 // // Check that the consistency of symbolic names for the activated subdetectors
652 // // in the geometry loaded by AliGeomManager
653 // AliRunLoader* runLoader = LoadRun("READ");
654 // if (!runLoader) return kFALSE;
656 // TString detsToBeChecked = "";
657 // TObjArray* detArray = runLoader->GetAliRun()->Detectors();
658 // for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
659 // AliModule* det = (AliModule*) detArray->At(iDet);
660 // if (!det || !det->IsActive()) continue;
661 // detsToBeChecked += det->GetName();
662 // detsToBeChecked += " ";
663 // } // end loop over detectors
664 // if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
665 if(!AliGeomManager::CheckSymNamesLUT("ALL"))
666 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
668 if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
670 if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
672 AliSysInfo::AddStamp("MissalignGeometry");
675 // hits -> summable digits
676 AliSysInfo::AddStamp("Start_sdigitization");
677 if (!fMakeSDigits.IsNull()) {
678 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
681 AliSysInfo::AddStamp("Stop_sdigitization");
683 AliSysInfo::AddStamp("Start_digitization");
684 // summable digits -> digits
685 if (!fMakeDigits.IsNull()) {
686 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
687 if (fStopOnError) return kFALSE;
690 AliSysInfo::AddStamp("Stop_digitization");
695 if (!fMakeDigitsFromHits.IsNull()) {
696 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
697 AliWarning(Form("Merging and direct creation of digits from hits "
698 "was selected for some detectors. "
699 "No merging will be done for the following detectors: %s",
700 fMakeDigitsFromHits.Data()));
702 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
703 if (fStopOnError) return kFALSE;
707 AliSysInfo::AddStamp("Hits2Digits");
711 if (!RunTrigger(fTriggerConfig,fMakeDigits)) {
712 if (fStopOnError) return kFALSE;
715 AliSysInfo::AddStamp("RunTrigger");
718 // digits -> raw data
719 if (!fWriteRawData.IsNull()) {
720 if (!WriteRawData(fWriteRawData, fRawDataFileName,
721 fDeleteIntermediateFiles,fWriteSelRawData)) {
722 if (fStopOnError) return kFALSE;
726 AliSysInfo::AddStamp("WriteRaw");
728 // run HLT simulation on simulated digit data if raw data is not
729 // simulated, otherwise its called as part of WriteRawData
730 if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
732 if (fStopOnError) return kFALSE;
736 AliSysInfo::AddStamp("RunHLT");
740 Bool_t rv = RunQA() ;
746 AliSysInfo::AddStamp("RunQA");
748 // Cleanup of CDB manager: cache and active storages!
749 AliCDBManager::Instance()->ClearCache();
754 //_______________________________________________________________________
755 Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
756 Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
757 Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
760 // Generates lego plots of:
761 // - radiation length map phi vs theta
762 // - radiation length map phi vs eta
763 // - interaction length map
764 // - g/cm2 length map
766 // ntheta bins in theta, eta
767 // themin minimum angle in theta (degrees)
768 // themax maximum angle in theta (degrees)
770 // phimin minimum angle in phi (degrees)
771 // phimax maximum angle in phi (degrees)
772 // rmin minimum radius
773 // rmax maximum radius
776 // The number of events generated = ntheta*nphi
777 // run input parameters in macro setup (default="Config.C")
779 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
782 <img src="picts/AliRunLego1.gif">
787 <img src="picts/AliRunLego2.gif">
792 <img src="picts/AliRunLego3.gif">
797 // run the generation and simulation
799 AliCodeTimerAuto("",0)
801 // initialize CDB storage and run number from external environment
802 // (either CDB manager or AliSimulation setters)
808 AliError("no gAlice object. Restart aliroot and try again.");
811 if (gAlice->Modules()->GetEntries() > 0) {
812 AliError("gAlice was already run. Restart aliroot and try again.");
816 AliInfo(Form("initializing gAlice with config file %s",
817 fConfigFileName.Data()));
820 if (nev == -1) nev = nc1 * nc2;
822 // check if initialisation has been done
823 // If runloader has been initialized, set the number of events per file to nc1 * nc2
826 if (!gener) gener = new AliLegoGenerator();
828 // Configure Generator
830 gener->SetRadiusRange(rmin, rmax);
831 gener->SetZMax(zmax);
832 gener->SetCoor1Range(nc1, c1min, c1max);
833 gener->SetCoor2Range(nc2, c2min, c2max);
837 fLego = new AliLego("lego",gener);
839 //__________________________________________________________________________
843 gROOT->LoadMacro(setup);
844 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
846 if(AliCDBManager::Instance()->GetRun() >= 0) {
847 SetRunNumber(AliCDBManager::Instance()->GetRun());
849 AliWarning("Run number not initialized!!");
852 AliRunLoader::Instance()->CdGAFile();
854 AliPDG::AddParticlesToPdgDataBase();
856 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
858 gAlice->GetMCApp()->Init();
861 //Must be here because some MCs (G4) adds detectors here and not in Config.C
862 gAlice->InitLoaders();
863 AliRunLoader::Instance()->MakeTree("E");
866 // Save stuff at the beginning of the file to avoid file corruption
867 AliRunLoader::Instance()->CdGAFile();
870 //Save current generator
871 AliGenerator *gen=gAlice->GetMCApp()->Generator();
872 gAlice->GetMCApp()->ResetGenerator(gener);
873 //Prepare MC for Lego Run
879 AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
880 gMC->ProcessRun(nev);
882 // End of this run, close files
884 // Restore current generator
885 gAlice->GetMCApp()->ResetGenerator(gen);
886 // Delete Lego Object
892 //_____________________________________________________________________________
893 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
897 AliCodeTimerAuto("",0)
899 // initialize CDB storage from external environment
900 // (either CDB manager or AliSimulation setters),
901 // if not already done in RunSimulation()
904 // Set run number in CDBManager from data
905 // From this point on the run number must be always loaded from data!
906 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
908 // Set CDB lock: from now on it is forbidden to reset the run number
909 // or the default storage or to activate any further storage!
912 AliRunLoader* runLoader = LoadRun("READ");
913 if (!runLoader) return kFALSE;
914 TString trconfiguration = config;
916 if (trconfiguration.IsNull()) {
917 if(!fTriggerConfig.IsNull()) {
918 trconfiguration = fTriggerConfig;
921 AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
924 runLoader->MakeTree( "GG" );
925 AliCentralTrigger* aCTP = runLoader->GetTrigger();
926 // Load Configuration
927 if (!aCTP->LoadConfiguration( trconfiguration ))
931 if( !aCTP->RunTrigger( runLoader , detectors ) ) {
943 //_____________________________________________________________________________
944 Bool_t AliSimulation::WriteTriggerRawData()
946 // Writes the CTP (trigger) DDL raw data
947 // Details of the format are given in the
948 // trigger TDR - pages 134 and 135.
949 AliCTPRawData writer;
955 //_____________________________________________________________________________
956 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
958 // run the generation and simulation
960 AliCodeTimerAuto("",0)
962 // initialize CDB storage and run number from external environment
963 // (either CDB manager or AliSimulation setters)
964 AliSysInfo::AddStamp("RunSimulation_Begin");
966 AliSysInfo::AddStamp("RunSimulation_InitCDB");
969 AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
972 AliError("no gAlice object. Restart aliroot and try again.");
975 if (gAlice->Modules()->GetEntries() > 0) {
976 AliError("gAlice was already run. Restart aliroot and try again.");
980 AliInfo(Form("initializing gAlice with config file %s",
981 fConfigFileName.Data()));
984 // Initialize ALICE Simulation run
989 gROOT->LoadMacro(fConfigFileName.Data());
990 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
991 AliSysInfo::AddStamp("RunSimulation_Config");
993 if(AliCDBManager::Instance()->GetRun() >= 0) {
994 AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
995 AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
997 AliWarning("Run number not initialized!!");
1000 AliRunLoader::Instance()->CdGAFile();
1002 AliPDG::AddParticlesToPdgDataBase();
1004 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1005 AliSysInfo::AddStamp("RunSimulation_GetField");
1007 gAlice->GetMCApp()->Init();
1008 AliSysInfo::AddStamp("RunSimulation_InitMCApp");
1010 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1011 gAlice->InitLoaders();
1012 AliRunLoader::Instance()->MakeTree("E");
1013 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1014 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1015 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1017 // Save stuff at the beginning of the file to avoid file corruption
1018 AliRunLoader::Instance()->CdGAFile();
1020 gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1021 AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1022 //___________________________________________________________________________________________
1024 AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
1026 // Set run number in CDBManager
1027 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
1029 AliRunLoader* runLoader = AliRunLoader::Instance();
1031 AliError(Form("gAlice has no run loader object. "
1032 "Check your config file: %s", fConfigFileName.Data()));
1035 SetGAliceFile(runLoader->GetFileName());
1037 // Misalign geometry
1038 #if ROOT_VERSION_CODE < 331527
1039 AliGeomManager::SetGeometry(gGeoManager);
1041 // Check that the consistency of symbolic names for the activated subdetectors
1042 // in the geometry loaded by AliGeomManager
1043 TString detsToBeChecked = "";
1044 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1045 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1046 AliModule* det = (AliModule*) detArray->At(iDet);
1047 if (!det || !det->IsActive()) continue;
1048 detsToBeChecked += det->GetName();
1049 detsToBeChecked += " ";
1050 } // end loop over detectors
1051 if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
1052 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
1053 MisalignGeometry(runLoader);
1054 AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
1057 // AliRunLoader* runLoader = AliRunLoader::Instance();
1058 // if (!runLoader) {
1059 // AliError(Form("gAlice has no run loader object. "
1060 // "Check your config file: %s", fConfigFileName.Data()));
1063 // SetGAliceFile(runLoader->GetFileName());
1065 if (!gAlice->GetMCApp()->Generator()) {
1066 AliError(Form("gAlice has no generator object. "
1067 "Check your config file: %s", fConfigFileName.Data()));
1071 // Write GRP entry corresponding to the setting found in Cofig.C
1074 AliSysInfo::AddStamp("RunSimulation_WriteGRP");
1076 if (nEvents <= 0) nEvents = fNEvents;
1078 // get vertex from background file in case of merging
1079 if (fUseBkgrdVertex &&
1080 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
1081 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
1082 const char* fileName = ((TObjString*)
1083 (fBkgrdFileNames->At(0)))->GetName();
1084 AliInfo(Form("The vertex will be taken from the background "
1085 "file %s with nSignalPerBackground = %d",
1086 fileName, signalPerBkgrd));
1087 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
1088 gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
1091 if (!fRunSimulation) {
1092 gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
1095 // set the number of events per file for given detectors and data types
1096 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
1097 if (!fEventsPerFile[i]) continue;
1098 const char* detName = fEventsPerFile[i]->GetName();
1099 const char* typeName = fEventsPerFile[i]->GetTitle();
1100 TString loaderName(detName);
1101 loaderName += "Loader";
1102 AliLoader* loader = runLoader->GetLoader(loaderName);
1104 AliError(Form("RunSimulation", "no loader for %s found\n"
1105 "Number of events per file not set for %s %s",
1106 detName, typeName, detName));
1109 AliDataLoader* dataLoader =
1110 loader->GetDataLoader(typeName);
1112 AliError(Form("no data loader for %s found\n"
1113 "Number of events per file not set for %s %s",
1114 typeName, detName, typeName));
1117 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
1118 AliDebug(1, Form("number of events per file set to %d for %s %s",
1119 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
1122 AliInfo("running gAlice");
1123 AliSysInfo::AddStamp("Start_ProcessRun");
1125 // Create the Root Tree with one branch per detector
1126 //Hits moved to begin event -> now we are crating separate tree for each event
1128 gMC->ProcessRun(nEvents);
1130 // End of this run, close files
1131 if(nEvents>0) FinishRun();
1133 AliSysInfo::AddStamp("Stop_ProcessRun");
1139 //_____________________________________________________________________________
1140 Bool_t AliSimulation::RunSDigitization(const char* detectors)
1142 // run the digitization and produce summable digits
1143 static Int_t eventNr=0;
1144 AliCodeTimerAuto("",0) ;
1146 // initialize CDB storage, run number, set CDB lock
1148 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1151 AliRunLoader* runLoader = LoadRun();
1152 if (!runLoader) return kFALSE;
1154 TString detStr = detectors;
1155 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1156 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1157 AliModule* det = (AliModule*) detArray->At(iDet);
1158 if (!det || !det->IsActive()) continue;
1159 if (IsSelected(det->GetName(), detStr)) {
1160 AliInfo(Form("creating summable digits for %s", det->GetName()));
1161 AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
1162 det->Hits2SDigits();
1163 AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
1164 AliSysInfo::AddStamp(Form("Digit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
1168 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1169 AliError(Form("the following detectors were not found: %s",
1171 if (fStopOnError) return kFALSE;
1180 //_____________________________________________________________________________
1181 Bool_t AliSimulation::RunDigitization(const char* detectors,
1182 const char* excludeDetectors)
1184 // run the digitization and produce digits from sdigits
1186 AliCodeTimerAuto("",0)
1188 // initialize CDB storage, run number, set CDB lock
1190 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1193 delete AliRunLoader::Instance();
1198 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
1199 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
1200 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
1201 // manager->SetEmbeddingFlag(fEmbeddingFlag);
1202 manager->SetInputStream(0, fGAliceFileName.Data());
1203 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
1204 const char* fileName = ((TObjString*)
1205 (fBkgrdFileNames->At(iStream-1)))->GetName();
1206 manager->SetInputStream(iStream, fileName);
1209 TString detStr = detectors;
1210 TString detExcl = excludeDetectors;
1211 manager->GetInputStream(0)->ImportgAlice();
1212 AliRunLoader* runLoader =
1213 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
1214 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1215 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1216 AliModule* det = (AliModule*) detArray->At(iDet);
1217 if (!det || !det->IsActive()) continue;
1218 if (IsSelected(det->GetName(), detStr) &&
1219 !IsSelected(det->GetName(), detExcl)) {
1220 AliDigitizer* digitizer = det->CreateDigitizer(manager);
1223 AliError(Form("no digitizer for %s", det->GetName()));
1224 if (fStopOnError) return kFALSE;
1226 digitizer->SetRegionOfInterest(fRegionOfInterest);
1231 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1232 AliError(Form("the following detectors were not found: %s",
1234 if (fStopOnError) return kFALSE;
1237 if (!manager->GetListOfTasks()->IsEmpty()) {
1238 AliInfo("executing digitization");
1247 //_____________________________________________________________________________
1248 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
1250 // run the digitization and produce digits from hits
1252 AliCodeTimerAuto("",0)
1254 // initialize CDB storage, run number, set CDB lock
1256 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1259 AliRunLoader* runLoader = LoadRun("READ");
1260 if (!runLoader) return kFALSE;
1262 TString detStr = detectors;
1263 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1264 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1265 AliModule* det = (AliModule*) detArray->At(iDet);
1266 if (!det || !det->IsActive()) continue;
1267 if (IsSelected(det->GetName(), detStr)) {
1268 AliInfo(Form("creating digits from hits for %s", det->GetName()));
1273 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1274 AliError(Form("the following detectors were not found: %s",
1276 if (fStopOnError) return kFALSE;
1282 //_____________________________________________________________________________
1283 Bool_t AliSimulation::WriteRawData(const char* detectors,
1284 const char* fileName,
1285 Bool_t deleteIntermediateFiles,
1288 // convert the digits to raw data
1289 // First DDL raw data files for the given detectors are created.
1290 // If a file name is given, the DDL files are then converted to a DATE file.
1291 // If deleteIntermediateFiles is true, the DDL raw files are deleted
1293 // If the file name has the extension ".root", the DATE file is converted
1295 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1296 // 'selrawdata' flag can be used to enable writing of detectors raw data
1297 // accoring to the trigger cluster.
1299 AliCodeTimerAuto("",0)
1300 AliSysInfo::AddStamp("WriteRawData_Start");
1302 TString detStr = detectors;
1303 if (!WriteRawFiles(detStr.Data())) {
1304 if (fStopOnError) return kFALSE;
1306 AliSysInfo::AddStamp("WriteRawFiles");
1308 // run HLT simulation on simulated DDL raw files
1309 // and produce HLT ddl raw files to be included in date/root file
1310 // bugfix 2009-06-26: the decision whether to write HLT raw data
1311 // is taken in RunHLT. Here HLT always needs to be run in order to
1312 // create HLT digits, unless its switched off. This is due to the
1313 // special placement of the HLT between the generation of DDL files
1314 // and conversion to DATE/Root file.
1315 detStr.ReplaceAll("HLT", "");
1316 if (!fRunHLT.IsNull()) {
1318 if (fStopOnError) return kFALSE;
1321 AliSysInfo::AddStamp("WriteRawData_RunHLT");
1323 TString dateFileName(fileName);
1324 if (!dateFileName.IsNull()) {
1325 Bool_t rootOutput = dateFileName.EndsWith(".root");
1326 if (rootOutput) dateFileName += ".date";
1327 TString selDateFileName;
1329 selDateFileName = "selected.";
1330 selDateFileName+= dateFileName;
1332 if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1333 if (fStopOnError) return kFALSE;
1335 AliSysInfo::AddStamp("ConvertRawFilesToDate");
1336 if (deleteIntermediateFiles) {
1337 AliRunLoader* runLoader = LoadRun("READ");
1338 if (runLoader) for (Int_t iEvent = 0;
1339 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1341 sprintf(command, "rm -r raw%d", iEvent);
1342 gSystem->Exec(command);
1348 if (!ConvertDateToRoot(dateFileName, fileName)) {
1349 if (fStopOnError) return kFALSE;
1351 AliSysInfo::AddStamp("ConvertDateToRoot");
1352 if (deleteIntermediateFiles) {
1353 gSystem->Unlink(dateFileName);
1356 TString selFileName = "selected.";
1357 selFileName += fileName;
1358 if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1359 if (fStopOnError) return kFALSE;
1361 if (deleteIntermediateFiles) {
1362 gSystem->Unlink(selDateFileName);
1371 //_____________________________________________________________________________
1372 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1374 // convert the digits to raw data DDL files
1376 AliCodeTimerAuto("",0)
1378 AliRunLoader* runLoader = LoadRun("READ");
1379 if (!runLoader) return kFALSE;
1381 // write raw data to DDL files
1382 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1383 AliInfo(Form("processing event %d", iEvent));
1384 runLoader->GetEvent(iEvent);
1385 TString baseDir = gSystem->WorkingDirectory();
1387 sprintf(dirName, "raw%d", iEvent);
1388 gSystem->MakeDirectory(dirName);
1389 if (!gSystem->ChangeDirectory(dirName)) {
1390 AliError(Form("couldn't change to directory %s", dirName));
1391 if (fStopOnError) return kFALSE; else continue;
1394 ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1397 TString detStr = detectors;
1398 if (IsSelected("HLT", detStr)) {
1399 // Do nothing. "HLT" will be removed from detStr and HLT raw
1400 // data files are generated in RunHLT.
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)) {
1408 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1413 if (!WriteTriggerRawData())
1414 if (fStopOnError) return kFALSE;
1416 gSystem->ChangeDirectory(baseDir);
1417 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1418 AliError(Form("the following detectors were not found: %s",
1420 if (fStopOnError) return kFALSE;
1429 //_____________________________________________________________________________
1430 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1431 const char* selDateFileName)
1433 // convert raw data DDL files to a DATE file with the program "dateStream"
1434 // The second argument is not empty when the user decides to write
1435 // the detectors raw data according to the trigger cluster.
1437 AliCodeTimerAuto("",0)
1439 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1441 AliError("the program dateStream was not found");
1442 if (fStopOnError) return kFALSE;
1447 AliRunLoader* runLoader = LoadRun("READ");
1448 if (!runLoader) return kFALSE;
1450 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1451 Bool_t selrawdata = kFALSE;
1452 if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1455 // Note the option -s. It is used in order to avoid
1456 // the generation of SOR/EOR events.
1457 sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d",
1458 dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1459 FILE* pipe = gSystem->OpenPipe(command, "w");
1462 AliError(Form("Cannot execute command: %s",command));
1466 Int_t selEvents = 0;
1467 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1469 UInt_t detectorPattern = 0;
1470 runLoader->GetEvent(iEvent);
1471 if (!runLoader->LoadTrigger()) {
1472 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1473 detectorPattern = aCTP->GetClusterMask();
1474 // Check if the event was triggered by CTP
1476 if (aCTP->GetClassMask()) selEvents++;
1480 AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1482 AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1483 selrawdata = kFALSE;
1487 fprintf(pipe, "GDC DetectorPattern %u\n", detectorPattern);
1491 // loop over detectors and DDLs
1492 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1493 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1495 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1496 Int_t ldcID = Int_t(ldc + 0.0001);
1497 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1499 char rawFileName[256];
1500 sprintf(rawFileName, "raw%d/%s",
1501 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1503 // check existence and size of raw data file
1504 FILE* file = fopen(rawFileName, "rb");
1505 if (!file) continue;
1506 fseek(file, 0, SEEK_END);
1507 unsigned long size = ftell(file);
1509 if (!size) continue;
1511 if (ldcID != prevLDC) {
1512 fprintf(pipe, " LDC Id %d\n", ldcID);
1515 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1520 Int_t result = gSystem->ClosePipe(pipe);
1522 if (!(selrawdata && selEvents > 0)) {
1524 return (result == 0);
1527 AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1529 sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d",
1530 selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1531 FILE* pipe2 = gSystem->OpenPipe(command, "w");
1533 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1535 // Get the trigger decision and cluster
1536 UInt_t detectorPattern = 0;
1538 runLoader->GetEvent(iEvent);
1539 if (!runLoader->LoadTrigger()) {
1540 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1541 if (aCTP->GetClassMask() == 0) continue;
1542 detectorPattern = aCTP->GetClusterMask();
1543 detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1544 AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1547 fprintf(pipe2, "GDC DetectorPattern %u\n", detectorPattern);
1551 // loop over detectors and DDLs
1552 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1553 // Write only raw data from detectors that
1554 // are contained in the trigger cluster(s)
1555 if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1557 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1559 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1560 Int_t ldcID = Int_t(ldc + 0.0001);
1561 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1563 char rawFileName[256];
1564 sprintf(rawFileName, "raw%d/%s",
1565 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1567 // check existence and size of raw data file
1568 FILE* file = fopen(rawFileName, "rb");
1569 if (!file) continue;
1570 fseek(file, 0, SEEK_END);
1571 unsigned long size = ftell(file);
1573 if (!size) continue;
1575 if (ldcID != prevLDC) {
1576 fprintf(pipe2, " LDC Id %d\n", ldcID);
1579 fprintf(pipe2, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1584 Int_t result2 = gSystem->ClosePipe(pipe2);
1587 return ((result == 0) && (result2 == 0));
1590 //_____________________________________________________________________________
1591 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1592 const char* rootFileName)
1594 // convert a DATE file to a root file with the program "alimdc"
1597 const Int_t kDBSize = 2000000000;
1598 const Int_t kTagDBSize = 1000000000;
1599 const Bool_t kFilter = kFALSE;
1600 const Int_t kCompression = 1;
1602 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1604 AliError("the program alimdc was not found");
1605 if (fStopOnError) return kFALSE;
1610 AliInfo(Form("converting DATE file %s to root file %s",
1611 dateFileName, rootFileName));
1613 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1614 const char* tagDBFS = "/tmp/mdc1/tags";
1616 // User defined file system locations
1617 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1618 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1619 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1620 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1621 if (gSystem->Getenv("ALIMDC_TAGDB"))
1622 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1624 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1625 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1626 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1628 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1629 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1630 gSystem->Exec(Form("mkdir %s",tagDBFS));
1632 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1633 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1634 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1636 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1637 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1638 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1640 return (result == 0);
1644 //_____________________________________________________________________________
1645 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1647 // delete existing run loaders, open a new one and load gAlice
1649 delete AliRunLoader::Instance();
1650 AliRunLoader* runLoader =
1651 AliRunLoader::Open(fGAliceFileName.Data(),
1652 AliConfig::GetDefaultEventFolderName(), mode);
1654 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1657 runLoader->LoadgAlice();
1658 runLoader->LoadHeader();
1659 gAlice = runLoader->GetAliRun();
1661 AliError(Form("no gAlice object found in file %s",
1662 fGAliceFileName.Data()));
1668 //_____________________________________________________________________________
1669 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1671 // get or calculate the number of signal events per background event
1673 if (!fBkgrdFileNames) return 1;
1674 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1675 if (nBkgrdFiles == 0) return 1;
1677 // get the number of signal events
1679 AliRunLoader* runLoader =
1680 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1681 if (!runLoader) return 1;
1683 nEvents = runLoader->GetNumberOfEvents();
1688 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1689 // get the number of background events
1690 const char* fileName = ((TObjString*)
1691 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1692 AliRunLoader* runLoader =
1693 AliRunLoader::Open(fileName, "BKGRD");
1694 if (!runLoader) continue;
1695 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1698 // get or calculate the number of signal per background events
1699 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1700 if (nSignalPerBkgrd <= 0) {
1701 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1702 } else if (result && (result != nSignalPerBkgrd)) {
1703 AliInfo(Form("the number of signal events per background event "
1704 "will be changed from %d to %d for stream %d",
1705 nSignalPerBkgrd, result, iBkgrdFile+1));
1706 nSignalPerBkgrd = result;
1709 if (!result) result = nSignalPerBkgrd;
1710 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1711 AliWarning(Form("not enough background events (%d) for %d signal events "
1712 "using %d signal per background events for stream %d",
1713 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1720 //_____________________________________________________________________________
1721 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1723 // check whether detName is contained in detectors
1724 // if yes, it is removed from detectors
1726 // check if all detectors are selected
1727 if ((detectors.CompareTo("ALL") == 0) ||
1728 detectors.BeginsWith("ALL ") ||
1729 detectors.EndsWith(" ALL") ||
1730 detectors.Contains(" ALL ")) {
1735 // search for the given detector
1736 Bool_t result = kFALSE;
1737 if ((detectors.CompareTo(detName) == 0) ||
1738 detectors.BeginsWith(detName+" ") ||
1739 detectors.EndsWith(" "+detName) ||
1740 detectors.Contains(" "+detName+" ")) {
1741 detectors.ReplaceAll(detName, "");
1745 // clean up the detectors string
1746 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1747 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1748 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1753 //_____________________________________________________________________________
1754 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName)
1757 // Steering routine to convert raw data in directory rawDirectory/ to fake SDigits.
1758 // These can be used for embedding of MC tracks into RAW data using the standard
1759 // merging procedure.
1761 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1764 AliError("no gAlice object. Restart aliroot and try again.");
1767 if (gAlice->Modules()->GetEntries() > 0) {
1768 AliError("gAlice was already run. Restart aliroot and try again.");
1772 AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1776 gROOT->LoadMacro(fConfigFileName.Data());
1777 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1779 if(AliCDBManager::Instance()->GetRun() >= 0) {
1780 SetRunNumber(AliCDBManager::Instance()->GetRun());
1782 AliWarning("Run number not initialized!!");
1785 AliRunLoader::Instance()->CdGAFile();
1787 AliPDG::AddParticlesToPdgDataBase();
1789 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1791 gAlice->GetMCApp()->Init();
1793 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1794 gAlice->InitLoaders();
1795 AliRunLoader::Instance()->MakeTree("E");
1796 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1797 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1798 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1800 // Save stuff at the beginning of the file to avoid file corruption
1801 AliRunLoader::Instance()->CdGAFile();
1806 //AliCDBManager* man = AliCDBManager::Instance();
1807 //man->SetRun(0); // Should this come from rawdata header ?
1811 // Get the runloader
1812 AliRunLoader* runLoader = AliRunLoader::Instance();
1814 // Open esd file if available
1815 TFile* esdFile = TFile::Open(esdFileName);
1817 AliESDEvent* esd = new AliESDEvent();
1818 esdFile->GetObject("esdTree", treeESD);
1819 if (treeESD) esd->ReadFromTree(treeESD);
1822 // Create the RawReader
1823 TString fileName(rawDirectory);
1824 AliRawReader* rawReader = 0x0;
1825 if (fileName.EndsWith("/")) {
1826 rawReader = new AliRawReaderFile(fileName);
1827 } else if (fileName.EndsWith(".root")) {
1828 rawReader = new AliRawReaderRoot(fileName);
1829 } else if (!fileName.IsNull()) {
1830 rawReader = new AliRawReaderDate(fileName);
1832 // if (!fEquipIdMap.IsNull() && fRawReader)
1833 // fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1835 // Get list of detectors
1836 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1839 AliHeader* header = runLoader->GetHeader();
1841 TString detStr = fMakeSDigits;
1845 if (!(rawReader->NextEvent())) break;
1848 for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1849 AliModule* det = (AliModule*) detArray->At(iDet);
1850 if (!det || !det->IsActive()) continue;
1851 if (IsSelected(det->GetName(), detStr)) {
1852 AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1853 det->Raw2SDigits(rawReader);
1860 // If ESD information available obtain reconstructed vertex and store in header.
1862 treeESD->GetEvent(nev);
1863 const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1864 Double_t position[3];
1865 esdVertex->GetXYZ(position);
1866 AliGenEventHeader* mcHeader = new AliGenEventHeader("ESD");
1869 for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1870 mcHeader->SetPrimaryVertex(mcV);
1871 header->Reset(0,nev);
1872 header->SetGenEventHeader(mcHeader);
1873 printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1878 runLoader->TreeE()->Fill();
1879 runLoader->SetNextEvent();
1885 runLoader->CdGAFile();
1886 runLoader->WriteHeader("OVERWRITE");
1887 runLoader->WriteRunLoader();
1892 //_____________________________________________________________________________
1893 void AliSimulation::FinishRun()
1896 // Called at the end of the run.
1901 AliDebug(1, "Finish Lego");
1902 AliRunLoader::Instance()->CdGAFile();
1906 // Clean detector information
1907 TIter next(gAlice->Modules());
1908 AliModule *detector;
1909 while((detector = dynamic_cast<AliModule*>(next()))) {
1910 AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
1911 detector->FinishRun();
1914 AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
1915 AliRunLoader::Instance()->WriteHeader("OVERWRITE");
1917 // Write AliRun info and all detectors parameters
1918 AliRunLoader::Instance()->CdGAFile();
1919 gAlice->Write(0,TObject::kOverwrite);//write AliRun
1920 AliRunLoader::Instance()->Write(0,TObject::kOverwrite);//write RunLoader itself
1922 if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();
1923 AliRunLoader::Instance()->Synchronize();
1926 //_____________________________________________________________________________
1927 Int_t AliSimulation::GetDetIndex(const char* detector)
1929 // return the detector index corresponding to detector
1931 for (index = 0; index < fgkNDetectors ; index++) {
1932 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
1938 //_____________________________________________________________________________
1939 Bool_t AliSimulation::CreateHLT()
1941 // Init the HLT simulation.
1942 // The function loads the library and creates the instance of AliHLTSimulation.
1943 // the main reason for the decoupled creation is to set the transient OCDB
1944 // objects before the OCDB is locked
1946 // load the library dynamically
1947 gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1949 // check for the library version
1950 AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
1952 AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
1955 if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
1956 AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
1959 // print compile info
1960 typedef void (*CompileInfo)( const char*& date, const char*& time);
1961 CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
1963 const char* date="";
1964 const char* time="";
1965 (*fctInfo)(date, time);
1966 if (!date) date="unknown";
1967 if (!time) time="unknown";
1968 AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
1970 AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
1973 // create instance of the HLT simulation
1974 AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
1975 if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
1976 AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
1980 TString specObjects;
1981 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
1982 if (specObjects.Length()>0) specObjects+=" ";
1983 specObjects+=fSpecCDBUri[i]->GetName();
1986 AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
1987 if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
1988 AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
1994 //_____________________________________________________________________________
1995 Bool_t AliSimulation::RunHLT()
1997 // Run the HLT simulation
1998 // HLT simulation is implemented in HLT/sim/AliHLTSimulation
1999 // Disabled if fRunHLT is empty, default vaule is "default".
2000 // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
2001 // The default simulation depends on the HLT component libraries and their
2002 // corresponding agents which define components and chains to run. See
2003 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
2004 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
2006 // The libraries to be loaded can be specified as an option.
2008 // AliSimulation sim;
2009 // sim.SetRunHLT("libAliHLTSample.so");
2011 // will only load <tt>libAliHLTSample.so</tt>
2013 // Other available options:
2014 // \li loglevel=<i>level</i> <br>
2015 // logging level for this processing
2017 // disable redirection of log messages to AliLog class
2018 // \li config=<i>macro</i>
2019 // configuration macro
2020 // \li chains=<i>configuration</i>
2021 // comma separated list of configurations to be run during simulation
2022 // \li rawfile=<i>file</i>
2023 // source for the RawReader to be created, the default is <i>./</i> if
2024 // raw data is simulated
2028 if (!fpHLT && !CreateHLT()) {
2031 AliHLTSimulation* pHLT=fpHLT;
2033 AliRunLoader* pRunLoader = LoadRun("READ");
2034 if (!pRunLoader) return kFALSE;
2036 // initialize CDB storage, run number, set CDB lock
2037 // thats for the case of running HLT simulation without all the other steps
2038 // multiple calls are handled by the function, so we can just call
2040 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
2043 // init the HLT simulation
2045 if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2046 TString detStr = fWriteRawData;
2047 if (!IsSelected("HLT", detStr)) {
2048 options+=" writerawfiles=";
2050 options+=" writerawfiles=HLT";
2053 if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2054 // as a matter of fact, HLT will run reconstruction and needs the RawReader
2055 // in order to get detector data. By default, RawReaderFile is used to read
2056 // the already simulated ddl files. Date and Root files from the raw data
2057 // are generated after the HLT simulation.
2058 options+=" rawfile=./";
2061 AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2062 if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2063 AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2065 // run the HLT simulation
2066 AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2067 if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2068 AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2072 // delete the instance
2073 AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2074 if (fctDelete==NULL || fctDelete(pHLT)<0) {
2075 AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2079 return iResult>=0?kTRUE:kFALSE;
2082 //_____________________________________________________________________________
2083 Bool_t AliSimulation::RunQA()
2085 // run the QA on summable hits, digits or digits
2087 if(!gAlice) return kFALSE;
2088 AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2090 TString detectorsw("") ;
2092 AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2093 detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ;
2094 if ( detectorsw.IsNull() )
2099 //_____________________________________________________________________________
2100 Bool_t AliSimulation::SetRunQA(TString detAndAction)
2102 // Allows to run QA for a selected set of detectors
2103 // and a selected set of tasks among HITS, SDIGITS and DIGITS
2104 // all selected detectors run the same selected tasks
2106 if (!detAndAction.Contains(":")) {
2107 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2111 Int_t colon = detAndAction.Index(":") ;
2112 fQADetectors = detAndAction(0, colon) ;
2113 if (fQADetectors.Contains("ALL") ){
2114 TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2115 Int_t minus = fQADetectors.Last('-') ;
2116 TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2117 TString toRemove("") ;
2118 while (minus >= 0) {
2119 toRemove = fQADetectors(minus+1, fQADetectors.Length()) ;
2120 toRemove = toRemove.Strip() ;
2121 toKeep.ReplaceAll(toRemove, "") ;
2122 fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ;
2123 minus = fQADetectors.Last('-') ;
2125 fQADetectors = toKeep ;
2127 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2128 if (fQATasks.Contains("ALL") ) {
2129 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
2131 fQATasks.ToUpper() ;
2133 if ( fQATasks.Contains("HIT") )
2134 tempo = Form("%d ", AliQAv1::kHITS) ;
2135 if ( fQATasks.Contains("SDIGIT") )
2136 tempo += Form("%d ", AliQAv1::kSDIGITS) ;
2137 if ( fQATasks.Contains("DIGIT") )
2138 tempo += Form("%d ", AliQAv1::kDIGITS) ;
2140 if (fQATasks.IsNull()) {
2141 AliInfo("No QA requested\n") ;
2146 TString tempo(fQATasks) ;
2147 tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS)) ;
2148 tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;
2149 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;
2150 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2152 AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ;
2153 AliQAManager::QAManager()->SetTasks(fQATasks) ;
2154 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++)
2155 AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2160 //_____________________________________________________________________________
2161 void AliSimulation::ProcessEnvironmentVars()
2163 // Extract run number and random generator seed from env variables
2165 AliInfo("Processing environment variables");
2167 // Random Number seed
2169 // first check that seed is not already set
2171 if (gSystem->Getenv("CONFIG_SEED")) {
2172 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2175 if (gSystem->Getenv("CONFIG_SEED")) {
2176 AliInfo(Form("Seed for random number generation already set (%d)"
2177 ": CONFIG_SEED variable ignored!", fSeed));
2181 AliInfo(Form("Seed for random number generation = %d ", fSeed));
2185 // first check that run number is not already set
2187 if (gSystem->Getenv("DC_RUN")) {
2188 fRun = atoi(gSystem->Getenv("DC_RUN"));
2191 if (gSystem->Getenv("DC_RUN")) {
2192 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2196 AliInfo(Form("Run number = %d", fRun));
2199 //---------------------------------------------------------------------
2200 void AliSimulation::WriteGRPEntry()
2202 // Get the necessary information from galice (generator, trigger etc) and
2203 // write a GRP entry corresponding to the settings in the Config.C used
2204 // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2207 AliInfo("Writing global run parameters entry into the OCDB");
2209 AliGRPObject* grpObj = new AliGRPObject();
2211 grpObj->SetRunType("PHYSICS");
2212 grpObj->SetTimeStart(0);
2214 grpObj->SetTimeStart(0);
2215 grpObj->SetTimeEnd(curtime.Convert());
2216 grpObj->SetBeamEnergyIsSqrtSHalfGeV(); // new format of GRP: store sqrt(s)/2 in GeV
2218 const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2220 grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
2223 gen->GetProjectile(projectile,a,z);
2225 gen->GetTarget(target,a,z);
2226 TString beamType = projectile + "-" + target;
2227 beamType.ReplaceAll(" ","");
2228 if (!beamType.CompareTo("-")) {
2229 grpObj->SetBeamType("UNKNOWN");
2232 grpObj->SetBeamType(beamType);
2233 // Heavy ion run, the event specie is set to kHighMult
2234 fEventSpecie = AliRecoParam::kHighMult;
2235 if ((strcmp(beamType,"p-p") == 0) ||
2236 (strcmp(beamType,"p-") == 0) ||
2237 (strcmp(beamType,"-p") == 0) ||
2238 (strcmp(beamType,"P-P") == 0) ||
2239 (strcmp(beamType,"P-") == 0) ||
2240 (strcmp(beamType,"-P") == 0)) {
2241 // Proton run, the event specie is set to kLowMult
2242 fEventSpecie = AliRecoParam::kLowMult;
2246 AliWarning("Unknown beam type and energy! Setting energy to 0");
2247 grpObj->SetBeamEnergy(0);
2248 grpObj->SetBeamType("UNKNOWN");
2251 UInt_t detectorPattern = 0;
2253 TObjArray *detArray = gAlice->Detectors();
2254 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2255 if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2256 detectorPattern |= (1 << iDet);
2261 if (!fTriggerConfig.IsNull())
2262 detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2265 if (!fRunHLT.IsNull())
2266 detectorPattern |= (1 << AliDAQ::kHLTId);
2268 grpObj->SetNumberOfDetectors((Char_t)nDets);
2269 grpObj->SetDetectorMask((Int_t)detectorPattern);
2270 grpObj->SetLHCPeriod("LHC08c");
2271 grpObj->SetLHCState("STABLE_BEAMS");
2273 AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2274 Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2276 Float_t factorSol = field ? field->GetFactorSol() : 0;
2277 Float_t currentSol = TMath::Abs(factorSol)>1E-6 ?
2278 TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
2280 Float_t factorDip = field ? field->GetFactorDip() : 0;
2281 Float_t currentDip = 6000.*TMath::Abs(factorDip);
2283 grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
2284 grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);
2285 grpObj->SetL3Polarity(factorSol>0 ? 0:1);
2286 grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
2287 grpObj->SetUniformBMap(field->IsUniform()); // for special MC with k5kGUniform map
2288 grpObj->SetPolarityConventionLHC(); // LHC convention +/+ current -> -/- field main components
2290 grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2292 //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2294 // Now store the entry in OCDB
2295 AliCDBManager* man = AliCDBManager::Instance();
2297 AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2298 AliCDBMetaData *metadata= new AliCDBMetaData();
2300 metadata->SetResponsible("alice-off@cern.ch");
2301 metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2303 man->Put(grpObj,id,metadata);