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(TObjArray* alObjArray)
317 // Read collection of alignment objects (AliAlignObj derived) saved
318 // in the TClonesArray ClArrayName and apply them to the geometry
319 // manager singleton.
322 Int_t nvols = alObjArray->GetEntriesFast();
326 for(Int_t j=0; j<nvols; j++)
328 AliAlignObj* alobj = (AliAlignObj*) alObjArray->UncheckedAt(j);
329 if (alobj->ApplyToGeometry() == kFALSE) flag = kFALSE;
332 if (AliDebugLevelClass() >= 1) {
333 gGeoManager->CheckOverlaps(20);
334 TObjArray* ovexlist = gGeoManager->GetListOfOverlaps();
335 if(ovexlist->GetEntriesFast()){
336 AliError("The application of alignment objects to the geometry caused huge overlaps/extrusions!");
344 //_____________________________________________________________________________
345 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* fileName, const char* clArrayName)
347 // read collection of alignment objects (AliAlignObj derived) saved
348 // in the TClonesArray ClArrayName in the file fileName and apply
349 // them to the TGeo geometry passed as argument
352 TFile* inFile = TFile::Open(fileName,"READ");
353 if (!inFile || !inFile->IsOpen()) {
354 AliErrorClass(Form("Could not open file %s !",fileName));
358 TClonesArray* alObjArray = ((TClonesArray*) inFile->Get(clArrayName));
361 AliErrorClass(Form("Could not get array (%s) from file (%s) !",clArrayName,fileName));
365 return ApplyAlignObjsToGeom(alObjArray);
369 //_____________________________________________________________________________
370 Bool_t AliSimulation::ApplyAlignObjsToGeom(AliCDBParam* param, AliCDBId& Id)
372 // read collection of alignment objects (AliAlignObj derived) saved
373 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
374 // param (to get the AliCDBStorage) and Id; apply the alignment objects
375 // to the TGeo geometry passed as argument
378 AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage(param);
379 AliCDBEntry* entry = storage->Get(Id);
380 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
382 return ApplyAlignObjsToGeom(AlObjArray);
386 //_____________________________________________________________________________
387 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* uri, const char* path, Int_t runnum, Int_t version, Int_t sversion)
389 // read collection of alignment objects (AliAlignObj derived) saved
390 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
391 // param (to get the AliCDBStorage) and Id; apply the alignment objects
392 // to the TGeo geometry passed as argument
395 AliCDBParam* param = AliCDBManager::Instance()->CreateParameter(uri);
396 AliCDBId id(path, runnum, runnum, version, sversion);
398 return ApplyAlignObjsToGeom(param, id);
402 //_____________________________________________________________________________
403 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* detName, Int_t runnum, Int_t version, Int_t sversion)
405 // read collection of alignment objects (AliAlignObj derived) saved
406 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
407 // param (to get the AliCDBStorage) and Id; apply the alignment objects
408 // to the TGeo geometry passed as argument
411 InitCDBStorage("local://$ALICE_ROOT");
412 AliCDBPath path(detName,"Align","Data");
413 AliCDBEntry* entry = AliCDBManager::Instance()->Get(path.GetPath(),runnum,version,sversion);
415 if(!entry) return kFALSE;
416 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
418 return ApplyAlignObjsToGeom(AlObjArray);
421 //_____________________________________________________________________________
422 Bool_t AliSimulation::SetAlignObjArraySingleDet(const char* detName)
424 // Fills array of single detector's alignable objects from CDB
426 AliDebug(2, Form("Loading alignment data for detector: %s",detName));
430 AliCDBPath path(detName,"Align","Data");
432 entry=AliCDBManager::Instance()->Get(path.GetPath());
434 AliDebug(2,Form("Couldn't load alignment data for detector %s",detName));
438 TClonesArray *alignArray = (TClonesArray*) entry->GetObject();
439 alignArray->SetOwner(0);
440 AliDebug(2,Form("Found %d alignment objects for %s",
441 alignArray->GetEntries(),detName));
443 AliAlignObj *alignObj=0;
444 TIter iter(alignArray);
446 // loop over align objects in detector
447 while( ( alignObj=(AliAlignObj *) iter.Next() ) ){
448 fAlignObjArray->Add(alignObj);
450 // delete entry --- Don't delete, it is cached!
452 AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
457 //_____________________________________________________________________________
458 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
460 // Read the alignment objects from CDB.
461 // Each detector is supposed to have the
462 // alignment objects in DET/Align/Data CDB path.
463 // All the detector objects are then collected,
464 // sorted by geometry level (starting from ALIC) and
465 // then applied to the TGeo geometry.
466 // Finally an overlaps check is performed.
468 Bool_t delRunLoader = kFALSE;
470 runLoader = LoadRun("READ");
471 if (!runLoader) return kFALSE;
472 delRunLoader = kTRUE;
475 // Load alignment data from CDB and fill fAlignObjArray
476 if(fLoadAlignFromCDB){
477 if(!fAlignObjArray) fAlignObjArray = new TObjArray();
479 //fAlignObjArray->RemoveAll();
480 fAlignObjArray->Clear();
481 fAlignObjArray->SetOwner(0);
483 TString detStr = fLoadAlignData;
484 TString dataNotLoaded="";
485 TString dataLoaded="";
487 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
488 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
489 AliModule* det = (AliModule*) detArray->At(iDet);
490 if (!det || !det->IsActive()) continue;
491 if (IsSelected(det->GetName(), detStr)) {
492 if(!SetAlignObjArraySingleDet(det->GetName())){
493 dataNotLoaded += det->GetName();
494 dataNotLoaded += " ";
496 dataLoaded += det->GetName();
500 } // end loop over detectors
502 if ((detStr.CompareTo("ALL") == 0)) detStr = "";
503 dataNotLoaded += detStr;
504 AliInfo(Form("Alignment data loaded for: %s",
506 AliInfo(Form("Didn't/couldn't load alignment data for: %s",
507 dataNotLoaded.Data()));
508 } // fLoadAlignFromCDB flag
510 // Check if the array with alignment objects was
511 // provided by the user. If yes, apply the objects
512 // to the present TGeo geometry
513 if (fAlignObjArray) {
514 if (gGeoManager && gGeoManager->IsClosed()) {
515 if (ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
516 AliError("The misalignment of one or more volumes failed!"
517 "Compare the list of simulated detectors and the list of detector alignment data!");
518 if (delRunLoader) delete runLoader;
523 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
524 if (delRunLoader) delete runLoader;
529 if (delRunLoader) delete runLoader;
535 //_____________________________________________________________________________
536 Bool_t AliSimulation::SetRunNumber()
538 // Set the CDB manager run number
539 // The run number is retrieved from gAlice
541 if(AliCDBManager::Instance()->GetRun() < 0) {
542 AliRunLoader* runLoader = LoadRun("READ");
543 if (!runLoader) return kFALSE;
545 AliCDBManager::Instance()->SetRun(runLoader->GetAliRun()->GetRunNumber());
546 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
553 //_____________________________________________________________________________
554 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
556 // add a file with background events for merging
558 TObjString* fileNameStr = new TObjString(fileName);
559 fileNameStr->SetUniqueID(nSignalPerBkgrd);
560 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
561 fBkgrdFileNames->Add(fileNameStr);
565 //_____________________________________________________________________________
566 Bool_t AliSimulation::Run(Int_t nEvents)
568 // run the generation, simulation and digitization
570 InitCDBStorage(fCDBUri);
572 if (nEvents > 0) fNEvents = nEvents;
574 // generation and simulation -> hits
575 if (fRunGeneration) {
576 if (!RunSimulation()) if (fStopOnError) return kFALSE;
579 // Set run number in CDBManager (if it is not already set in RunSimulation)
580 if (!SetRunNumber()) if (fStopOnError) return kFALSE;
582 // Load and misalign the geometry
584 TGeoManager::Import("geometry.root");
585 if (!gGeoManager) if (fStopOnError) return kFALSE;
586 if (!MisalignGeometry()) if (fStopOnError) return kFALSE;
589 // hits -> summable digits
590 if (!fMakeSDigits.IsNull()) {
591 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
594 // summable digits -> digits
595 if (!fMakeDigits.IsNull()) {
596 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
597 if (fStopOnError) return kFALSE;
602 if (!fMakeDigitsFromHits.IsNull()) {
603 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
604 AliWarning(Form("Merging and direct creation of digits from hits "
605 "was selected for some detectors. "
606 "No merging will be done for the following detectors: %s",
607 fMakeDigitsFromHits.Data()));
609 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
610 if (fStopOnError) return kFALSE;
615 if (!fMakeTrigger.IsNull()) {
616 if (!RunTrigger(fMakeTrigger)) {
617 if (fStopOnError) return kFALSE;
621 // digits -> raw data
622 if (!fWriteRawData.IsNull()) {
623 if (!WriteRawData(fWriteRawData, fRawDataFileName,
624 fDeleteIntermediateFiles)) {
625 if (fStopOnError) return kFALSE;
632 //_____________________________________________________________________________
633 Bool_t AliSimulation::RunTrigger(const char* descriptors)
637 TStopwatch stopwatch;
640 AliRunLoader* runLoader = LoadRun("READ");
641 if (!runLoader) return kFALSE;
642 TString des = descriptors;
644 runLoader->MakeTree( "CT" );
645 AliCentralTrigger* aCTP = runLoader->GetTrigger();
647 aCTP->LoadDescriptor( des );
650 if( !aCTP->RunTrigger( runLoader ) ) {
657 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
658 stopwatch.RealTime(),stopwatch.CpuTime()));
667 //_____________________________________________________________________________
668 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
670 // run the generation and simulation
672 TStopwatch stopwatch;
676 AliError("no gAlice object. Restart aliroot and try again.");
679 if (gAlice->Modules()->GetEntries() > 0) {
680 AliError("gAlice was already run. Restart aliroot and try again.");
684 AliInfo(Form("initializing gAlice with config file %s",
685 fConfigFileName.Data()));
686 StdoutToAliInfo(StderrToAliError(
687 gAlice->Init(fConfigFileName.Data());
690 // Set run number in CDBManager
691 AliCDBManager::Instance()->SetRun(gAlice->GetRunNumber());
692 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
694 AliRunLoader* runLoader = gAlice->GetRunLoader();
696 AliError(Form("gAlice has no run loader object. "
697 "Check your config file: %s", fConfigFileName.Data()));
700 SetGAliceFile(runLoader->GetFileName());
702 // Export ideal geometry
703 if (gGeoManager) gGeoManager->Export("geometry.root");
706 if (!MisalignGeometry(runLoader)) return kFALSE;
708 // Temporary fix by A.Gheata
709 // Could be removed with the next Root version (>5.11)
711 TIter next(gGeoManager->GetListOfVolumes());
713 while ((vol = (TGeoVolume *)next())) {
714 if (vol->GetVoxels()) {
715 if (vol->GetVoxels()->NeedRebuild()) {
716 vol->GetVoxels()->Voxelize();
723 // Export (mis)aligned geometry
724 if (gGeoManager) gGeoManager->Export("misaligned_geometry.root");
726 // AliRunLoader* runLoader = gAlice->GetRunLoader();
728 // AliError(Form("gAlice has no run loader object. "
729 // "Check your config file: %s", fConfigFileName.Data()));
732 // SetGAliceFile(runLoader->GetFileName());
734 if (!gAlice->Generator()) {
735 AliError(Form("gAlice has no generator object. "
736 "Check your config file: %s", fConfigFileName.Data()));
739 if (nEvents <= 0) nEvents = fNEvents;
741 // get vertex from background file in case of merging
742 if (fUseBkgrdVertex &&
743 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
744 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
745 const char* fileName = ((TObjString*)
746 (fBkgrdFileNames->At(0)))->GetName();
747 AliInfo(Form("The vertex will be taken from the background "
748 "file %s with nSignalPerBackground = %d",
749 fileName, signalPerBkgrd));
750 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
751 gAlice->Generator()->SetVertexGenerator(vtxGen);
754 if (!fRunSimulation) {
755 gAlice->Generator()->SetTrackingFlag(0);
758 // set the number of events per file for given detectors and data types
759 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
760 if (!fEventsPerFile[i]) continue;
761 const char* detName = fEventsPerFile[i]->GetName();
762 const char* typeName = fEventsPerFile[i]->GetTitle();
763 TString loaderName(detName);
764 loaderName += "Loader";
765 AliLoader* loader = runLoader->GetLoader(loaderName);
767 AliError(Form("RunSimulation", "no loader for %s found\n"
768 "Number of events per file not set for %s %s",
769 detName, typeName, detName));
772 AliDataLoader* dataLoader =
773 loader->GetDataLoader(typeName);
775 AliError(Form("no data loader for %s found\n"
776 "Number of events per file not set for %s %s",
777 typeName, detName, typeName));
780 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
781 AliDebug(1, Form("number of events per file set to %d for %s %s",
782 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
785 AliInfo("running gAlice");
786 StdoutToAliInfo(StderrToAliError(
787 gAlice->Run(nEvents);
792 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
793 stopwatch.RealTime(),stopwatch.CpuTime()));
798 //_____________________________________________________________________________
799 Bool_t AliSimulation::RunSDigitization(const char* detectors)
801 // run the digitization and produce summable digits
803 TStopwatch stopwatch;
806 AliRunLoader* runLoader = LoadRun();
807 if (!runLoader) return kFALSE;
809 TString detStr = detectors;
810 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
811 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
812 AliModule* det = (AliModule*) detArray->At(iDet);
813 if (!det || !det->IsActive()) continue;
814 if (IsSelected(det->GetName(), detStr)) {
815 AliInfo(Form("creating summable digits for %s", det->GetName()));
816 TStopwatch stopwatchDet;
817 stopwatchDet.Start();
819 AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
820 det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
824 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
825 AliError(Form("the following detectors were not found: %s",
827 if (fStopOnError) return kFALSE;
832 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
833 stopwatch.RealTime(),stopwatch.CpuTime()));
839 //_____________________________________________________________________________
840 Bool_t AliSimulation::RunDigitization(const char* detectors,
841 const char* excludeDetectors)
843 // run the digitization and produce digits from sdigits
845 TStopwatch stopwatch;
848 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
849 if (gAlice) delete gAlice;
853 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
854 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
855 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
856 manager->SetInputStream(0, fGAliceFileName.Data());
857 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
858 const char* fileName = ((TObjString*)
859 (fBkgrdFileNames->At(iStream-1)))->GetName();
860 manager->SetInputStream(iStream, fileName);
863 TString detStr = detectors;
864 TString detExcl = excludeDetectors;
865 manager->GetInputStream(0)->ImportgAlice();
866 AliRunLoader* runLoader =
867 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
868 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
869 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
870 AliModule* det = (AliModule*) detArray->At(iDet);
871 if (!det || !det->IsActive()) continue;
872 if (IsSelected(det->GetName(), detStr) &&
873 !IsSelected(det->GetName(), detExcl)) {
874 AliDigitizer* digitizer = det->CreateDigitizer(manager);
876 AliError(Form("no digitizer for %s", det->GetName()));
877 if (fStopOnError) return kFALSE;
879 digitizer->SetRegionOfInterest(fRegionOfInterest);
884 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
885 AliError(Form("the following detectors were not found: %s",
887 if (fStopOnError) return kFALSE;
890 if (!manager->GetListOfTasks()->IsEmpty()) {
891 AliInfo("executing digitization");
897 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
898 stopwatch.RealTime(),stopwatch.CpuTime()));
903 //_____________________________________________________________________________
904 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
906 // run the digitization and produce digits from hits
908 TStopwatch stopwatch;
911 AliRunLoader* runLoader = LoadRun("READ");
912 if (!runLoader) return kFALSE;
914 TString detStr = detectors;
915 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
916 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
917 AliModule* det = (AliModule*) detArray->At(iDet);
918 if (!det || !det->IsActive()) continue;
919 if (IsSelected(det->GetName(), detStr)) {
920 AliInfo(Form("creating digits from hits for %s", det->GetName()));
925 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
926 AliError(Form("the following detectors were not found: %s",
928 if (fStopOnError) return kFALSE;
932 //PH Temporary fix to avoid interference with the PHOS loder/getter
933 //PH The problem has to be solved in more general way 09/06/05
935 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
936 stopwatch.RealTime(),stopwatch.CpuTime()));
941 //_____________________________________________________________________________
942 Bool_t AliSimulation::WriteRawData(const char* detectors,
943 const char* fileName,
944 Bool_t deleteIntermediateFiles)
946 // convert the digits to raw data
947 // First DDL raw data files for the given detectors are created.
948 // If a file name is given, the DDL files are then converted to a DATE file.
949 // If deleteIntermediateFiles is true, the DDL raw files are deleted
951 // If the file name has the extension ".root", the DATE file is converted
953 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
955 TStopwatch stopwatch;
958 if (!WriteRawFiles(detectors)) {
959 if (fStopOnError) return kFALSE;
962 TString dateFileName(fileName);
963 if (!dateFileName.IsNull()) {
964 Bool_t rootOutput = dateFileName.EndsWith(".root");
965 if (rootOutput) dateFileName += ".date";
966 if (!ConvertRawFilesToDate(dateFileName)) {
967 if (fStopOnError) return kFALSE;
969 if (deleteIntermediateFiles) {
970 AliRunLoader* runLoader = LoadRun("READ");
971 if (runLoader) for (Int_t iEvent = 0;
972 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
974 sprintf(command, "rm -r raw%d", iEvent);
975 gSystem->Exec(command);
980 if (!ConvertDateToRoot(dateFileName, fileName)) {
981 if (fStopOnError) return kFALSE;
983 if (deleteIntermediateFiles) {
984 gSystem->Unlink(dateFileName);
989 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
990 stopwatch.RealTime(),stopwatch.CpuTime()));
995 //_____________________________________________________________________________
996 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
998 // convert the digits to raw data DDL files
1000 AliRunLoader* runLoader = LoadRun("READ");
1001 if (!runLoader) return kFALSE;
1003 // write raw data to DDL files
1004 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1005 AliInfo(Form("processing event %d", iEvent));
1006 runLoader->GetEvent(iEvent);
1007 TString baseDir = gSystem->WorkingDirectory();
1009 sprintf(dirName, "raw%d", iEvent);
1010 gSystem->MakeDirectory(dirName);
1011 if (!gSystem->ChangeDirectory(dirName)) {
1012 AliError(Form("couldn't change to directory %s", dirName));
1013 if (fStopOnError) return kFALSE; else continue;
1016 TString detStr = detectors;
1017 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1018 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1019 AliModule* det = (AliModule*) detArray->At(iDet);
1020 if (!det || !det->IsActive()) continue;
1021 if (IsSelected(det->GetName(), detStr)) {
1022 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1027 gSystem->ChangeDirectory(baseDir);
1028 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1029 AliError(Form("the following detectors were not found: %s",
1031 if (fStopOnError) return kFALSE;
1039 //_____________________________________________________________________________
1040 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1042 // convert raw data DDL files to a DATE file with the program "dateStream"
1044 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1046 AliError("the program dateStream was not found");
1047 if (fStopOnError) return kFALSE;
1052 AliRunLoader* runLoader = LoadRun("READ");
1053 if (!runLoader) return kFALSE;
1055 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1057 sprintf(command, "dateStream -o %s -# %d -C",
1058 dateFileName, runLoader->GetNumberOfEvents());
1059 FILE* pipe = gSystem->OpenPipe(command, "w");
1061 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1062 fprintf(pipe, "GDC\n");
1066 // loop over detectors and DDLs
1067 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1068 for (Int_t iDDL = 0; iDDL < kDetectorDDLs[iDet]; iDDL++) {
1070 Int_t ddlID = 0x100*iDet + iDDL;
1071 Int_t ldcID = Int_t(ldc + 0.0001);
1072 ldc += kDetectorLDCs[iDet] / kDetectorDDLs[iDet];
1074 char rawFileName[256];
1075 sprintf(rawFileName, "raw%d/%s_%d.ddl",
1076 iEvent, kDetectors[iDet], ddlID);
1078 // check existence and size of raw data file
1079 FILE* file = fopen(rawFileName, "rb");
1080 if (!file) continue;
1081 fseek(file, 0, SEEK_END);
1082 unsigned long size = ftell(file);
1084 if (!size) continue;
1086 if (ldcID != prevLDC) {
1087 fprintf(pipe, " LDC Id %d\n", ldcID);
1090 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1095 Int_t result = gSystem->ClosePipe(pipe);
1098 return (result == 0);
1101 //_____________________________________________________________________________
1102 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1103 const char* rootFileName)
1105 // convert a DATE file to a root file with the program "alimdc"
1108 const Int_t kDBSize = 1000000000;
1109 const Int_t kTagDBSize = 1000000000;
1110 const Bool_t kFilter = kFALSE;
1111 const Int_t kCompression = 1;
1113 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1115 AliError("the program alimdc was not found");
1116 if (fStopOnError) return kFALSE;
1121 AliInfo(Form("converting DATE file %s to root file %s",
1122 dateFileName, rootFileName));
1124 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1125 const char* tagDBFS = "/tmp/mdc1/tags";
1126 const char* runDBFS = "/tmp/mdc1/meta";
1128 // User defined file system locations
1129 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1130 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1131 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1132 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1133 if (gSystem->Getenv("ALIMDC_TAGDB"))
1134 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1135 if (gSystem->Getenv("ALIMDC_RUNDB"))
1136 runDBFS = gSystem->Getenv("ALIMDC_RUNDB");
1138 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1139 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1140 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1141 gSystem->Exec(Form("rm -rf %s",runDBFS));
1143 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1144 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1145 gSystem->Exec(Form("mkdir %s",tagDBFS));
1146 gSystem->Exec(Form("mkdir %s",runDBFS));
1148 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1149 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1150 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1152 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1153 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1154 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1155 gSystem->Exec(Form("rm -rf %s",runDBFS));
1157 return (result == 0);
1161 //_____________________________________________________________________________
1162 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1164 // delete existing run loaders, open a new one and load gAlice
1166 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1167 AliRunLoader* runLoader =
1168 AliRunLoader::Open(fGAliceFileName.Data(),
1169 AliConfig::GetDefaultEventFolderName(), mode);
1171 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1174 runLoader->LoadgAlice();
1175 gAlice = runLoader->GetAliRun();
1177 AliError(Form("no gAlice object found in file %s",
1178 fGAliceFileName.Data()));
1184 //_____________________________________________________________________________
1185 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1187 // get or calculate the number of signal events per background event
1189 if (!fBkgrdFileNames) return 1;
1190 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1191 if (nBkgrdFiles == 0) return 1;
1193 // get the number of signal events
1195 AliRunLoader* runLoader =
1196 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1197 if (!runLoader) return 1;
1198 nEvents = runLoader->GetNumberOfEvents();
1203 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1204 // get the number of background events
1205 const char* fileName = ((TObjString*)
1206 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1207 AliRunLoader* runLoader =
1208 AliRunLoader::Open(fileName, "BKGRD");
1209 if (!runLoader) continue;
1210 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1213 // get or calculate the number of signal per background events
1214 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1215 if (nSignalPerBkgrd <= 0) {
1216 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1217 } else if (result && (result != nSignalPerBkgrd)) {
1218 AliInfo(Form("the number of signal events per background event "
1219 "will be changed from %d to %d for stream %d",
1220 nSignalPerBkgrd, result, iBkgrdFile+1));
1221 nSignalPerBkgrd = result;
1224 if (!result) result = nSignalPerBkgrd;
1225 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1226 AliWarning(Form("not enough background events (%d) for %d signal events "
1227 "using %d signal per background events for stream %d",
1228 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1235 //_____________________________________________________________________________
1236 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1238 // check whether detName is contained in detectors
1239 // if yes, it is removed from detectors
1241 // check if all detectors are selected
1242 if ((detectors.CompareTo("ALL") == 0) ||
1243 detectors.BeginsWith("ALL ") ||
1244 detectors.EndsWith(" ALL") ||
1245 detectors.Contains(" ALL ")) {
1250 // search for the given detector
1251 Bool_t result = kFALSE;
1252 if ((detectors.CompareTo(detName) == 0) ||
1253 detectors.BeginsWith(detName+" ") ||
1254 detectors.EndsWith(" "+detName) ||
1255 detectors.Contains(" "+detName+" ")) {
1256 detectors.ReplaceAll(detName, "");
1260 // clean up the detectors string
1261 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1262 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1263 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);