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();
324 for(Int_t j=0; j<nvols; j++)
326 AliAlignObj* alobj = (AliAlignObj*) alObjArray->UncheckedAt(j);
327 if (alobj->ApplyToGeometry() == kFALSE)
331 if (AliDebugLevelClass() >= 1) {
332 gGeoManager->CheckOverlaps(20);
333 TObjArray* ovexlist = gGeoManager->GetListOfOverlaps();
334 if(ovexlist->GetEntriesFast()){
335 AliError("The application of alignment objects to the geometry caused huge overlaps/extrusions!");
343 //_____________________________________________________________________________
344 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* fileName, const char* clArrayName)
346 // read collection of alignment objects (AliAlignObj derived) saved
347 // in the TClonesArray ClArrayName in the file fileName and apply
348 // them to the TGeo geometry passed as argument
351 TFile* inFile = TFile::Open(fileName,"READ");
352 if (!inFile || !inFile->IsOpen()) {
353 AliErrorClass(Form("Could not open file %s !",fileName));
357 TClonesArray* alObjArray = ((TClonesArray*) inFile->Get(clArrayName));
360 AliErrorClass(Form("Could not get array (%s) from file (%s) !",clArrayName,fileName));
364 return ApplyAlignObjsToGeom(alObjArray);
368 //_____________________________________________________________________________
369 Bool_t AliSimulation::ApplyAlignObjsToGeom(AliCDBParam* param, AliCDBId& Id)
371 // read collection of alignment objects (AliAlignObj derived) saved
372 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
373 // param (to get the AliCDBStorage) and Id; apply the alignment objects
374 // to the TGeo geometry passed as argument
377 AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage(param);
378 AliCDBEntry* entry = storage->Get(Id);
379 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
381 return ApplyAlignObjsToGeom(AlObjArray);
385 //_____________________________________________________________________________
386 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* uri, const char* path, Int_t runnum, Int_t version, Int_t sversion)
388 // read collection of alignment objects (AliAlignObj derived) saved
389 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
390 // param (to get the AliCDBStorage) and Id; apply the alignment objects
391 // to the TGeo geometry passed as argument
394 AliCDBParam* param = AliCDBManager::Instance()->CreateParameter(uri);
395 AliCDBId id(path, runnum, runnum, version, sversion);
397 return ApplyAlignObjsToGeom(param, id);
401 //_____________________________________________________________________________
402 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* detName, Int_t runnum, Int_t version, Int_t sversion)
404 // read collection of alignment objects (AliAlignObj derived) saved
405 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
406 // param (to get the AliCDBStorage) and Id; apply the alignment objects
407 // to the TGeo geometry passed as argument
410 InitCDBStorage("local://$ALICE_ROOT");
411 AliCDBPath path(detName,"Align","Data");
412 AliCDBEntry* entry = AliCDBManager::Instance()->Get(path.GetPath(),runnum,version,sversion);
414 if(!entry) return kFALSE;
415 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
417 return ApplyAlignObjsToGeom(AlObjArray);
420 //_____________________________________________________________________________
421 Bool_t AliSimulation::SetAlignObjArraySingleDet(const char* detName)
423 // Fills array of single detector's alignable objects from CDB
425 AliDebug(2, Form("Loading alignment data for detector: %s",detName));
429 AliCDBPath path(detName,"Align","Data");
431 entry=AliCDBManager::Instance()->Get(path.GetPath());
433 AliDebug(2,Form("Couldn't load alignment data for detector %s",detName));
437 TClonesArray *alignArray = (TClonesArray*) entry->GetObject();
438 alignArray->SetOwner(0);
439 AliDebug(2,Form("Found %d alignment objects for %s",
440 alignArray->GetEntries(),detName));
442 AliAlignObj *alignObj=0;
443 TIter iter(alignArray);
445 // loop over align objects in detector
446 while( ( alignObj=(AliAlignObj *) iter.Next() ) ){
447 fAlignObjArray->Add(alignObj);
449 // delete entry --- Don't delete, it is cached!
451 AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
456 //_____________________________________________________________________________
457 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
459 // Read the alignment objects from CDB.
460 // Each detector is supposed to have the
461 // alignment objects in DET/Align/Data CDB path.
462 // All the detector objects are then collected,
463 // sorted by geometry level (starting from ALIC) and
464 // then applied to the TGeo geometry.
465 // Finally an overlaps check is performed.
467 Bool_t delRunLoader = kFALSE;
469 runLoader = LoadRun("READ");
470 if (!runLoader) return kFALSE;
471 delRunLoader = kTRUE;
474 // Load alignment data from CDB and fill fAlignObjArray
475 if(fLoadAlignFromCDB){
476 if(!fAlignObjArray) fAlignObjArray = new TObjArray();
478 //fAlignObjArray->RemoveAll();
479 fAlignObjArray->Clear();
480 fAlignObjArray->SetOwner(0);
482 TString detStr = fLoadAlignData;
483 TString dataNotLoaded="";
484 TString dataLoaded="";
486 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
487 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
488 AliModule* det = (AliModule*) detArray->At(iDet);
489 if (!det || !det->IsActive()) continue;
490 if (IsSelected(det->GetName(), detStr)) {
491 if(!SetAlignObjArraySingleDet(det->GetName())){
492 dataNotLoaded += det->GetName();
493 dataNotLoaded += " ";
495 dataLoaded += det->GetName();
499 } // end loop over detectors
501 if ((detStr.CompareTo("ALL") == 0)) detStr = "";
502 dataNotLoaded += detStr;
503 AliInfo(Form("Alignment data loaded for: %s",
505 AliInfo(Form("Didn't/couldn't load alignment data for: %s",
506 dataNotLoaded.Data()));
507 } // fLoadAlignFromCDB flag
509 // Check if the array with alignment objects was
510 // provided by the user. If yes, apply the objects
511 // to the present TGeo geometry
512 if (fAlignObjArray) {
513 if (gGeoManager && gGeoManager->IsClosed()) {
514 if (ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
515 AliError("The application of misalignment failed! Restart aliroot and try again. ");
520 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
525 if (delRunLoader) delete runLoader;
531 //_____________________________________________________________________________
532 Bool_t AliSimulation::SetRunNumber()
534 // Set the CDB manager run number
535 // The run number is retrieved from gAlice
537 if(AliCDBManager::Instance()->GetRun() < 0) {
538 AliRunLoader* runLoader = LoadRun("READ");
539 if (!runLoader) return kFALSE;
541 AliCDBManager::Instance()->SetRun(runLoader->GetAliRun()->GetRunNumber());
542 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
549 //_____________________________________________________________________________
550 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
552 // add a file with background events for merging
554 TObjString* fileNameStr = new TObjString(fileName);
555 fileNameStr->SetUniqueID(nSignalPerBkgrd);
556 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
557 fBkgrdFileNames->Add(fileNameStr);
561 //_____________________________________________________________________________
562 Bool_t AliSimulation::Run(Int_t nEvents)
564 // run the generation, simulation and digitization
566 InitCDBStorage(fCDBUri);
568 if (nEvents > 0) fNEvents = nEvents;
570 // generation and simulation -> hits
571 if (fRunGeneration) {
572 if (!RunSimulation()) if (fStopOnError) return kFALSE;
575 // Set run number in CDBManager (if it is not already set in RunSimulation)
576 if (!SetRunNumber()) if (fStopOnError) return kFALSE;
578 // Load and misalign the geometry
580 TGeoManager::Import("geometry.root");
581 if (!gGeoManager) if (fStopOnError) return kFALSE;
582 if (!MisalignGeometry()) if (fStopOnError) return kFALSE;
585 // hits -> summable digits
586 if (!fMakeSDigits.IsNull()) {
587 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
590 // summable digits -> digits
591 if (!fMakeDigits.IsNull()) {
592 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
593 if (fStopOnError) return kFALSE;
598 if (!fMakeDigitsFromHits.IsNull()) {
599 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
600 AliWarning(Form("Merging and direct creation of digits from hits "
601 "was selected for some detectors. "
602 "No merging will be done for the following detectors: %s",
603 fMakeDigitsFromHits.Data()));
605 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
606 if (fStopOnError) return kFALSE;
611 if (!fMakeTrigger.IsNull()) {
612 if (!RunTrigger(fMakeTrigger)) {
613 if (fStopOnError) return kFALSE;
617 // digits -> raw data
618 if (!fWriteRawData.IsNull()) {
619 if (!WriteRawData(fWriteRawData, fRawDataFileName,
620 fDeleteIntermediateFiles)) {
621 if (fStopOnError) return kFALSE;
628 //_____________________________________________________________________________
629 Bool_t AliSimulation::RunTrigger(const char* descriptors)
633 TStopwatch stopwatch;
636 AliRunLoader* runLoader = LoadRun("READ");
637 if (!runLoader) return kFALSE;
638 TString des = descriptors;
640 runLoader->MakeTree( "CT" );
641 AliCentralTrigger* aCTP = runLoader->GetTrigger();
643 aCTP->LoadDescriptor( des );
646 if( !aCTP->RunTrigger( runLoader ) ) {
653 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
654 stopwatch.RealTime(),stopwatch.CpuTime()));
663 //_____________________________________________________________________________
664 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
666 // run the generation and simulation
668 TStopwatch stopwatch;
672 AliError("no gAlice object. Restart aliroot and try again.");
675 if (gAlice->Modules()->GetEntries() > 0) {
676 AliError("gAlice was already run. Restart aliroot and try again.");
680 AliInfo(Form("initializing gAlice with config file %s",
681 fConfigFileName.Data()));
682 StdoutToAliInfo(StderrToAliError(
683 gAlice->Init(fConfigFileName.Data());
686 // Set run number in CDBManager
687 AliCDBManager::Instance()->SetRun(gAlice->GetRunNumber());
688 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
690 AliRunLoader* runLoader = gAlice->GetRunLoader();
692 AliError(Form("gAlice has no run loader object. "
693 "Check your config file: %s", fConfigFileName.Data()));
696 SetGAliceFile(runLoader->GetFileName());
698 // Export ideal geometry
699 if (gGeoManager) gGeoManager->Export("geometry.root");
702 if (!MisalignGeometry(runLoader)) return kFALSE;
704 // Temporary fix by A.Gheata
705 // Could be removed with the next Root version (>5.11)
707 TIter next(gGeoManager->GetListOfVolumes());
709 while ((vol = (TGeoVolume *)next())) {
710 if (vol->GetVoxels()) {
711 if (vol->GetVoxels()->NeedRebuild()) {
712 vol->GetVoxels()->Voxelize();
719 // Export (mis)aligned geometry
720 if (gGeoManager) gGeoManager->Export("misaligned_geometry.root");
722 // AliRunLoader* runLoader = gAlice->GetRunLoader();
724 // AliError(Form("gAlice has no run loader object. "
725 // "Check your config file: %s", fConfigFileName.Data()));
728 // SetGAliceFile(runLoader->GetFileName());
730 if (!gAlice->Generator()) {
731 AliError(Form("gAlice has no generator object. "
732 "Check your config file: %s", fConfigFileName.Data()));
735 if (nEvents <= 0) nEvents = fNEvents;
737 // get vertex from background file in case of merging
738 if (fUseBkgrdVertex &&
739 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
740 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
741 const char* fileName = ((TObjString*)
742 (fBkgrdFileNames->At(0)))->GetName();
743 AliInfo(Form("The vertex will be taken from the background "
744 "file %s with nSignalPerBackground = %d",
745 fileName, signalPerBkgrd));
746 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
747 gAlice->Generator()->SetVertexGenerator(vtxGen);
750 if (!fRunSimulation) {
751 gAlice->Generator()->SetTrackingFlag(0);
754 // set the number of events per file for given detectors and data types
755 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
756 if (!fEventsPerFile[i]) continue;
757 const char* detName = fEventsPerFile[i]->GetName();
758 const char* typeName = fEventsPerFile[i]->GetTitle();
759 TString loaderName(detName);
760 loaderName += "Loader";
761 AliLoader* loader = runLoader->GetLoader(loaderName);
763 AliError(Form("RunSimulation", "no loader for %s found\n"
764 "Number of events per file not set for %s %s",
765 detName, typeName, detName));
768 AliDataLoader* dataLoader =
769 loader->GetDataLoader(typeName);
771 AliError(Form("no data loader for %s found\n"
772 "Number of events per file not set for %s %s",
773 typeName, detName, typeName));
776 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
777 AliDebug(1, Form("number of events per file set to %d for %s %s",
778 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
781 AliInfo("running gAlice");
782 StdoutToAliInfo(StderrToAliError(
783 gAlice->Run(nEvents);
788 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
789 stopwatch.RealTime(),stopwatch.CpuTime()));
794 //_____________________________________________________________________________
795 Bool_t AliSimulation::RunSDigitization(const char* detectors)
797 // run the digitization and produce summable digits
799 TStopwatch stopwatch;
802 AliRunLoader* runLoader = LoadRun();
803 if (!runLoader) return kFALSE;
805 TString detStr = detectors;
806 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
807 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
808 AliModule* det = (AliModule*) detArray->At(iDet);
809 if (!det || !det->IsActive()) continue;
810 if (IsSelected(det->GetName(), detStr)) {
811 AliInfo(Form("creating summable digits for %s", det->GetName()));
812 TStopwatch stopwatchDet;
813 stopwatchDet.Start();
815 AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
816 det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
820 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
821 AliError(Form("the following detectors were not found: %s",
823 if (fStopOnError) return kFALSE;
828 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
829 stopwatch.RealTime(),stopwatch.CpuTime()));
835 //_____________________________________________________________________________
836 Bool_t AliSimulation::RunDigitization(const char* detectors,
837 const char* excludeDetectors)
839 // run the digitization and produce digits from sdigits
841 TStopwatch stopwatch;
844 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
845 if (gAlice) delete gAlice;
849 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
850 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
851 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
852 manager->SetInputStream(0, fGAliceFileName.Data());
853 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
854 const char* fileName = ((TObjString*)
855 (fBkgrdFileNames->At(iStream-1)))->GetName();
856 manager->SetInputStream(iStream, fileName);
859 TString detStr = detectors;
860 TString detExcl = excludeDetectors;
861 manager->GetInputStream(0)->ImportgAlice();
862 AliRunLoader* runLoader =
863 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
864 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
865 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
866 AliModule* det = (AliModule*) detArray->At(iDet);
867 if (!det || !det->IsActive()) continue;
868 if (IsSelected(det->GetName(), detStr) &&
869 !IsSelected(det->GetName(), detExcl)) {
870 AliDigitizer* digitizer = det->CreateDigitizer(manager);
872 AliError(Form("no digitizer for %s", det->GetName()));
873 if (fStopOnError) return kFALSE;
875 digitizer->SetRegionOfInterest(fRegionOfInterest);
880 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
881 AliError(Form("the following detectors were not found: %s",
883 if (fStopOnError) return kFALSE;
886 if (!manager->GetListOfTasks()->IsEmpty()) {
887 AliInfo("executing digitization");
893 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
894 stopwatch.RealTime(),stopwatch.CpuTime()));
899 //_____________________________________________________________________________
900 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
902 // run the digitization and produce digits from hits
904 TStopwatch stopwatch;
907 AliRunLoader* runLoader = LoadRun("READ");
908 if (!runLoader) return kFALSE;
910 TString detStr = detectors;
911 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
912 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
913 AliModule* det = (AliModule*) detArray->At(iDet);
914 if (!det || !det->IsActive()) continue;
915 if (IsSelected(det->GetName(), detStr)) {
916 AliInfo(Form("creating digits from hits for %s", det->GetName()));
921 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
922 AliError(Form("the following detectors were not found: %s",
924 if (fStopOnError) return kFALSE;
928 //PH Temporary fix to avoid interference with the PHOS loder/getter
929 //PH The problem has to be solved in more general way 09/06/05
931 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
932 stopwatch.RealTime(),stopwatch.CpuTime()));
937 //_____________________________________________________________________________
938 Bool_t AliSimulation::WriteRawData(const char* detectors,
939 const char* fileName,
940 Bool_t deleteIntermediateFiles)
942 // convert the digits to raw data
943 // First DDL raw data files for the given detectors are created.
944 // If a file name is given, the DDL files are then converted to a DATE file.
945 // If deleteIntermediateFiles is true, the DDL raw files are deleted
947 // If the file name has the extension ".root", the DATE file is converted
949 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
951 TStopwatch stopwatch;
954 if (!WriteRawFiles(detectors)) {
955 if (fStopOnError) return kFALSE;
958 TString dateFileName(fileName);
959 if (!dateFileName.IsNull()) {
960 Bool_t rootOutput = dateFileName.EndsWith(".root");
961 if (rootOutput) dateFileName += ".date";
962 if (!ConvertRawFilesToDate(dateFileName)) {
963 if (fStopOnError) return kFALSE;
965 if (deleteIntermediateFiles) {
966 AliRunLoader* runLoader = LoadRun("READ");
967 if (runLoader) for (Int_t iEvent = 0;
968 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
970 sprintf(command, "rm -r raw%d", iEvent);
971 gSystem->Exec(command);
976 if (!ConvertDateToRoot(dateFileName, fileName)) {
977 if (fStopOnError) return kFALSE;
979 if (deleteIntermediateFiles) {
980 gSystem->Unlink(dateFileName);
985 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
986 stopwatch.RealTime(),stopwatch.CpuTime()));
991 //_____________________________________________________________________________
992 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
994 // convert the digits to raw data DDL files
996 AliRunLoader* runLoader = LoadRun("READ");
997 if (!runLoader) return kFALSE;
999 // write raw data to DDL files
1000 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1001 AliInfo(Form("processing event %d", iEvent));
1002 runLoader->GetEvent(iEvent);
1003 TString baseDir = gSystem->WorkingDirectory();
1005 sprintf(dirName, "raw%d", iEvent);
1006 gSystem->MakeDirectory(dirName);
1007 if (!gSystem->ChangeDirectory(dirName)) {
1008 AliError(Form("couldn't change to directory %s", dirName));
1009 if (fStopOnError) return kFALSE; else continue;
1012 TString detStr = detectors;
1013 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1014 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1015 AliModule* det = (AliModule*) detArray->At(iDet);
1016 if (!det || !det->IsActive()) continue;
1017 if (IsSelected(det->GetName(), detStr)) {
1018 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1023 gSystem->ChangeDirectory(baseDir);
1024 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1025 AliError(Form("the following detectors were not found: %s",
1027 if (fStopOnError) return kFALSE;
1035 //_____________________________________________________________________________
1036 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1038 // convert raw data DDL files to a DATE file with the program "dateStream"
1040 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1042 AliError("the program dateStream was not found");
1043 if (fStopOnError) return kFALSE;
1048 AliRunLoader* runLoader = LoadRun("READ");
1049 if (!runLoader) return kFALSE;
1051 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1053 sprintf(command, "dateStream -o %s -# %d -C",
1054 dateFileName, runLoader->GetNumberOfEvents());
1055 FILE* pipe = gSystem->OpenPipe(command, "w");
1057 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1058 fprintf(pipe, "GDC\n");
1062 // loop over detectors and DDLs
1063 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1064 for (Int_t iDDL = 0; iDDL < kDetectorDDLs[iDet]; iDDL++) {
1066 Int_t ddlID = 0x100*iDet + iDDL;
1067 Int_t ldcID = Int_t(ldc + 0.0001);
1068 ldc += kDetectorLDCs[iDet] / kDetectorDDLs[iDet];
1070 char rawFileName[256];
1071 sprintf(rawFileName, "raw%d/%s_%d.ddl",
1072 iEvent, kDetectors[iDet], ddlID);
1074 // check existence and size of raw data file
1075 FILE* file = fopen(rawFileName, "rb");
1076 if (!file) continue;
1077 fseek(file, 0, SEEK_END);
1078 unsigned long size = ftell(file);
1080 if (!size) continue;
1082 if (ldcID != prevLDC) {
1083 fprintf(pipe, " LDC Id %d\n", ldcID);
1086 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1091 Int_t result = gSystem->ClosePipe(pipe);
1094 return (result == 0);
1097 //_____________________________________________________________________________
1098 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1099 const char* rootFileName)
1101 // convert a DATE file to a root file with the program "alimdc"
1104 const Int_t kDBSize = 1000000000;
1105 const Int_t kTagDBSize = 1000000000;
1106 const Bool_t kFilter = kFALSE;
1107 const Int_t kCompression = 1;
1109 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1111 AliError("the program alimdc was not found");
1112 if (fStopOnError) return kFALSE;
1117 AliInfo(Form("converting DATE file %s to root file %s",
1118 dateFileName, rootFileName));
1120 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1121 const char* tagDBFS = "/tmp/mdc1/tags";
1122 const char* runDBFS = "/tmp/mdc1/meta";
1124 // User defined file system locations
1125 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1126 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1127 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1128 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1129 if (gSystem->Getenv("ALIMDC_TAGDB"))
1130 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1131 if (gSystem->Getenv("ALIMDC_RUNDB"))
1132 runDBFS = gSystem->Getenv("ALIMDC_RUNDB");
1134 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1135 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1136 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1137 gSystem->Exec(Form("rm -rf %s",runDBFS));
1139 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1140 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1141 gSystem->Exec(Form("mkdir %s",tagDBFS));
1142 gSystem->Exec(Form("mkdir %s",runDBFS));
1144 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1145 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1146 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1148 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1149 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1150 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1151 gSystem->Exec(Form("rm -rf %s",runDBFS));
1153 return (result == 0);
1157 //_____________________________________________________________________________
1158 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1160 // delete existing run loaders, open a new one and load gAlice
1162 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1163 AliRunLoader* runLoader =
1164 AliRunLoader::Open(fGAliceFileName.Data(),
1165 AliConfig::GetDefaultEventFolderName(), mode);
1167 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1170 runLoader->LoadgAlice();
1171 gAlice = runLoader->GetAliRun();
1173 AliError(Form("no gAlice object found in file %s",
1174 fGAliceFileName.Data()));
1180 //_____________________________________________________________________________
1181 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1183 // get or calculate the number of signal events per background event
1185 if (!fBkgrdFileNames) return 1;
1186 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1187 if (nBkgrdFiles == 0) return 1;
1189 // get the number of signal events
1191 AliRunLoader* runLoader =
1192 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1193 if (!runLoader) return 1;
1194 nEvents = runLoader->GetNumberOfEvents();
1199 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1200 // get the number of background events
1201 const char* fileName = ((TObjString*)
1202 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1203 AliRunLoader* runLoader =
1204 AliRunLoader::Open(fileName, "BKGRD");
1205 if (!runLoader) continue;
1206 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1209 // get or calculate the number of signal per background events
1210 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1211 if (nSignalPerBkgrd <= 0) {
1212 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1213 } else if (result && (result != nSignalPerBkgrd)) {
1214 AliInfo(Form("the number of signal events per background event "
1215 "will be changed from %d to %d for stream %d",
1216 nSignalPerBkgrd, result, iBkgrdFile+1));
1217 nSignalPerBkgrd = result;
1220 if (!result) result = nSignalPerBkgrd;
1221 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1222 AliWarning(Form("not enough background events (%d) for %d signal events "
1223 "using %d signal per background events for stream %d",
1224 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1231 //_____________________________________________________________________________
1232 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1234 // check whether detName is contained in detectors
1235 // if yes, it is removed from detectors
1237 // check if all detectors are selected
1238 if ((detectors.CompareTo("ALL") == 0) ||
1239 detectors.BeginsWith("ALL ") ||
1240 detectors.EndsWith(" ALL") ||
1241 detectors.Contains(" ALL ")) {
1246 // search for the given detector
1247 Bool_t result = kFALSE;
1248 if ((detectors.CompareTo(detName) == 0) ||
1249 detectors.BeginsWith(detName+" ") ||
1250 detectors.EndsWith(" "+detName) ||
1251 detectors.Contains(" "+detName+" ")) {
1252 detectors.ReplaceAll(detName, "");
1256 // clean up the detectors string
1257 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1258 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1259 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);