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>
118 #include "AliAlignObj.h"
119 #include "AliCDBEntry.h"
120 #include "AliCDBManager.h"
121 #include "AliCDBStorage.h"
122 #include "AliCTPRawData.h"
123 #include "AliCentralTrigger.h"
124 #include "AliCentralTrigger.h"
125 #include "AliCodeTimer.h"
127 #include "AliDigitizer.h"
129 #include "AliGRPObject.h"
130 #include "AliGenEventHeader.h"
131 #include "AliGenerator.h"
132 #include "AliGeomManager.h"
133 #include "AliHLTSimulation.h"
134 #include "AliHeader.h"
136 #include "AliLegoGenerator.h"
140 #include "AliModule.h"
142 #include "AliRawReaderDate.h"
143 #include "AliRawReaderFile.h"
144 #include "AliRawReaderRoot.h"
146 #include "AliRunDigitizer.h"
147 #include "AliRunLoader.h"
148 #include "AliSimulation.h"
149 #include "AliSysInfo.h"
150 #include "AliVertexGenFile.h"
152 ClassImp(AliSimulation)
154 AliSimulation *AliSimulation::fgInstance = 0;
155 const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
157 //_____________________________________________________________________________
158 AliSimulation::AliSimulation(const char* configFileName,
159 const char* name, const char* title) :
162 fRunGeneration(kTRUE),
163 fRunSimulation(kTRUE),
164 fLoadAlignFromCDB(kTRUE),
165 fLoadAlObjsListOfDets("ALL"),
169 fMakeDigitsFromHits(""),
171 fRawDataFileName(""),
172 fDeleteIntermediateFiles(kFALSE),
173 fWriteSelRawData(kFALSE),
174 fStopOnError(kFALSE),
176 fConfigFileName(configFileName),
177 fGAliceFileName("galice.root"),
179 fBkgrdFileNames(NULL),
180 fAlignObjArray(NULL),
181 fUseBkgrdVertex(kTRUE),
182 fRegionOfInterest(kFALSE),
188 fInitCDBCalled(kFALSE),
189 fInitRunNumberCalled(kFALSE),
190 fSetRunNumberFromDataCalled(kFALSE),
191 fEmbeddingFlag(kFALSE),
196 fEventSpecie(AliRecoParam::kDefault),
197 fWriteQAExpertData(kTRUE),
199 fWriteGRPEntry(kTRUE)
201 // create simulation object with default parameters
203 SetGAliceFile("galice.root");
206 AliQAManager * qam = AliQAManager::QAManager("sim") ;
207 qam->SetActiveDetectors(fQADetectors) ;
208 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
209 qam->SetTasks(fQATasks) ;
212 //_____________________________________________________________________________
213 AliSimulation::~AliSimulation()
217 fEventsPerFile.Delete();
218 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
219 // delete fAlignObjArray; fAlignObjArray=0;
221 if (fBkgrdFileNames) {
222 fBkgrdFileNames->Delete();
223 delete fBkgrdFileNames;
226 fSpecCDBUri.Delete();
227 if (fgInstance==this) fgInstance = 0;
229 AliQAManager::QAManager()->ShowQA() ;
230 AliQAManager::Destroy() ;
231 AliCodeTimer::Instance()->Print();
235 //_____________________________________________________________________________
236 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
238 // set the number of events for one run
243 //_____________________________________________________________________________
244 void AliSimulation::InitQA()
246 // activate a default CDB storage
247 // First check if we have any CDB storage set, because it is used
248 // to retrieve the calibration and alignment constants
250 if (fInitCDBCalled) return;
251 fInitCDBCalled = kTRUE;
253 AliQAManager * qam = AliQAManager::QAManager("sim") ;
254 qam->SetActiveDetectors(fQADetectors) ;
255 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
256 qam->SetTasks(fQATasks) ;
257 if (fWriteQAExpertData)
258 qam->SetWriteExpert() ;
260 if (qam->IsDefaultStorageSet()) {
261 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
262 AliWarning("Default QA reference storage has been already set !");
263 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fQARefUri.Data()));
264 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
265 fQARefUri = qam->GetDefaultStorage()->GetURI();
267 if (fQARefUri.Length() > 0) {
268 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
269 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
270 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
272 fQARefUri="local://$ALICE_ROOT/QARef";
273 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
274 AliWarning("Default QA reference storage not yet set !!!!");
275 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
276 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
278 qam->SetDefaultStorage(fQARefUri);
282 //_____________________________________________________________________________
283 void AliSimulation::InitCDB()
285 // activate a default CDB storage
286 // First check if we have any CDB storage set, because it is used
287 // to retrieve the calibration and alignment constants
289 if (fInitCDBCalled) return;
290 fInitCDBCalled = kTRUE;
292 AliCDBManager* man = AliCDBManager::Instance();
293 if (man->IsDefaultStorageSet())
295 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
296 AliWarning("Default CDB storage has been already set !");
297 AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
298 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
299 fCDBUri = man->GetDefaultStorage()->GetURI();
302 if (fCDBUri.Length() > 0)
304 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
305 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
306 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
308 fCDBUri="local://$ALICE_ROOT/OCDB";
309 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
310 AliWarning("Default CDB storage not yet set !!!!");
311 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
312 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
315 man->SetDefaultStorage(fCDBUri);
318 // Now activate the detector specific CDB storage locations
319 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
320 TObject* obj = fSpecCDBUri[i];
322 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
323 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
324 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
325 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
330 //_____________________________________________________________________________
331 void AliSimulation::InitRunNumber(){
332 // check run number. If not set, set it to 0 !!!!
334 if (fInitRunNumberCalled) return;
335 fInitRunNumberCalled = kTRUE;
337 AliCDBManager* man = AliCDBManager::Instance();
338 if (man->GetRun() >= 0)
340 AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
341 "Use external variable DC_RUN or AliSimulation::SetRun()!"));
345 AliDebug(2, Form("Setting CDB run number to: %d",fRun));
348 AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
357 //_____________________________________________________________________________
358 void AliSimulation::SetCDBLock() {
359 // Set CDB lock: from now on it is forbidden to reset the run number
360 // or the default storage or to activate any further storage!
362 AliCDBManager::Instance()->SetLock(1);
365 //_____________________________________________________________________________
366 void AliSimulation::SetDefaultStorage(const char* uri) {
367 // Store the desired default CDB storage location
368 // Activate it later within the Run() method
374 //_____________________________________________________________________________
375 void AliSimulation::SetQARefDefaultStorage(const char* uri) {
376 // Store the desired default CDB storage location
377 // Activate it later within the Run() method
380 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
383 //_____________________________________________________________________________
384 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
385 // Store a detector-specific CDB storage location
386 // Activate it later within the Run() method
388 AliCDBPath aPath(calibType);
389 if(!aPath.IsValid()){
390 AliError(Form("Not a valid path: %s", calibType));
394 TObject* obj = fSpecCDBUri.FindObject(calibType);
395 if (obj) fSpecCDBUri.Remove(obj);
396 fSpecCDBUri.Add(new TNamed(calibType, uri));
400 //_____________________________________________________________________________
401 void AliSimulation::SetRunNumber(Int_t run)
404 // Activate it later within the Run() method
409 //_____________________________________________________________________________
410 void AliSimulation::SetSeed(Int_t seed)
413 // Activate it later within the Run() method
418 //_____________________________________________________________________________
419 Bool_t AliSimulation::SetRunNumberFromData()
421 // Set the CDB manager run number
422 // The run number is retrieved from gAlice
424 if (fSetRunNumberFromDataCalled) return kTRUE;
425 fSetRunNumberFromDataCalled = kTRUE;
427 AliCDBManager* man = AliCDBManager::Instance();
428 Int_t runData = -1, runCDB = -1;
430 AliRunLoader* runLoader = LoadRun("READ");
431 if (!runLoader) return kFALSE;
433 runData = runLoader->GetHeader()->GetRun();
437 runCDB = man->GetRun();
439 if (runCDB != runData) {
440 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
441 AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
442 AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
443 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
448 man->SetRun(runData);
451 if(man->GetRun() < 0) {
452 AliError("Run number not properly initalized!");
461 //_____________________________________________________________________________
462 void AliSimulation::SetConfigFile(const char* fileName)
464 // set the name of the config file
466 fConfigFileName = fileName;
469 //_____________________________________________________________________________
470 void AliSimulation::SetGAliceFile(const char* fileName)
472 // set the name of the galice file
473 // the path is converted to an absolute one if it is relative
475 fGAliceFileName = fileName;
476 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
477 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
479 fGAliceFileName = absFileName;
480 delete[] absFileName;
483 AliDebug(2, Form("galice file name set to %s", fileName));
486 //_____________________________________________________________________________
487 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
490 // set the number of events per file for the given detector and data type
491 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
493 TNamed* obj = new TNamed(detector, type);
494 obj->SetUniqueID(nEvents);
495 fEventsPerFile.Add(obj);
498 //_____________________________________________________________________________
499 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
501 // Read the alignment objects from CDB.
502 // Each detector is supposed to have the
503 // alignment objects in DET/Align/Data CDB path.
504 // All the detector objects are then collected,
505 // sorted by geometry level (starting from ALIC) and
506 // then applied to the TGeo geometry.
507 // Finally an overlaps check is performed.
509 if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
510 AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
514 // initialize CDB storage, run number, set CDB lock
516 // if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
519 Bool_t delRunLoader = kFALSE;
521 runLoader = LoadRun("READ");
522 if (!runLoader) return kFALSE;
523 delRunLoader = kTRUE;
526 // Export ideal geometry
527 if(!gAlice->IsRootGeometry()) AliGeomManager::GetGeometry()->Export("geometry.root");
529 // Load alignment data from CDB and apply to geometry through AliGeomManager
530 if(fLoadAlignFromCDB){
532 TString detStr = fLoadAlObjsListOfDets;
533 TString loadAlObjsListOfDets = "";
535 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
536 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
537 AliModule* det = (AliModule*) detArray->At(iDet);
538 if (!det || !det->IsActive()) continue;
539 if (IsSelected(det->GetName(), detStr)) {
540 //add det to list of dets to be aligned from CDB
541 loadAlObjsListOfDets += det->GetName();
542 loadAlObjsListOfDets += " ";
544 } // end loop over detectors
545 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
546 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
548 // Check if the array with alignment objects was
549 // provided by the user. If yes, apply the objects
550 // to the present TGeo geometry
551 if (fAlignObjArray) {
552 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
553 AliError("The misalignment of one or more volumes failed!"
554 "Compare the list of simulated detectors and the list of detector alignment data!");
555 if (delRunLoader) delete runLoader;
561 // Update the internal geometry of modules (ITS needs it)
562 TString detStr = fLoadAlObjsListOfDets;
563 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
564 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
566 AliModule* det = (AliModule*) detArray->At(iDet);
567 if (!det || !det->IsActive()) continue;
568 if (IsSelected(det->GetName(), detStr)) {
569 det->UpdateInternalGeometry();
571 } // end loop over detectors
574 if (delRunLoader) delete runLoader;
579 //_____________________________________________________________________________
580 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
582 // add a file with background events for merging
584 TObjString* fileNameStr = new TObjString(fileName);
585 fileNameStr->SetUniqueID(nSignalPerBkgrd);
586 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
587 fBkgrdFileNames->Add(fileNameStr);
590 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
592 // add a file with background events for embeddin
593 MergeWith(fileName, nSignalPerBkgrd);
594 fEmbeddingFlag = kTRUE;
597 //_____________________________________________________________________________
598 Bool_t AliSimulation::Run(Int_t nEvents)
600 // run the generation, simulation and digitization
604 AliSysInfo::AddStamp("Start_Run");
606 // Load run number and seed from environmental vars
607 ProcessEnvironmentVars();
608 AliSysInfo::AddStamp("ProcessEnvironmentVars");
610 gRandom->SetSeed(fSeed);
612 if (nEvents > 0) fNEvents = nEvents;
615 // generation and simulation -> hits
616 if (fRunGeneration) {
617 if (!RunSimulation()) if (fStopOnError) return kFALSE;
619 AliSysInfo::AddStamp("RunSimulation");
621 // initialize CDB storage from external environment
622 // (either CDB manager or AliSimulation setters),
623 // if not already done in RunSimulation()
625 AliSysInfo::AddStamp("InitCDB");
627 // Set run number in CDBManager from data
628 // From this point on the run number must be always loaded from data!
629 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
631 // Set CDB lock: from now on it is forbidden to reset the run number
632 // or the default storage or to activate any further storage!
635 // If RunSimulation was not called, load the geometry and misalign it
636 if (!AliGeomManager::GetGeometry()) {
637 // Initialize the geometry manager
638 AliGeomManager::LoadGeometry("geometry.root");
639 AliSysInfo::AddStamp("GetGeometry");
642 // // Check that the consistency of symbolic names for the activated subdetectors
643 // // in the geometry loaded by AliGeomManager
644 // AliRunLoader* runLoader = LoadRun("READ");
645 // if (!runLoader) return kFALSE;
647 // TString detsToBeChecked = "";
648 // TObjArray* detArray = runLoader->GetAliRun()->Detectors();
649 // for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
650 // AliModule* det = (AliModule*) detArray->At(iDet);
651 // if (!det || !det->IsActive()) continue;
652 // detsToBeChecked += det->GetName();
653 // detsToBeChecked += " ";
654 // } // end loop over detectors
655 // if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
656 if(!AliGeomManager::CheckSymNamesLUT("ALL"))
657 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
659 if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
661 if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
663 AliSysInfo::AddStamp("MissalignGeometry");
666 // hits -> summable digits
667 AliSysInfo::AddStamp("Start_sdigitization");
668 if (!fMakeSDigits.IsNull()) {
669 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
672 AliSysInfo::AddStamp("Stop_sdigitization");
674 AliSysInfo::AddStamp("Start_digitization");
675 // summable digits -> digits
676 if (!fMakeDigits.IsNull()) {
677 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
678 if (fStopOnError) return kFALSE;
681 AliSysInfo::AddStamp("Stop_digitization");
686 if (!fMakeDigitsFromHits.IsNull()) {
687 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
688 AliWarning(Form("Merging and direct creation of digits from hits "
689 "was selected for some detectors. "
690 "No merging will be done for the following detectors: %s",
691 fMakeDigitsFromHits.Data()));
693 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
694 if (fStopOnError) return kFALSE;
698 AliSysInfo::AddStamp("Hits2Digits");
702 if (!RunTrigger(fMakeTrigger,fMakeDigits)) {
703 if (fStopOnError) return kFALSE;
706 AliSysInfo::AddStamp("RunTrigger");
709 // digits -> raw data
710 if (!fWriteRawData.IsNull()) {
711 if (!WriteRawData(fWriteRawData, fRawDataFileName,
712 fDeleteIntermediateFiles,fWriteSelRawData)) {
713 if (fStopOnError) return kFALSE;
717 AliSysInfo::AddStamp("WriteRaw");
719 // run HLT simulation on simulated digit data if raw data is not
720 // simulated, otherwise its called as part of WriteRawData
721 if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
723 if (fStopOnError) return kFALSE;
727 AliSysInfo::AddStamp("RunHLT");
731 Bool_t rv = RunQA() ;
737 AliSysInfo::AddStamp("RunQA");
739 // Cleanup of CDB manager: cache and active storages!
740 AliCDBManager::Instance()->ClearCache();
745 //_______________________________________________________________________
746 Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
747 Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
748 Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
751 // Generates lego plots of:
752 // - radiation length map phi vs theta
753 // - radiation length map phi vs eta
754 // - interaction length map
755 // - g/cm2 length map
757 // ntheta bins in theta, eta
758 // themin minimum angle in theta (degrees)
759 // themax maximum angle in theta (degrees)
761 // phimin minimum angle in phi (degrees)
762 // phimax maximum angle in phi (degrees)
763 // rmin minimum radius
764 // rmax maximum radius
767 // The number of events generated = ntheta*nphi
768 // run input parameters in macro setup (default="Config.C")
770 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
773 <img src="picts/AliRunLego1.gif">
778 <img src="picts/AliRunLego2.gif">
783 <img src="picts/AliRunLego3.gif">
788 // run the generation and simulation
792 // initialize CDB storage and run number from external environment
793 // (either CDB manager or AliSimulation setters)
799 AliError("no gAlice object. Restart aliroot and try again.");
802 if (gAlice->Modules()->GetEntries() > 0) {
803 AliError("gAlice was already run. Restart aliroot and try again.");
807 AliInfo(Form("initializing gAlice with config file %s",
808 fConfigFileName.Data()));
811 if (nev == -1) nev = nc1 * nc2;
813 // check if initialisation has been done
814 // If runloader has been initialized, set the number of events per file to nc1 * nc2
817 if (!gener) gener = new AliLegoGenerator();
819 // Configure Generator
821 gener->SetRadiusRange(rmin, rmax);
822 gener->SetZMax(zmax);
823 gener->SetCoor1Range(nc1, c1min, c1max);
824 gener->SetCoor2Range(nc2, c2min, c2max);
828 fLego = new AliLego("lego",gener);
830 //__________________________________________________________________________
834 gROOT->LoadMacro(setup);
835 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
837 if(AliCDBManager::Instance()->GetRun() >= 0) {
838 SetRunNumber(AliCDBManager::Instance()->GetRun());
840 AliWarning("Run number not initialized!!");
843 AliRunLoader::Instance()->CdGAFile();
845 AliPDG::AddParticlesToPdgDataBase();
847 gAlice->GetMCApp()->Init();
849 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
852 //Must be here because some MCs (G4) adds detectors here and not in Config.C
853 gAlice->InitLoaders();
854 AliRunLoader::Instance()->MakeTree("E");
857 // Save stuff at the beginning of the file to avoid file corruption
858 AliRunLoader::Instance()->CdGAFile();
861 //Save current generator
862 AliGenerator *gen=gAlice->GetMCApp()->Generator();
863 gAlice->GetMCApp()->ResetGenerator(gener);
864 //Prepare MC for Lego Run
870 AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
871 gMC->ProcessRun(nev);
873 // End of this run, close files
875 // Restore current generator
876 gAlice->GetMCApp()->ResetGenerator(gen);
877 // Delete Lego Object
883 //_____________________________________________________________________________
884 Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
890 // initialize CDB storage from external environment
891 // (either CDB manager or AliSimulation setters),
892 // if not already done in RunSimulation()
895 // Set run number in CDBManager from data
896 // From this point on the run number must be always loaded from data!
897 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
899 // Set CDB lock: from now on it is forbidden to reset the run number
900 // or the default storage or to activate any further storage!
903 AliRunLoader* runLoader = LoadRun("READ");
904 if (!runLoader) return kFALSE;
905 TString trconfiguration = config;
907 if (trconfiguration.IsNull()) {
908 if (strcmp(gAlice->GetTriggerDescriptor(),"")) {
909 trconfiguration = gAlice->GetTriggerDescriptor();
912 AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
915 runLoader->MakeTree( "GG" );
916 AliCentralTrigger* aCTP = runLoader->GetTrigger();
917 // Load Configuration
918 if (!aCTP->LoadConfiguration( trconfiguration ))
922 if( !aCTP->RunTrigger( runLoader , detectors ) ) {
934 //_____________________________________________________________________________
935 Bool_t AliSimulation::WriteTriggerRawData()
937 // Writes the CTP (trigger) DDL raw data
938 // Details of the format are given in the
939 // trigger TDR - pages 134 and 135.
940 AliCTPRawData writer;
946 //_____________________________________________________________________________
947 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
949 // run the generation and simulation
953 // initialize CDB storage and run number from external environment
954 // (either CDB manager or AliSimulation setters)
955 AliSysInfo::AddStamp("RunSimulation_Begin");
957 AliSysInfo::AddStamp("RunSimulation_InitCDB");
960 AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
963 AliError("no gAlice object. Restart aliroot and try again.");
966 if (gAlice->Modules()->GetEntries() > 0) {
967 AliError("gAlice was already run. Restart aliroot and try again.");
971 AliInfo(Form("initializing gAlice with config file %s",
972 fConfigFileName.Data()));
975 // Initialize ALICE Simulation run
980 gROOT->LoadMacro(fConfigFileName.Data());
981 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
982 AliSysInfo::AddStamp("RunSimulation_Config");
984 if(AliCDBManager::Instance()->GetRun() >= 0) {
985 AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
986 AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
988 AliWarning("Run number not initialized!!");
991 AliRunLoader::Instance()->CdGAFile();
993 AliPDG::AddParticlesToPdgDataBase();
995 gAlice->GetMCApp()->Init();
996 AliSysInfo::AddStamp("RunSimulation_InitMCApp");
998 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
999 AliSysInfo::AddStamp("RunSimulation_GetField");
1001 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1002 gAlice->InitLoaders();
1003 AliRunLoader::Instance()->MakeTree("E");
1004 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1005 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1006 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1008 // Save stuff at the beginning of the file to avoid file corruption
1009 AliRunLoader::Instance()->CdGAFile();
1011 gAlice->SetEventNrInRun(-1); //important - we start Begin event from increasing current number in run
1012 AliSysInfo::AddStamp("RunSimulation_InitLoaders");
1013 //___________________________________________________________________________________________
1015 // Get the trigger descriptor string
1016 // Either from AliSimulation or from
1018 if (fMakeTrigger.IsNull()) {
1019 if (strcmp(gAlice->GetTriggerDescriptor(),""))
1020 fMakeTrigger = gAlice->GetTriggerDescriptor();
1023 gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
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("") ;
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("")
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("")
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("")
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 if (IsSelected("HLT", detStr) && !fRunHLT.IsNull()) {
1312 if (fStopOnError) return kFALSE;
1315 AliSysInfo::AddStamp("WriteRawData_RunHLT");
1317 TString dateFileName(fileName);
1318 if (!dateFileName.IsNull()) {
1319 Bool_t rootOutput = dateFileName.EndsWith(".root");
1320 if (rootOutput) dateFileName += ".date";
1321 TString selDateFileName;
1323 selDateFileName = "selected.";
1324 selDateFileName+= dateFileName;
1326 if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
1327 if (fStopOnError) return kFALSE;
1329 AliSysInfo::AddStamp("ConvertRawFilesToDate");
1330 if (deleteIntermediateFiles) {
1331 AliRunLoader* runLoader = LoadRun("READ");
1332 if (runLoader) for (Int_t iEvent = 0;
1333 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1335 sprintf(command, "rm -r raw%d", iEvent);
1336 gSystem->Exec(command);
1341 if (!ConvertDateToRoot(dateFileName, fileName)) {
1342 if (fStopOnError) return kFALSE;
1344 AliSysInfo::AddStamp("ConvertDateToRoot");
1345 if (deleteIntermediateFiles) {
1346 gSystem->Unlink(dateFileName);
1349 TString selFileName = "selected.";
1350 selFileName += fileName;
1351 if (!ConvertDateToRoot(selDateFileName, selFileName)) {
1352 if (fStopOnError) return kFALSE;
1354 if (deleteIntermediateFiles) {
1355 gSystem->Unlink(selDateFileName);
1364 //_____________________________________________________________________________
1365 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1367 // convert the digits to raw data DDL files
1369 AliCodeTimerAuto("")
1371 AliRunLoader* runLoader = LoadRun("READ");
1372 if (!runLoader) return kFALSE;
1374 // write raw data to DDL files
1375 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1376 AliInfo(Form("processing event %d", iEvent));
1377 runLoader->GetEvent(iEvent);
1378 TString baseDir = gSystem->WorkingDirectory();
1380 sprintf(dirName, "raw%d", iEvent);
1381 gSystem->MakeDirectory(dirName);
1382 if (!gSystem->ChangeDirectory(dirName)) {
1383 AliError(Form("couldn't change to directory %s", dirName));
1384 if (fStopOnError) return kFALSE; else continue;
1387 ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
1390 TString detStr = detectors;
1391 if (IsSelected("HLT", detStr)) {
1392 // Do nothing. "HLT" will be removed from detStr and HLT raw
1393 // data files are generated in RunHLT.
1396 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1397 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1398 AliModule* det = (AliModule*) detArray->At(iDet);
1399 if (!det || !det->IsActive()) continue;
1400 if (IsSelected(det->GetName(), detStr)) {
1401 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1406 if (!WriteTriggerRawData())
1407 if (fStopOnError) return kFALSE;
1409 gSystem->ChangeDirectory(baseDir);
1410 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1411 AliError(Form("the following detectors were not found: %s",
1413 if (fStopOnError) return kFALSE;
1422 //_____________________________________________________________________________
1423 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
1424 const char* selDateFileName)
1426 // convert raw data DDL files to a DATE file with the program "dateStream"
1427 // The second argument is not empty when the user decides to write
1428 // the detectors raw data according to the trigger cluster.
1430 AliCodeTimerAuto("")
1432 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1434 AliError("the program dateStream was not found");
1435 if (fStopOnError) return kFALSE;
1440 AliRunLoader* runLoader = LoadRun("READ");
1441 if (!runLoader) return kFALSE;
1443 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1444 Bool_t selrawdata = kFALSE;
1445 if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
1448 // Note the option -s. It is used in order to avoid
1449 // the generation of SOR/EOR events.
1450 sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d",
1451 dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
1452 FILE* pipe = gSystem->OpenPipe(command, "w");
1455 AliError(Form("Cannot execute command: %s",command));
1459 Int_t selEvents = 0;
1460 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1462 UInt_t detectorPattern = 0;
1463 runLoader->GetEvent(iEvent);
1464 if (!runLoader->LoadTrigger()) {
1465 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1466 detectorPattern = aCTP->GetClusterMask();
1467 // Check if the event was triggered by CTP
1469 if (aCTP->GetClassMask()) selEvents++;
1473 AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
1475 AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
1476 selrawdata = kFALSE;
1480 fprintf(pipe, "GDC DetectorPattern %u\n", detectorPattern);
1484 // loop over detectors and DDLs
1485 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1486 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1488 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1489 Int_t ldcID = Int_t(ldc + 0.0001);
1490 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1492 char rawFileName[256];
1493 sprintf(rawFileName, "raw%d/%s",
1494 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1496 // check existence and size of raw data file
1497 FILE* file = fopen(rawFileName, "rb");
1498 if (!file) continue;
1499 fseek(file, 0, SEEK_END);
1500 unsigned long size = ftell(file);
1502 if (!size) continue;
1504 if (ldcID != prevLDC) {
1505 fprintf(pipe, " LDC Id %d\n", ldcID);
1508 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1513 Int_t result = gSystem->ClosePipe(pipe);
1515 if (!(selrawdata && selEvents > 0)) {
1517 return (result == 0);
1520 AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
1522 sprintf(command, "dateStream -c -s -D -o %s -# %d -C -run %d",
1523 selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
1524 FILE* pipe2 = gSystem->OpenPipe(command, "w");
1526 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1528 // Get the trigger decision and cluster
1529 UInt_t detectorPattern = 0;
1531 runLoader->GetEvent(iEvent);
1532 if (!runLoader->LoadTrigger()) {
1533 AliCentralTrigger *aCTP = runLoader->GetTrigger();
1534 if (aCTP->GetClassMask() == 0) continue;
1535 detectorPattern = aCTP->GetClusterMask();
1536 detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
1537 AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
1540 fprintf(pipe2, "GDC DetectorPattern %u\n", detectorPattern);
1544 // loop over detectors and DDLs
1545 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1546 // Write only raw data from detectors that
1547 // are contained in the trigger cluster(s)
1548 if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
1550 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1552 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1553 Int_t ldcID = Int_t(ldc + 0.0001);
1554 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1556 char rawFileName[256];
1557 sprintf(rawFileName, "raw%d/%s",
1558 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1560 // check existence and size of raw data file
1561 FILE* file = fopen(rawFileName, "rb");
1562 if (!file) continue;
1563 fseek(file, 0, SEEK_END);
1564 unsigned long size = ftell(file);
1566 if (!size) continue;
1568 if (ldcID != prevLDC) {
1569 fprintf(pipe2, " LDC Id %d\n", ldcID);
1572 fprintf(pipe2, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1577 Int_t result2 = gSystem->ClosePipe(pipe2);
1580 return ((result == 0) && (result2 == 0));
1583 //_____________________________________________________________________________
1584 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1585 const char* rootFileName)
1587 // convert a DATE file to a root file with the program "alimdc"
1590 const Int_t kDBSize = 2000000000;
1591 const Int_t kTagDBSize = 1000000000;
1592 const Bool_t kFilter = kFALSE;
1593 const Int_t kCompression = 1;
1595 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1597 AliError("the program alimdc was not found");
1598 if (fStopOnError) return kFALSE;
1603 AliInfo(Form("converting DATE file %s to root file %s",
1604 dateFileName, rootFileName));
1606 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1607 const char* tagDBFS = "/tmp/mdc1/tags";
1609 // User defined file system locations
1610 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1611 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1612 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1613 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1614 if (gSystem->Getenv("ALIMDC_TAGDB"))
1615 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1617 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1618 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1619 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1621 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1622 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1623 gSystem->Exec(Form("mkdir %s",tagDBFS));
1625 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1626 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1627 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1629 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1630 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1631 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1633 return (result == 0);
1637 //_____________________________________________________________________________
1638 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1640 // delete existing run loaders, open a new one and load gAlice
1642 delete AliRunLoader::Instance();
1643 AliRunLoader* runLoader =
1644 AliRunLoader::Open(fGAliceFileName.Data(),
1645 AliConfig::GetDefaultEventFolderName(), mode);
1647 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1650 runLoader->LoadgAlice();
1651 runLoader->LoadHeader();
1652 gAlice = runLoader->GetAliRun();
1654 AliError(Form("no gAlice object found in file %s",
1655 fGAliceFileName.Data()));
1661 //_____________________________________________________________________________
1662 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1664 // get or calculate the number of signal events per background event
1666 if (!fBkgrdFileNames) return 1;
1667 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1668 if (nBkgrdFiles == 0) return 1;
1670 // get the number of signal events
1672 AliRunLoader* runLoader =
1673 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1674 if (!runLoader) return 1;
1676 nEvents = runLoader->GetNumberOfEvents();
1681 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1682 // get the number of background events
1683 const char* fileName = ((TObjString*)
1684 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1685 AliRunLoader* runLoader =
1686 AliRunLoader::Open(fileName, "BKGRD");
1687 if (!runLoader) continue;
1688 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1691 // get or calculate the number of signal per background events
1692 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1693 if (nSignalPerBkgrd <= 0) {
1694 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1695 } else if (result && (result != nSignalPerBkgrd)) {
1696 AliInfo(Form("the number of signal events per background event "
1697 "will be changed from %d to %d for stream %d",
1698 nSignalPerBkgrd, result, iBkgrdFile+1));
1699 nSignalPerBkgrd = result;
1702 if (!result) result = nSignalPerBkgrd;
1703 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1704 AliWarning(Form("not enough background events (%d) for %d signal events "
1705 "using %d signal per background events for stream %d",
1706 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1713 //_____________________________________________________________________________
1714 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1716 // check whether detName is contained in detectors
1717 // if yes, it is removed from detectors
1719 // check if all detectors are selected
1720 if ((detectors.CompareTo("ALL") == 0) ||
1721 detectors.BeginsWith("ALL ") ||
1722 detectors.EndsWith(" ALL") ||
1723 detectors.Contains(" ALL ")) {
1728 // search for the given detector
1729 Bool_t result = kFALSE;
1730 if ((detectors.CompareTo(detName) == 0) ||
1731 detectors.BeginsWith(detName+" ") ||
1732 detectors.EndsWith(" "+detName) ||
1733 detectors.Contains(" "+detName+" ")) {
1734 detectors.ReplaceAll(detName, "");
1738 // clean up the detectors string
1739 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1740 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1741 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1746 //_____________________________________________________________________________
1747 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName)
1750 // Steering routine to convert raw data in directory rawDirectory/ to fake SDigits.
1751 // These can be used for embedding of MC tracks into RAW data using the standard
1752 // merging procedure.
1754 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1757 AliError("no gAlice object. Restart aliroot and try again.");
1760 if (gAlice->Modules()->GetEntries() > 0) {
1761 AliError("gAlice was already run. Restart aliroot and try again.");
1765 AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1769 gROOT->LoadMacro(fConfigFileName.Data());
1770 gInterpreter->ProcessLine(gAlice->GetConfigFunction());
1772 if(AliCDBManager::Instance()->GetRun() >= 0) {
1773 SetRunNumber(AliCDBManager::Instance()->GetRun());
1775 AliWarning("Run number not initialized!!");
1778 AliRunLoader::Instance()->CdGAFile();
1780 AliPDG::AddParticlesToPdgDataBase();
1782 gAlice->GetMCApp()->Init();
1784 gMC->SetMagField(TGeoGlobalMagField::Instance()->GetField());
1786 //Must be here because some MCs (G4) adds detectors here and not in Config.C
1787 gAlice->InitLoaders();
1788 AliRunLoader::Instance()->MakeTree("E");
1789 AliRunLoader::Instance()->LoadKinematics("RECREATE");
1790 AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
1791 AliRunLoader::Instance()->LoadHits("all","RECREATE");
1793 // Save stuff at the beginning of the file to avoid file corruption
1794 AliRunLoader::Instance()->CdGAFile();
1799 //AliCDBManager* man = AliCDBManager::Instance();
1800 //man->SetRun(0); // Should this come from rawdata header ?
1804 // Get the runloader
1805 AliRunLoader* runLoader = AliRunLoader::Instance();
1807 // Open esd file if available
1808 TFile* esdFile = TFile::Open(esdFileName);
1809 Bool_t esdOK = (esdFile != 0);
1810 AliESD* esd = new AliESD;
1813 treeESD = (TTree*) esdFile->Get("esdTree");
1815 AliWarning("No ESD tree found");
1818 treeESD->SetBranchAddress("ESD", &esd);
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::RunHLT()
1941 // Run the HLT simulation
1942 // HLT simulation is implemented in HLT/sim/AliHLTSimulation
1943 // Disabled if fRunHLT is empty, default vaule is "default".
1944 // AliSimulation::SetRunHLT can be used to set the options for HLT simulation
1945 // The default simulation depends on the HLT component libraries and their
1946 // corresponding agents which define components and chains to run. See
1947 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/
1948 // http://web.ift.uib.no/~kjeks/doc/alice-hlt-current/classAliHLTModuleAgent.html
1950 // The libraries to be loaded can be specified as an option.
1952 // AliSimulation sim;
1953 // sim.SetRunHLT("libAliHLTSample.so");
1955 // will only load <tt>libAliHLTSample.so</tt>
1957 // Other available options:
1958 // \li loglevel=<i>level</i> <br>
1959 // logging level for this processing
1961 // disable redirection of log messages to AliLog class
1962 // \li config=<i>macro</i>
1963 // configuration macro
1964 // \li chains=<i>configuration</i>
1965 // comma separated list of configurations to be run during simulation
1966 // \li rawfile=<i>file</i>
1967 // source for the RawReader to be created, the default is <i>./</i> if
1968 // raw data is simulated
1971 AliRunLoader* pRunLoader = LoadRun("READ");
1972 if (!pRunLoader) return kFALSE;
1974 // initialize CDB storage, run number, set CDB lock
1976 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
1979 // load the library dynamically
1980 gSystem->Load(ALIHLTSIMULATION_LIBRARY);
1982 // check for the library version
1983 AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
1985 AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
1988 if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
1989 AliError(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
1993 // print compile info
1994 typedef void (*CompileInfo)( const char*& date, const char*& time);
1995 CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
1997 const char* date="";
1998 const char* time="";
1999 (*fctInfo)(date, time);
2000 if (!date) date="unknown";
2001 if (!time) time="unknown";
2002 AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
2004 AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
2007 // create instance of the HLT simulation
2008 AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
2009 AliHLTSimulation* pHLT=NULL;
2010 if (fctCreate==NULL || (pHLT=(fctCreate()))==NULL) {
2011 AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
2015 // init the HLT simulation
2017 if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
2018 TString detStr = fWriteRawData;
2019 if (!IsSelected("HLT", detStr)) {
2020 options+=" writerawfiles=";
2022 options+=" writerawfiles=HLT";
2025 if (!detStr.IsNull() && !options.Contains("rawfile=")) {
2026 // as a matter of fact, HLT will run reconstruction and needs the RawReader
2027 // in order to get detector data. By default, RawReaderFile is used to read
2028 // the already simulated ddl files. Date and Root files from the raw data
2029 // are generated after the HLT simulation.
2030 options+=" rawfile=./";
2033 AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
2034 if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
2035 AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
2037 // run the HLT simulation
2038 AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
2039 if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
2040 AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
2044 // delete the instance
2045 AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
2046 if (fctDelete==NULL || fctDelete(pHLT)<0) {
2047 AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
2051 return iResult>=0?kTRUE:kFALSE;
2054 //_____________________________________________________________________________
2055 Bool_t AliSimulation::RunQA()
2057 // run the QA on summable hits, digits or digits
2059 if(!gAlice) return kFALSE;
2060 AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
2062 TString detectorsw("") ;
2064 AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
2065 detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ;
2066 if ( detectorsw.IsNull() )
2071 //_____________________________________________________________________________
2072 Bool_t AliSimulation::SetRunQA(TString detAndAction)
2074 // Allows to run QA for a selected set of detectors
2075 // and a selected set of tasks among HITS, SDIGITS and DIGITS
2076 // all selected detectors run the same selected tasks
2078 if (!detAndAction.Contains(":")) {
2079 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2083 Int_t colon = detAndAction.Index(":") ;
2084 fQADetectors = detAndAction(0, colon) ;
2085 if (fQADetectors.Contains("ALL") )
2086 fQADetectors = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
2087 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2088 if (fQATasks.Contains("ALL") ) {
2089 fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
2091 fQATasks.ToUpper() ;
2093 if ( fQATasks.Contains("HIT") )
2094 tempo = Form("%d ", AliQAv1::kHITS) ;
2095 if ( fQATasks.Contains("SDIGIT") )
2096 tempo += Form("%d ", AliQAv1::kSDIGITS) ;
2097 if ( fQATasks.Contains("DIGIT") )
2098 tempo += Form("%d ", AliQAv1::kDIGITS) ;
2100 if (fQATasks.IsNull()) {
2101 AliInfo("No QA requested\n") ;
2106 TString tempo(fQATasks) ;
2107 tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS)) ;
2108 tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;
2109 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;
2110 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2112 AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ;
2113 AliQAManager::QAManager()->SetTasks(fQATasks) ;
2114 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++)
2115 AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
2120 //_____________________________________________________________________________
2121 void AliSimulation::ProcessEnvironmentVars()
2123 // Extract run number and random generator seed from env variables
2125 AliInfo("Processing environment variables");
2127 // Random Number seed
2129 // first check that seed is not already set
2131 if (gSystem->Getenv("CONFIG_SEED")) {
2132 fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
2135 if (gSystem->Getenv("CONFIG_SEED")) {
2136 AliInfo(Form("Seed for random number generation already set (%d)"
2137 ": CONFIG_SEED variable ignored!", fSeed));
2141 AliInfo(Form("Seed for random number generation = %d ", fSeed));
2145 // first check that run number is not already set
2147 if (gSystem->Getenv("DC_RUN")) {
2148 fRun = atoi(gSystem->Getenv("DC_RUN"));
2151 if (gSystem->Getenv("DC_RUN")) {
2152 AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
2156 AliInfo(Form("Run number = %d", fRun));
2159 //---------------------------------------------------------------------
2160 void AliSimulation::WriteGRPEntry()
2162 // Get the necessary information from galice (generator, trigger etc) and
2163 // write a GRP entry corresponding to the settings in the Config.C used
2164 // note that Hall probes and Cavern and Surface Atmos pressures are not simulated.
2167 AliInfo("Writing global run parameters entry into the OCDB");
2169 AliGRPObject* grpObj = new AliGRPObject();
2171 grpObj->SetRunType("PHYSICS");
2172 grpObj->SetTimeStart(0);
2173 grpObj->SetTimeEnd(9999);
2175 const AliGenerator *gen = gAlice->GetMCApp()->Generator();
2177 grpObj->SetBeamEnergy(gen->GetEnergyCMS());
2180 gen->GetProjectile(projectile,a,z);
2182 gen->GetTarget(target,a,z);
2183 TString beamType = projectile + "-" + target;
2184 beamType.ReplaceAll(" ","");
2185 if (!beamType.CompareTo("-")) {
2186 grpObj->SetBeamType("UNKNOWN");
2189 grpObj->SetBeamType(beamType);
2190 // Heavy ion run, the event specie is set to kHighMult
2191 fEventSpecie = AliRecoParam::kHighMult;
2192 if ((strcmp(beamType,"p-p") == 0) ||
2193 (strcmp(beamType,"p-") == 0) ||
2194 (strcmp(beamType,"-p") == 0) ||
2195 (strcmp(beamType,"P-P") == 0) ||
2196 (strcmp(beamType,"P-") == 0) ||
2197 (strcmp(beamType,"-P") == 0)) {
2198 // Proton run, the event specie is set to kLowMult
2199 fEventSpecie = AliRecoParam::kLowMult;
2203 AliWarning("Unknown beam type and energy! Setting energy to 0");
2204 grpObj->SetBeamEnergy(0);
2205 grpObj->SetBeamType("UNKNOWN");
2208 UInt_t detectorPattern = 0;
2210 TObjArray *detArray = gAlice->Detectors();
2211 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
2212 if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
2213 detectorPattern |= (1 << iDet);
2218 if (!fMakeTrigger.IsNull() || strcmp(gAlice->GetTriggerDescriptor(),""))
2219 detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
2222 if (!fRunHLT.IsNull())
2223 detectorPattern |= (1 << AliDAQ::kHLTId);
2225 grpObj->SetNumberOfDetectors((Char_t)nDets);
2226 grpObj->SetDetectorMask((Int_t)detectorPattern);
2227 grpObj->SetLHCPeriod("LHC08c");
2228 grpObj->SetLHCState("STABLE_BEAMS");
2229 grpObj->SetLHCLuminosity(0,(AliGRPObject::Stats)0);
2230 grpObj->SetBeamIntensity(0,(AliGRPObject::Stats)0);
2232 AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
2233 Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
2234 Float_t factor = field ? field->Factor() : 0;
2235 Float_t l3current = TMath::Abs(factor)*solenoidField*30000./5.;
2236 grpObj->SetL3Current(l3current,(AliGRPObject::Stats)0);
2239 grpObj->SetL3Polarity(0);
2240 grpObj->SetDipolePolarity(0);
2243 grpObj->SetL3Polarity(1);
2244 grpObj->SetDipolePolarity(1);
2247 if (TMath::Abs(factor) != 0)
2248 grpObj->SetDipoleCurrent(6000,(AliGRPObject::Stats)0);
2250 grpObj->SetDipoleCurrent(0,(AliGRPObject::Stats)0);
2252 grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
2254 //grpMap->Add(new TObjString("fCavernPressure"),new TObjString("0")); ---> not inserted in simulation with the new object, since it is now an AliDCSSensor
2256 // Now store the entry in OCDB
2257 AliCDBManager* man = AliCDBManager::Instance();
2259 AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
2260 AliCDBMetaData *metadata= new AliCDBMetaData();
2262 metadata->SetResponsible("alice-off@cern.ch");
2263 metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
2265 man->Put(grpObj,id,metadata);