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 ///////////////////////////////////////////////////////////////////////////////
108 #include <TGeoManager.h>
109 #include <TObjString.h>
110 #include <TStopwatch.h>
114 #include "AliCDBStorage.h"
115 #include "AliCDBEntry.h"
116 #include "AliCDBManager.h"
117 #include "AliAlignObj.h"
118 #include "AliCentralTrigger.h"
119 #include "AliDAQConfig.h"
120 #include "AliDigitizer.h"
121 #include "AliGenerator.h"
123 #include "AliModule.h"
125 #include "AliRunDigitizer.h"
126 #include "AliRunLoader.h"
127 #include "AliSimulation.h"
128 #include "AliVertexGenFile.h"
129 #include "AliCentralTrigger.h"
131 ClassImp(AliSimulation)
134 //_____________________________________________________________________________
135 AliSimulation::AliSimulation(const char* configFileName, const char* cdbUri,
136 const char* name, const char* title) :
139 fRunGeneration(kTRUE),
140 fRunSimulation(kTRUE),
141 fLoadAlignFromCDB(kTRUE),
142 fLoadAlignData("ALL"),
146 fMakeDigitsFromHits(""),
148 fRawDataFileName(""),
149 fDeleteIntermediateFiles(kFALSE),
150 fStopOnError(kFALSE),
153 fConfigFileName(configFileName),
154 fGAliceFileName("galice.root"),
156 fBkgrdFileNames(NULL),
157 fAlignObjArray(NULL),
158 fUseBkgrdVertex(kTRUE),
159 fRegionOfInterest(kFALSE),
162 // create simulation object with default parameters
164 SetGAliceFile("galice.root");
167 //_____________________________________________________________________________
168 AliSimulation::AliSimulation(const AliSimulation& sim) :
171 fRunGeneration(sim.fRunGeneration),
172 fRunSimulation(sim.fRunSimulation),
173 fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
174 fLoadAlignData(sim.fLoadAlignData),
175 fMakeSDigits(sim.fMakeSDigits),
176 fMakeDigits(sim.fMakeDigits),
177 fMakeTrigger(sim.fMakeTrigger),
178 fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
179 fWriteRawData(sim.fWriteRawData),
180 fRawDataFileName(""),
181 fDeleteIntermediateFiles(kFALSE),
182 fStopOnError(sim.fStopOnError),
184 fNEvents(sim.fNEvents),
185 fConfigFileName(sim.fConfigFileName),
186 fGAliceFileName(sim.fGAliceFileName),
188 fBkgrdFileNames(NULL),
189 fAlignObjArray(NULL),
190 fUseBkgrdVertex(sim.fUseBkgrdVertex),
191 fRegionOfInterest(sim.fRegionOfInterest),
196 for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
197 if (!sim.fEventsPerFile[i]) continue;
198 fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
201 fBkgrdFileNames = new TObjArray;
202 for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
203 if (!sim.fBkgrdFileNames->At(i)) continue;
204 fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
208 //_____________________________________________________________________________
209 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
211 // assignment operator
213 this->~AliSimulation();
214 new(this) AliSimulation(sim);
218 //_____________________________________________________________________________
219 AliSimulation::~AliSimulation()
223 fEventsPerFile.Delete();
224 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
225 // delete fAlignObjArray; fAlignObjArray=0;
227 if (fBkgrdFileNames) {
228 fBkgrdFileNames->Delete();
229 delete fBkgrdFileNames;
234 //_____________________________________________________________________________
235 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
237 // set the number of events for one run
242 //_____________________________________________________________________________
243 void AliSimulation::InitCDBStorage(const char* uri)
245 // activate a default CDB storage
246 // First check if we have any CDB storage set, because it is used
247 // to retrieve the calibration and alignment constants
249 AliCDBManager* man = AliCDBManager::Instance();
250 if (!man->IsDefaultStorageSet())
252 AliWarningClass("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
253 AliWarningClass("Default CDB storage not yet set");
254 AliWarningClass(Form("Using default storage declared in AliSimulation: %s",uri));
255 AliWarningClass("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
256 SetDefaultStorage(uri);
261 //_____________________________________________________________________________
262 void AliSimulation::SetDefaultStorage(const char* uri) {
263 // activate a default CDB storage
265 AliCDBManager::Instance()->SetDefaultStorage(uri);
269 //_____________________________________________________________________________
270 void AliSimulation::SetSpecificStorage(const char* detName, const char* uri) {
271 // activate a detector-specific CDB storage
273 AliCDBManager::Instance()->SetSpecificStorage(detName, uri);
277 //_____________________________________________________________________________
278 void AliSimulation::SetConfigFile(const char* fileName)
280 // set the name of the config file
282 fConfigFileName = fileName;
285 //_____________________________________________________________________________
286 void AliSimulation::SetGAliceFile(const char* fileName)
288 // set the name of the galice file
289 // the path is converted to an absolute one if it is relative
291 fGAliceFileName = fileName;
292 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
293 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
295 fGAliceFileName = absFileName;
296 delete[] absFileName;
299 AliDebug(2, Form("galice file name set to %s", fileName));
302 //_____________________________________________________________________________
303 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
306 // set the number of events per file for the given detector and data type
307 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
309 TNamed* obj = new TNamed(detector, type);
310 obj->SetUniqueID(nEvents);
311 fEventsPerFile.Add(obj);
314 //_____________________________________________________________________________
315 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* fileName, const char* clArrayName)
317 // read collection of alignment objects (AliAlignObj derived) saved
318 // in the TClonesArray ClArrayName in the file fileName and apply
319 // them to the TGeo geometry passed as argument
322 TFile* inFile = TFile::Open(fileName,"READ");
323 if (!inFile || !inFile->IsOpen()) {
324 AliErrorClass(Form("Could not open file %s !",fileName));
328 TClonesArray* alObjArray = ((TClonesArray*) inFile->Get(clArrayName));
331 AliErrorClass(Form("Could not get array (%s) from file (%s) !",clArrayName,fileName));
335 return gAlice->ApplyAlignObjsToGeom(alObjArray);
339 //_____________________________________________________________________________
340 Bool_t AliSimulation::ApplyAlignObjsToGeom(AliCDBParam* param, AliCDBId& Id)
342 // read collection of alignment objects (AliAlignObj derived) saved
343 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
344 // param (to get the AliCDBStorage) and Id; apply the alignment objects
345 // to the TGeo geometry passed as argument
348 AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage(param);
349 AliCDBEntry* entry = storage->Get(Id);
350 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
352 return gAlice->ApplyAlignObjsToGeom(AlObjArray);
356 //_____________________________________________________________________________
357 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* uri, const char* path, Int_t runnum, Int_t version, Int_t sversion)
359 // read collection of alignment objects (AliAlignObj derived) saved
360 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
361 // param (to get the AliCDBStorage) and Id; apply the alignment objects
362 // to the TGeo geometry passed as argument
365 AliCDBParam* param = AliCDBManager::Instance()->CreateParameter(uri);
366 AliCDBId id(path, runnum, runnum, version, sversion);
368 return ApplyAlignObjsToGeom(param, id);
372 //_____________________________________________________________________________
373 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* detName, Int_t runnum, Int_t version, Int_t sversion)
375 // read collection of alignment objects (AliAlignObj derived) saved
376 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
377 // param (to get the AliCDBStorage) and Id; apply the alignment objects
378 // to the TGeo geometry passed as argument
381 InitCDBStorage("local://$ALICE_ROOT");
382 AliCDBPath path(detName,"Align","Data");
383 AliCDBEntry* entry = AliCDBManager::Instance()->Get(path.GetPath(),runnum,version,sversion);
385 if(!entry) return kFALSE;
386 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
388 return gAlice->ApplyAlignObjsToGeom(AlObjArray);
392 //_____________________________________________________________________________
393 void AliSimulation::SetAlignObjArray(const char* detectors)
395 // Fills array of detectors' alignable objects from CDB
396 // detectors can be "ALL" or "ITS TPC ..."
398 AliRunLoader* runLoader = LoadRun();
399 if (!runLoader) return;
401 InitCDBStorage(fCDBUri);
403 if(!fAlignObjArray) fAlignObjArray = new TObjArray();
404 fAlignObjArray->SetOwner(0); // AliCDBEntry is owner of the align objects!
405 fAlignObjArray->Clear();
407 TString detStr = detectors;
408 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
409 TString dataNotLoaded="";
410 TString dataLoaded="";
412 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
413 AliModule* det = (AliModule*) detArray->At(iDet);
414 if (!det || !det->IsActive()) continue;
415 if (IsSelected(det->GetName(), detStr)) {
417 if(!SetAlignObjArraySingleDet(det->GetName())){
418 dataNotLoaded += det->GetName();
419 dataNotLoaded += " ";
421 dataLoaded += det->GetName();
425 } // end loop over all detectors
427 if ((detStr.CompareTo("ALL") == 0)) detStr = "";
428 dataNotLoaded += detStr;
429 AliInfo(Form("Alignment data loaded for: %s",
431 AliInfo(Form("Didn't/couldn't load alignment data for: %s",
432 dataNotLoaded.Data()));
434 AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
439 //_____________________________________________________________________________
440 Bool_t AliSimulation::SetAlignObjArraySingleDet(const char* detName)
442 // Fills array of single detector's alignable objects from CDB
444 AliDebug(2, Form("Loading alignment data for detector: %s",detName));
448 AliCDBPath path(detName,"Align","Data");
450 entry=AliCDBManager::Instance()->Get(path.GetPath());
452 AliDebug(2,Form("Couldn't load alignment data for detector %s",detName));
456 TClonesArray *alignArray = (TClonesArray*) entry->GetObject();
457 alignArray->SetOwner(0);
458 AliDebug(2,Form("Found %d alignment objects for %s",
459 alignArray->GetEntries(),detName));
461 AliAlignObj *alignObj=0;
462 TIter iter(alignArray);
464 // loop over align objects in detector
465 while( ( alignObj=(AliAlignObj *) iter.Next() ) ){
466 fAlignObjArray->Add(alignObj);
468 // delete entry --- Don't delete, it is cached!
470 AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
475 //_____________________________________________________________________________
476 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
478 // add a file with background events for merging
480 TObjString* fileNameStr = new TObjString(fileName);
481 fileNameStr->SetUniqueID(nSignalPerBkgrd);
482 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
483 fBkgrdFileNames->Add(fileNameStr);
487 //_____________________________________________________________________________
488 Bool_t AliSimulation::Run(Int_t nEvents)
490 // run the generation, simulation and digitization
492 InitCDBStorage(fCDBUri);
494 if (nEvents > 0) fNEvents = nEvents;
496 // generation and simulation -> hits
497 if (fRunGeneration) {
498 if (!RunSimulation()) if (fStopOnError) return kFALSE;
501 // hits -> summable digits
502 if (!fMakeSDigits.IsNull()) {
503 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
506 // summable digits -> digits
507 if (!fMakeDigits.IsNull()) {
508 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
509 if (fStopOnError) return kFALSE;
514 if (!fMakeDigitsFromHits.IsNull()) {
515 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
516 AliWarning(Form("Merging and direct creation of digits from hits "
517 "was selected for some detectors. "
518 "No merging will be done for the following detectors: %s",
519 fMakeDigitsFromHits.Data()));
521 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
522 if (fStopOnError) return kFALSE;
527 if (!fMakeTrigger.IsNull()) {
528 if (!RunTrigger(fMakeTrigger)) {
529 if (fStopOnError) return kFALSE;
533 // digits -> raw data
534 if (!fWriteRawData.IsNull()) {
535 if (!WriteRawData(fWriteRawData, fRawDataFileName,
536 fDeleteIntermediateFiles)) {
537 if (fStopOnError) return kFALSE;
544 //_____________________________________________________________________________
545 Bool_t AliSimulation::RunTrigger(const char* descriptors)
549 TStopwatch stopwatch;
552 AliRunLoader* runLoader = LoadRun("READ");
553 if (!runLoader) return kFALSE;
554 TString des = descriptors;
556 AliCentralTrigger* aCTP = new AliCentralTrigger( des );
559 if( !aCTP->RunTrigger( runLoader ) ) {
567 // Process each event
568 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
569 AliInfo(Form("processing event %d", iEvent));
570 runLoader->GetEvent(iEvent);
572 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
573 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
574 AliModule* det = (AliModule*) detArray->At(iDet);
575 if (!det || !det->IsActive()) continue;
576 if (IsSelected(det->GetName(), detStr)) {
577 AliInfo(Form("triggering from digits for %s", det->GetName()));
579 // AliLoader* loader = fLoader[iDet];
580 // loader->LoadDigits("read");
581 // TTree* digitsTree = loader->TreeD();
582 // det->Trigger( digitsTree );
584 AliTriggerDetector* tdet = det->CreateTriggerDetector();
585 TObjArray* detInp = dtrg->GetTriggerInputs();
586 for( Int_t i=0; i<detInp->GetEntriesFast(); i++ )
587 fInputs.AddLast( detInp->At(i) );
589 AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
590 det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
594 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
595 AliError(Form("the following detectors were not found: %s",
603 // Check trigger conditions
604 centralTP->TriggerConditions();
606 // Write trigger ????
611 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
612 stopwatch.RealTime(),stopwatch.CpuTime()));
622 //_____________________________________________________________________________
623 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
625 // run the generation and simulation
627 TStopwatch stopwatch;
631 AliError("no gAlice object. Restart aliroot and try again.");
634 if (gAlice->Modules()->GetEntries() > 0) {
635 AliError("gAlice was already run. Restart aliroot and try again.");
639 AliInfo(Form("initializing gAlice with config file %s",
640 fConfigFileName.Data()));
641 StdoutToAliInfo(StderrToAliError(
642 gAlice->Init(fConfigFileName.Data());
645 // Set run number in CDBManager (here????)
646 AliCDBManager::Instance()->SetRun(gAlice->GetRunNumber());
647 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
649 AliRunLoader* runLoader = gAlice->GetRunLoader();
651 AliError(Form("gAlice has no run loader object. "
652 "Check your config file: %s", fConfigFileName.Data()));
655 SetGAliceFile(runLoader->GetFileName());
657 // Export ideal geometry
658 if (gGeoManager) gGeoManager->Export("geometry.root");
660 // Load alignment data from CDB and fill fAlignObjArray
661 if(fLoadAlignFromCDB){
662 if(!fAlignObjArray) fAlignObjArray = new TObjArray();
664 //fAlignObjArray->RemoveAll();
665 fAlignObjArray->Clear();
666 fAlignObjArray->SetOwner(0);
668 TString detStr = fLoadAlignData;
669 TString dataNotLoaded="";
670 TString dataLoaded="";
672 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
673 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
674 AliModule* det = (AliModule*) detArray->At(iDet);
675 if (!det || !det->IsActive()) continue;
676 if (IsSelected(det->GetName(), detStr)) {
677 if(!SetAlignObjArraySingleDet(det->GetName())){
678 dataNotLoaded += det->GetName();
679 dataNotLoaded += " ";
681 dataLoaded += det->GetName();
685 } // end loop over detectors
687 if ((detStr.CompareTo("ALL") == 0)) detStr = "";
688 dataNotLoaded += detStr;
689 AliInfo(Form("Alignment data loaded for: %s",
691 AliInfo(Form("Didn't/couldn't load alignment data for: %s",
692 dataNotLoaded.Data()));
693 } // fLoadAlignFromCDB flag
695 // Check if the array with alignment objects was
696 // provided by the user. If yes, apply the objects
697 // to the present TGeo geometry
698 if (fAlignObjArray) {
699 if (gGeoManager && gGeoManager->IsClosed()) {
700 if (gAlice->ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
701 AliError("The application of misalignment failed! Restart aliroot and try again. ");
706 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
711 // Export (mis)aligned geometry
712 if (gGeoManager) gGeoManager->Export("misaligned_geometry.root");
714 // AliRunLoader* runLoader = gAlice->GetRunLoader();
716 // AliError(Form("gAlice has no run loader object. "
717 // "Check your config file: %s", fConfigFileName.Data()));
720 // SetGAliceFile(runLoader->GetFileName());
722 if (!gAlice->Generator()) {
723 AliError(Form("gAlice has no generator object. "
724 "Check your config file: %s", fConfigFileName.Data()));
727 if (nEvents <= 0) nEvents = fNEvents;
729 // get vertex from background file in case of merging
730 if (fUseBkgrdVertex &&
731 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
732 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
733 const char* fileName = ((TObjString*)
734 (fBkgrdFileNames->At(0)))->GetName();
735 AliInfo(Form("The vertex will be taken from the background "
736 "file %s with nSignalPerBackground = %d",
737 fileName, signalPerBkgrd));
738 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
739 gAlice->Generator()->SetVertexGenerator(vtxGen);
742 if (!fRunSimulation) {
743 gAlice->Generator()->SetTrackingFlag(0);
746 // set the number of events per file for given detectors and data types
747 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
748 if (!fEventsPerFile[i]) continue;
749 const char* detName = fEventsPerFile[i]->GetName();
750 const char* typeName = fEventsPerFile[i]->GetTitle();
751 TString loaderName(detName);
752 loaderName += "Loader";
753 AliLoader* loader = runLoader->GetLoader(loaderName);
755 AliError(Form("RunSimulation", "no loader for %s found\n"
756 "Number of events per file not set for %s %s",
757 detName, typeName, detName));
760 AliDataLoader* dataLoader =
761 loader->GetDataLoader(typeName);
763 AliError(Form("no data loader for %s found\n"
764 "Number of events per file not set for %s %s",
765 typeName, detName, typeName));
768 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
769 AliDebug(1, Form("number of events per file set to %d for %s %s",
770 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
773 AliInfo("running gAlice");
774 StdoutToAliInfo(StderrToAliError(
775 gAlice->Run(nEvents);
780 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
781 stopwatch.RealTime(),stopwatch.CpuTime()));
786 //_____________________________________________________________________________
787 Bool_t AliSimulation::RunSDigitization(const char* detectors)
789 // run the digitization and produce summable digits
791 TStopwatch stopwatch;
794 AliRunLoader* runLoader = LoadRun();
795 if (!runLoader) return kFALSE;
797 TString detStr = detectors;
798 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
799 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
800 AliModule* det = (AliModule*) detArray->At(iDet);
801 if (!det || !det->IsActive()) continue;
802 if (IsSelected(det->GetName(), detStr)) {
803 AliInfo(Form("creating summable digits for %s", det->GetName()));
804 TStopwatch stopwatchDet;
805 stopwatchDet.Start();
807 AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
808 det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
812 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
813 AliError(Form("the following detectors were not found: %s",
815 if (fStopOnError) return kFALSE;
820 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
821 stopwatch.RealTime(),stopwatch.CpuTime()));
827 //_____________________________________________________________________________
828 Bool_t AliSimulation::RunDigitization(const char* detectors,
829 const char* excludeDetectors)
831 // run the digitization and produce digits from sdigits
833 TStopwatch stopwatch;
836 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
837 if (gAlice) delete gAlice;
841 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
842 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
843 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
844 manager->SetInputStream(0, fGAliceFileName.Data());
845 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
846 const char* fileName = ((TObjString*)
847 (fBkgrdFileNames->At(iStream-1)))->GetName();
848 manager->SetInputStream(iStream, fileName);
851 TString detStr = detectors;
852 TString detExcl = excludeDetectors;
853 manager->GetInputStream(0)->ImportgAlice();
854 AliRunLoader* runLoader =
855 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
856 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
857 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
858 AliModule* det = (AliModule*) detArray->At(iDet);
859 if (!det || !det->IsActive()) continue;
860 if (IsSelected(det->GetName(), detStr) &&
861 !IsSelected(det->GetName(), detExcl)) {
862 AliDigitizer* digitizer = det->CreateDigitizer(manager);
864 AliError(Form("no digitizer for %s", det->GetName()));
865 if (fStopOnError) return kFALSE;
867 digitizer->SetRegionOfInterest(fRegionOfInterest);
872 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
873 AliError(Form("the following detectors were not found: %s",
875 if (fStopOnError) return kFALSE;
878 if (!manager->GetListOfTasks()->IsEmpty()) {
879 AliInfo("executing digitization");
885 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
886 stopwatch.RealTime(),stopwatch.CpuTime()));
891 //_____________________________________________________________________________
892 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
894 // run the digitization and produce digits from hits
896 TStopwatch stopwatch;
899 AliRunLoader* runLoader = LoadRun("READ");
900 if (!runLoader) return kFALSE;
902 TString detStr = detectors;
903 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
904 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
905 AliModule* det = (AliModule*) detArray->At(iDet);
906 if (!det || !det->IsActive()) continue;
907 if (IsSelected(det->GetName(), detStr)) {
908 AliInfo(Form("creating digits from hits for %s", det->GetName()));
913 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
914 AliError(Form("the following detectors were not found: %s",
916 if (fStopOnError) return kFALSE;
920 //PH Temporary fix to avoid interference with the PHOS loder/getter
921 //PH The problem has to be solved in more general way 09/06/05
923 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
924 stopwatch.RealTime(),stopwatch.CpuTime()));
929 //_____________________________________________________________________________
930 Bool_t AliSimulation::WriteRawData(const char* detectors,
931 const char* fileName,
932 Bool_t deleteIntermediateFiles)
934 // convert the digits to raw data
935 // First DDL raw data files for the given detectors are created.
936 // If a file name is given, the DDL files are then converted to a DATE file.
937 // If deleteIntermediateFiles is true, the DDL raw files are deleted
939 // If the file name has the extension ".root", the DATE file is converted
941 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
943 TStopwatch stopwatch;
946 if (!WriteRawFiles(detectors)) {
947 if (fStopOnError) return kFALSE;
950 TString dateFileName(fileName);
951 if (!dateFileName.IsNull()) {
952 Bool_t rootOutput = dateFileName.EndsWith(".root");
953 if (rootOutput) dateFileName += ".date";
954 if (!ConvertRawFilesToDate(dateFileName)) {
955 if (fStopOnError) return kFALSE;
957 if (deleteIntermediateFiles) {
958 AliRunLoader* runLoader = LoadRun("READ");
959 if (runLoader) for (Int_t iEvent = 0;
960 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
962 sprintf(command, "rm -r raw%d", iEvent);
963 gSystem->Exec(command);
968 if (!ConvertDateToRoot(dateFileName, fileName)) {
969 if (fStopOnError) return kFALSE;
971 if (deleteIntermediateFiles) {
972 gSystem->Unlink(dateFileName);
977 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
978 stopwatch.RealTime(),stopwatch.CpuTime()));
983 //_____________________________________________________________________________
984 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
986 // convert the digits to raw data DDL files
988 AliRunLoader* runLoader = LoadRun("READ");
989 if (!runLoader) return kFALSE;
991 // write raw data to DDL files
992 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
993 AliInfo(Form("processing event %d", iEvent));
994 runLoader->GetEvent(iEvent);
995 TString baseDir = gSystem->WorkingDirectory();
997 sprintf(dirName, "raw%d", iEvent);
998 gSystem->MakeDirectory(dirName);
999 if (!gSystem->ChangeDirectory(dirName)) {
1000 AliError(Form("couldn't change to directory %s", dirName));
1001 if (fStopOnError) return kFALSE; else continue;
1004 TString detStr = detectors;
1005 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1006 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1007 AliModule* det = (AliModule*) detArray->At(iDet);
1008 if (!det || !det->IsActive()) continue;
1009 if (IsSelected(det->GetName(), detStr)) {
1010 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1015 gSystem->ChangeDirectory(baseDir);
1016 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1017 AliError(Form("the following detectors were not found: %s",
1019 if (fStopOnError) return kFALSE;
1027 //_____________________________________________________________________________
1028 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1030 // convert raw data DDL files to a DATE file with the program "dateStream"
1032 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1034 AliError("the program dateStream was not found");
1035 if (fStopOnError) return kFALSE;
1040 AliRunLoader* runLoader = LoadRun("READ");
1041 if (!runLoader) return kFALSE;
1043 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1045 sprintf(command, "dateStream -o %s -# %d -C",
1046 dateFileName, runLoader->GetNumberOfEvents());
1047 FILE* pipe = gSystem->OpenPipe(command, "w");
1049 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1050 fprintf(pipe, "GDC\n");
1054 // loop over detectors and DDLs
1055 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1056 for (Int_t iDDL = 0; iDDL < kDetectorDDLs[iDet]; iDDL++) {
1058 Int_t ddlID = 0x100*iDet + iDDL;
1059 Int_t ldcID = Int_t(ldc + 0.0001);
1060 ldc += kDetectorLDCs[iDet] / kDetectorDDLs[iDet];
1062 char rawFileName[256];
1063 sprintf(rawFileName, "raw%d/%s_%d.ddl",
1064 iEvent, kDetectors[iDet], ddlID);
1066 // check existence and size of raw data file
1067 FILE* file = fopen(rawFileName, "rb");
1068 if (!file) continue;
1069 fseek(file, 0, SEEK_END);
1070 unsigned long size = ftell(file);
1072 if (!size) continue;
1074 if (ldcID != prevLDC) {
1075 fprintf(pipe, " LDC Id %d\n", ldcID);
1078 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1083 Int_t result = gSystem->ClosePipe(pipe);
1086 return (result == 0);
1089 //_____________________________________________________________________________
1090 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1091 const char* rootFileName)
1093 // convert a DATE file to a root file with the program "alimdc"
1096 const Int_t kDBSize = 1000000000;
1097 const Int_t kTagDBSize = 1000000000;
1098 const Bool_t kFilter = kFALSE;
1099 const Int_t kCompression = 1;
1101 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1103 AliError("the program alimdc was not found");
1104 if (fStopOnError) return kFALSE;
1109 AliInfo(Form("converting DATE file %s to root file %s",
1110 dateFileName, rootFileName));
1112 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1113 const char* tagDBFS = "/tmp/mdc1/tags";
1114 const char* runDBFS = "/tmp/mdc1/meta";
1116 // User defined file system locations
1117 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1118 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1119 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1120 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1121 if (gSystem->Getenv("ALIMDC_TAGDB"))
1122 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1123 if (gSystem->Getenv("ALIMDC_RUNDB"))
1124 runDBFS = gSystem->Getenv("ALIMDC_RUNDB");
1126 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1127 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1128 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1129 gSystem->Exec(Form("rm -rf %s",runDBFS));
1131 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1132 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1133 gSystem->Exec(Form("mkdir %s",tagDBFS));
1134 gSystem->Exec(Form("mkdir %s",runDBFS));
1136 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1137 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1138 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1140 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1141 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1142 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1143 gSystem->Exec(Form("rm -rf %s",runDBFS));
1145 return (result == 0);
1149 //_____________________________________________________________________________
1150 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1152 // delete existing run loaders, open a new one and load gAlice
1154 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1155 AliRunLoader* runLoader =
1156 AliRunLoader::Open(fGAliceFileName.Data(),
1157 AliConfig::GetDefaultEventFolderName(), mode);
1159 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1162 runLoader->LoadgAlice();
1163 gAlice = runLoader->GetAliRun();
1165 AliError(Form("no gAlice object found in file %s",
1166 fGAliceFileName.Data()));
1172 //_____________________________________________________________________________
1173 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1175 // get or calculate the number of signal events per background event
1177 if (!fBkgrdFileNames) return 1;
1178 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1179 if (nBkgrdFiles == 0) return 1;
1181 // get the number of signal events
1183 AliRunLoader* runLoader =
1184 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1185 if (!runLoader) return 1;
1186 nEvents = runLoader->GetNumberOfEvents();
1191 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1192 // get the number of background events
1193 const char* fileName = ((TObjString*)
1194 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1195 AliRunLoader* runLoader =
1196 AliRunLoader::Open(fileName, "BKGRD");
1197 if (!runLoader) continue;
1198 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1201 // get or calculate the number of signal per background events
1202 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1203 if (nSignalPerBkgrd <= 0) {
1204 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1205 } else if (result && (result != nSignalPerBkgrd)) {
1206 AliInfo(Form("the number of signal events per background event "
1207 "will be changed from %d to %d for stream %d",
1208 nSignalPerBkgrd, result, iBkgrdFile+1));
1209 nSignalPerBkgrd = result;
1212 if (!result) result = nSignalPerBkgrd;
1213 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1214 AliWarning(Form("not enough background events (%d) for %d signal events "
1215 "using %d signal per background events for stream %d",
1216 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1223 //_____________________________________________________________________________
1224 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1226 // check whether detName is contained in detectors
1227 // if yes, it is removed from detectors
1229 // check if all detectors are selected
1230 if ((detectors.CompareTo("ALL") == 0) ||
1231 detectors.BeginsWith("ALL ") ||
1232 detectors.EndsWith(" ALL") ||
1233 detectors.Contains(" ALL ")) {
1238 // search for the given detector
1239 Bool_t result = kFALSE;
1240 if ((detectors.CompareTo(detName) == 0) ||
1241 detectors.BeginsWith(detName+" ") ||
1242 detectors.EndsWith(" "+detName) ||
1243 detectors.Contains(" "+detName+" ")) {
1244 detectors.ReplaceAll(detName, "");
1248 // clean up the detectors string
1249 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1250 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1251 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);