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"
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"
130 #include "AliCTPRawData.h"
132 ClassImp(AliSimulation)
135 //_____________________________________________________________________________
136 AliSimulation::AliSimulation(const char* configFileName, const char* cdbUri,
137 const char* name, const char* title) :
140 fRunGeneration(kTRUE),
141 fRunSimulation(kTRUE),
142 fLoadAlignFromCDB(kTRUE),
143 fLoadAlignData("ALL"),
147 fMakeDigitsFromHits(""),
149 fRawDataFileName(""),
150 fDeleteIntermediateFiles(kFALSE),
151 fStopOnError(kFALSE),
154 fConfigFileName(configFileName),
155 fGAliceFileName("galice.root"),
157 fBkgrdFileNames(NULL),
158 fAlignObjArray(NULL),
159 fUseBkgrdVertex(kTRUE),
160 fRegionOfInterest(kFALSE),
164 // create simulation object with default parameters
166 SetGAliceFile("galice.root");
169 //_____________________________________________________________________________
170 AliSimulation::AliSimulation(const AliSimulation& sim) :
173 fRunGeneration(sim.fRunGeneration),
174 fRunSimulation(sim.fRunSimulation),
175 fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
176 fLoadAlignData(sim.fLoadAlignData),
177 fMakeSDigits(sim.fMakeSDigits),
178 fMakeDigits(sim.fMakeDigits),
179 fMakeTrigger(sim.fMakeTrigger),
180 fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
181 fWriteRawData(sim.fWriteRawData),
182 fRawDataFileName(""),
183 fDeleteIntermediateFiles(kFALSE),
184 fStopOnError(sim.fStopOnError),
186 fNEvents(sim.fNEvents),
187 fConfigFileName(sim.fConfigFileName),
188 fGAliceFileName(sim.fGAliceFileName),
190 fBkgrdFileNames(NULL),
191 fAlignObjArray(NULL),
192 fUseBkgrdVertex(sim.fUseBkgrdVertex),
193 fRegionOfInterest(sim.fRegionOfInterest),
194 fCDBUri(sim.fCDBUri),
199 for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
200 if (!sim.fEventsPerFile[i]) continue;
201 fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
204 fBkgrdFileNames = new TObjArray;
205 for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
206 if (!sim.fBkgrdFileNames->At(i)) continue;
207 fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
210 for (Int_t i = 0; i < sim.fSpecCDBUri.GetEntriesFast(); i++) {
211 if (sim.fSpecCDBUri[i]) fSpecCDBUri.Add(sim.fSpecCDBUri[i]->Clone());
216 //_____________________________________________________________________________
217 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
219 // assignment operator
221 this->~AliSimulation();
222 new(this) AliSimulation(sim);
226 //_____________________________________________________________________________
227 AliSimulation::~AliSimulation()
231 fEventsPerFile.Delete();
232 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
233 // delete fAlignObjArray; fAlignObjArray=0;
235 if (fBkgrdFileNames) {
236 fBkgrdFileNames->Delete();
237 delete fBkgrdFileNames;
240 fSpecCDBUri.Delete();
244 //_____________________________________________________________________________
245 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
247 // set the number of events for one run
252 //_____________________________________________________________________________
253 void AliSimulation::InitCDBStorage()
255 // activate a default CDB storage
256 // First check if we have any CDB storage set, because it is used
257 // to retrieve the calibration and alignment constants
259 AliCDBManager* man = AliCDBManager::Instance();
260 if (man->IsDefaultStorageSet())
262 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
263 AliWarning("Default CDB storage has been already set !");
264 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
265 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
269 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
270 AliWarning(Form("Default CDB storage is set to: %s",fCDBUri.Data()));
271 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
272 man->SetDefaultStorage(fCDBUri);
275 // Now activate the detector specific CDB storage locations
276 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
277 TObject* obj = fSpecCDBUri[i];
279 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
280 AliWarning(Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
281 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
282 man->SetSpecificStorage(obj->GetName(),obj->GetTitle());
287 //_____________________________________________________________________________
288 void AliSimulation::SetDefaultStorage(const char* uri) {
289 // Store the desired default CDB storage location
290 // Activate it later within the Run() method
296 //_____________________________________________________________________________
297 void AliSimulation::SetSpecificStorage(const char* detName, const char* uri) {
298 // Store a detector-specific CDB storage location
299 // Activate it later within the Run() method
301 TObject* obj = fSpecCDBUri.FindObject(detName);
302 if (obj) fSpecCDBUri.Remove(obj);
303 fSpecCDBUri.Add(new TNamed(detName, uri));
307 //_____________________________________________________________________________
308 void AliSimulation::SetConfigFile(const char* fileName)
310 // set the name of the config file
312 fConfigFileName = fileName;
315 //_____________________________________________________________________________
316 void AliSimulation::SetGAliceFile(const char* fileName)
318 // set the name of the galice file
319 // the path is converted to an absolute one if it is relative
321 fGAliceFileName = fileName;
322 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
323 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
325 fGAliceFileName = absFileName;
326 delete[] absFileName;
329 AliDebug(2, Form("galice file name set to %s", fileName));
332 //_____________________________________________________________________________
333 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
336 // set the number of events per file for the given detector and data type
337 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
339 TNamed* obj = new TNamed(detector, type);
340 obj->SetUniqueID(nEvents);
341 fEventsPerFile.Add(obj);
344 //_____________________________________________________________________________
345 Bool_t AliSimulation::ApplyAlignObjsToGeom(TObjArray* alObjArray)
347 // Read collection of alignment objects (AliAlignObj derived) saved
348 // in the TClonesArray ClArrayName and apply them to the geometry
349 // manager singleton.
352 Int_t nvols = alObjArray->GetEntriesFast();
356 for(Int_t j=0; j<nvols; j++)
358 AliAlignObj* alobj = (AliAlignObj*) alObjArray->UncheckedAt(j);
359 if (alobj->ApplyToGeometry() == kFALSE) flag = kFALSE;
362 if (AliDebugLevelClass() >= 1) {
363 gGeoManager->GetTopNode()->CheckOverlaps(20);
364 TObjArray* ovexlist = gGeoManager->GetListOfOverlaps();
365 if(ovexlist->GetEntriesFast()){
366 AliError("The application of alignment objects to the geometry caused huge overlaps/extrusions!");
374 //_____________________________________________________________________________
375 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* fileName, const char* clArrayName)
377 // read collection of alignment objects (AliAlignObj derived) saved
378 // in the TClonesArray ClArrayName in the file fileName and apply
379 // them to the TGeo geometry passed as argument
382 TFile* inFile = TFile::Open(fileName,"READ");
383 if (!inFile || !inFile->IsOpen()) {
384 AliErrorClass(Form("Could not open file %s !",fileName));
388 TClonesArray* alObjArray = ((TClonesArray*) inFile->Get(clArrayName));
391 AliErrorClass(Form("Could not get array (%s) from file (%s) !",clArrayName,fileName));
395 return ApplyAlignObjsToGeom(alObjArray);
399 //_____________________________________________________________________________
400 Bool_t AliSimulation::ApplyAlignObjsToGeom(AliCDBParam* param, AliCDBId& Id)
402 // read collection of alignment objects (AliAlignObj derived) saved
403 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
404 // param (to get the AliCDBStorage) and Id; apply the alignment objects
405 // to the TGeo geometry passed as argument
408 AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage(param);
409 AliCDBEntry* entry = storage->Get(Id);
410 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
412 return ApplyAlignObjsToGeom(AlObjArray);
416 //_____________________________________________________________________________
417 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* uri, const char* path, Int_t runnum, Int_t version, Int_t sversion)
419 // read collection of alignment objects (AliAlignObj derived) saved
420 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
421 // param (to get the AliCDBStorage) and Id; apply the alignment objects
422 // to the TGeo geometry passed as argument
425 AliCDBParam* param = AliCDBManager::Instance()->CreateParameter(uri);
426 AliCDBId id(path, runnum, runnum, version, sversion);
428 return ApplyAlignObjsToGeom(param, id);
432 //_____________________________________________________________________________
433 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* detName, Int_t runnum, Int_t version, Int_t sversion)
435 // read collection of alignment objects (AliAlignObj derived) saved
436 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
437 // param (to get the AliCDBStorage) and Id; apply the alignment objects
438 // to the TGeo geometry passed as argument
441 AliCDBPath path(detName,"Align","Data");
442 AliCDBEntry* entry = AliCDBManager::Instance()->Get(path.GetPath(),runnum,version,sversion);
444 if(!entry) return kFALSE;
445 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
447 return ApplyAlignObjsToGeom(AlObjArray);
450 //_____________________________________________________________________________
451 Bool_t AliSimulation::SetAlignObjArraySingleDet(const char* detName)
453 // Fills array of single detector's alignable objects from CDB
455 AliDebug(2, Form("Loading alignment data for detector: %s",detName));
459 AliCDBPath path(detName,"Align","Data");
461 entry=AliCDBManager::Instance()->Get(path.GetPath());
463 AliDebug(2,Form("Couldn't load alignment data for detector %s",detName));
467 TClonesArray *alignArray = (TClonesArray*) entry->GetObject();
468 alignArray->SetOwner(0);
469 AliDebug(2,Form("Found %d alignment objects for %s",
470 alignArray->GetEntries(),detName));
472 AliAlignObj *alignObj=0;
473 TIter iter(alignArray);
475 // loop over align objects in detector
476 while( ( alignObj=(AliAlignObj *) iter.Next() ) ){
477 fAlignObjArray->Add(alignObj);
479 // delete entry --- Don't delete, it is cached!
481 AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
486 //_____________________________________________________________________________
487 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
489 // Read the alignment objects from CDB.
490 // Each detector is supposed to have the
491 // alignment objects in DET/Align/Data CDB path.
492 // All the detector objects are then collected,
493 // sorted by geometry level (starting from ALIC) and
494 // then applied to the TGeo geometry.
495 // Finally an overlaps check is performed.
497 Bool_t delRunLoader = kFALSE;
499 runLoader = LoadRun("READ");
500 if (!runLoader) return kFALSE;
501 delRunLoader = kTRUE;
504 // Load alignment data from CDB and fill fAlignObjArray
505 if(fLoadAlignFromCDB){
506 if(!fAlignObjArray) fAlignObjArray = new TObjArray();
508 //fAlignObjArray->RemoveAll();
509 fAlignObjArray->Clear();
510 fAlignObjArray->SetOwner(0);
512 TString detStr = fLoadAlignData;
513 TString dataNotLoaded="";
514 TString dataLoaded="";
516 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
517 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
518 AliModule* det = (AliModule*) detArray->At(iDet);
519 if (!det || !det->IsActive()) continue;
520 if (IsSelected(det->GetName(), detStr)) {
521 if(!SetAlignObjArraySingleDet(det->GetName())){
522 dataNotLoaded += det->GetName();
523 dataNotLoaded += " ";
525 dataLoaded += det->GetName();
529 } // end loop over detectors
531 if ((detStr.CompareTo("ALL") == 0)) detStr = "";
532 dataNotLoaded += detStr;
533 AliInfo(Form("Alignment data loaded for: %s",
535 AliInfo(Form("Didn't/couldn't load alignment data for: %s",
536 dataNotLoaded.Data()));
537 } // fLoadAlignFromCDB flag
539 // Check if the array with alignment objects was
540 // provided by the user. If yes, apply the objects
541 // to the present TGeo geometry
542 if (fAlignObjArray) {
543 if (gGeoManager && gGeoManager->IsClosed()) {
544 if (ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
545 AliError("The misalignment of one or more volumes failed!"
546 "Compare the list of simulated detectors and the list of detector alignment data!");
547 if (delRunLoader) delete runLoader;
552 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
553 if (delRunLoader) delete runLoader;
558 if (delRunLoader) delete runLoader;
564 //_____________________________________________________________________________
565 Bool_t AliSimulation::SetRunNumber()
567 // Set the CDB manager run number
568 // The run number is retrieved from gAlice
570 if(AliCDBManager::Instance()->GetRun() < 0) {
571 AliRunLoader* runLoader = LoadRun("READ");
572 if (!runLoader) return kFALSE;
574 AliCDBManager::Instance()->SetRun(runLoader->GetAliRun()->GetRunNumber());
575 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
582 //_____________________________________________________________________________
583 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
585 // add a file with background events for merging
587 TObjString* fileNameStr = new TObjString(fileName);
588 fileNameStr->SetUniqueID(nSignalPerBkgrd);
589 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
590 fBkgrdFileNames->Add(fileNameStr);
594 //_____________________________________________________________________________
595 Bool_t AliSimulation::Run(Int_t nEvents)
597 // run the generation, simulation and digitization
601 if (nEvents > 0) fNEvents = nEvents;
603 // generation and simulation -> hits
604 if (fRunGeneration) {
605 if (!RunSimulation()) if (fStopOnError) return kFALSE;
608 // Set run number in CDBManager (if it is not already set in RunSimulation)
609 if (!SetRunNumber()) if (fStopOnError) return kFALSE;
611 // Load and misalign the geometry
613 TGeoManager::Import("geometry.root");
614 if (!gGeoManager) if (fStopOnError) return kFALSE;
615 if (!MisalignGeometry()) if (fStopOnError) return kFALSE;
618 // hits -> summable digits
619 if (!fMakeSDigits.IsNull()) {
620 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
623 // summable digits -> digits
624 if (!fMakeDigits.IsNull()) {
625 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
626 if (fStopOnError) return kFALSE;
631 if (!fMakeDigitsFromHits.IsNull()) {
632 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
633 AliWarning(Form("Merging and direct creation of digits from hits "
634 "was selected for some detectors. "
635 "No merging will be done for the following detectors: %s",
636 fMakeDigitsFromHits.Data()));
638 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
639 if (fStopOnError) return kFALSE;
644 if (!RunTrigger(fMakeTrigger)) {
645 if (fStopOnError) return kFALSE;
648 // digits -> raw data
649 if (!fWriteRawData.IsNull()) {
650 if (!WriteRawData(fWriteRawData, fRawDataFileName,
651 fDeleteIntermediateFiles)) {
652 if (fStopOnError) return kFALSE;
659 //_____________________________________________________________________________
660 Bool_t AliSimulation::RunTrigger(const char* descriptors)
664 TStopwatch stopwatch;
667 AliRunLoader* runLoader = LoadRun("READ");
668 if (!runLoader) return kFALSE;
669 TString des = descriptors;
672 if (gAlice->GetTriggerDescriptor() != "") {
673 des = gAlice->GetTriggerDescriptor();
676 AliWarning("No trigger descriptor is specified. Skipping the trigger simulation...");
681 runLoader->MakeTree( "CT" );
682 AliCentralTrigger* aCTP = runLoader->GetTrigger();
684 aCTP->LoadDescriptor( des );
687 if( !aCTP->RunTrigger( runLoader ) ) {
694 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
695 stopwatch.RealTime(),stopwatch.CpuTime()));
702 //_____________________________________________________________________________
703 Bool_t AliSimulation::WriteTriggerRawData()
705 // Writes the CTP (trigger) DDL raw data
706 // Details of the format are given in the
707 // trigger TDR - pages 134 and 135.
708 AliCTPRawData writer;
714 //_____________________________________________________________________________
715 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
717 // run the generation and simulation
719 TStopwatch stopwatch;
723 AliError("no gAlice object. Restart aliroot and try again.");
726 if (gAlice->Modules()->GetEntries() > 0) {
727 AliError("gAlice was already run. Restart aliroot and try again.");
731 AliInfo(Form("initializing gAlice with config file %s",
732 fConfigFileName.Data()));
733 StdoutToAliInfo(StderrToAliError(
734 gAlice->Init(fConfigFileName.Data());
737 // Get the trigger descriptor string
738 // Either from AliSimulation or from
740 if (fMakeTrigger.IsNull()) {
741 if (gAlice->GetTriggerDescriptor() != "")
742 fMakeTrigger = gAlice->GetTriggerDescriptor();
745 gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
747 // Set run number in CDBManager
748 AliCDBManager::Instance()->SetRun(gAlice->GetRunNumber());
749 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
751 AliRunLoader* runLoader = gAlice->GetRunLoader();
753 AliError(Form("gAlice has no run loader object. "
754 "Check your config file: %s", fConfigFileName.Data()));
757 SetGAliceFile(runLoader->GetFileName());
759 // Export ideal geometry
760 if (gGeoManager) gGeoManager->Export("geometry.root");
763 // if (!MisalignGeometry(runLoader)) {
767 MisalignGeometry(runLoader);
769 // Temporary fix by A.Gheata
770 // Could be removed with the next Root version (>5.11)
772 TIter next(gGeoManager->GetListOfVolumes());
774 while ((vol = (TGeoVolume *)next())) {
775 if (vol->GetVoxels()) {
776 if (vol->GetVoxels()->NeedRebuild()) {
777 vol->GetVoxels()->Voxelize();
784 // Export (mis)aligned geometry
785 if (gGeoManager) gGeoManager->Export("misaligned_geometry.root");
787 // AliRunLoader* runLoader = gAlice->GetRunLoader();
789 // AliError(Form("gAlice has no run loader object. "
790 // "Check your config file: %s", fConfigFileName.Data()));
793 // SetGAliceFile(runLoader->GetFileName());
795 if (!gAlice->Generator()) {
796 AliError(Form("gAlice has no generator object. "
797 "Check your config file: %s", fConfigFileName.Data()));
800 if (nEvents <= 0) nEvents = fNEvents;
802 // get vertex from background file in case of merging
803 if (fUseBkgrdVertex &&
804 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
805 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
806 const char* fileName = ((TObjString*)
807 (fBkgrdFileNames->At(0)))->GetName();
808 AliInfo(Form("The vertex will be taken from the background "
809 "file %s with nSignalPerBackground = %d",
810 fileName, signalPerBkgrd));
811 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
812 gAlice->Generator()->SetVertexGenerator(vtxGen);
815 if (!fRunSimulation) {
816 gAlice->Generator()->SetTrackingFlag(0);
819 // set the number of events per file for given detectors and data types
820 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
821 if (!fEventsPerFile[i]) continue;
822 const char* detName = fEventsPerFile[i]->GetName();
823 const char* typeName = fEventsPerFile[i]->GetTitle();
824 TString loaderName(detName);
825 loaderName += "Loader";
826 AliLoader* loader = runLoader->GetLoader(loaderName);
828 AliError(Form("RunSimulation", "no loader for %s found\n"
829 "Number of events per file not set for %s %s",
830 detName, typeName, detName));
833 AliDataLoader* dataLoader =
834 loader->GetDataLoader(typeName);
836 AliError(Form("no data loader for %s found\n"
837 "Number of events per file not set for %s %s",
838 typeName, detName, typeName));
841 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
842 AliDebug(1, Form("number of events per file set to %d for %s %s",
843 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
846 AliInfo("running gAlice");
847 StdoutToAliInfo(StderrToAliError(
848 gAlice->Run(nEvents);
853 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
854 stopwatch.RealTime(),stopwatch.CpuTime()));
859 //_____________________________________________________________________________
860 Bool_t AliSimulation::RunSDigitization(const char* detectors)
862 // run the digitization and produce summable digits
864 TStopwatch stopwatch;
867 AliRunLoader* runLoader = LoadRun();
868 if (!runLoader) return kFALSE;
870 TString detStr = detectors;
871 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
872 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
873 AliModule* det = (AliModule*) detArray->At(iDet);
874 if (!det || !det->IsActive()) continue;
875 if (IsSelected(det->GetName(), detStr)) {
876 AliInfo(Form("creating summable digits for %s", det->GetName()));
877 TStopwatch stopwatchDet;
878 stopwatchDet.Start();
880 AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
881 det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
885 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
886 AliError(Form("the following detectors were not found: %s",
888 if (fStopOnError) return kFALSE;
893 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
894 stopwatch.RealTime(),stopwatch.CpuTime()));
900 //_____________________________________________________________________________
901 Bool_t AliSimulation::RunDigitization(const char* detectors,
902 const char* excludeDetectors)
904 // run the digitization and produce digits from sdigits
906 TStopwatch stopwatch;
909 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
910 if (gAlice) delete gAlice;
914 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
915 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
916 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
917 manager->SetInputStream(0, fGAliceFileName.Data());
918 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
919 const char* fileName = ((TObjString*)
920 (fBkgrdFileNames->At(iStream-1)))->GetName();
921 manager->SetInputStream(iStream, fileName);
924 TString detStr = detectors;
925 TString detExcl = excludeDetectors;
926 manager->GetInputStream(0)->ImportgAlice();
927 AliRunLoader* runLoader =
928 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
929 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
930 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
931 AliModule* det = (AliModule*) detArray->At(iDet);
932 if (!det || !det->IsActive()) continue;
933 if (IsSelected(det->GetName(), detStr) &&
934 !IsSelected(det->GetName(), detExcl)) {
935 AliDigitizer* digitizer = det->CreateDigitizer(manager);
937 AliError(Form("no digitizer for %s", det->GetName()));
938 if (fStopOnError) return kFALSE;
940 digitizer->SetRegionOfInterest(fRegionOfInterest);
945 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
946 AliError(Form("the following detectors were not found: %s",
948 if (fStopOnError) return kFALSE;
951 if (!manager->GetListOfTasks()->IsEmpty()) {
952 AliInfo("executing digitization");
958 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
959 stopwatch.RealTime(),stopwatch.CpuTime()));
964 //_____________________________________________________________________________
965 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
967 // run the digitization and produce digits from hits
969 TStopwatch stopwatch;
972 AliRunLoader* runLoader = LoadRun("READ");
973 if (!runLoader) return kFALSE;
975 TString detStr = detectors;
976 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
977 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
978 AliModule* det = (AliModule*) detArray->At(iDet);
979 if (!det || !det->IsActive()) continue;
980 if (IsSelected(det->GetName(), detStr)) {
981 AliInfo(Form("creating digits from hits for %s", det->GetName()));
986 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
987 AliError(Form("the following detectors were not found: %s",
989 if (fStopOnError) return kFALSE;
993 //PH Temporary fix to avoid interference with the PHOS loder/getter
994 //PH The problem has to be solved in more general way 09/06/05
996 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
997 stopwatch.RealTime(),stopwatch.CpuTime()));
1002 //_____________________________________________________________________________
1003 Bool_t AliSimulation::WriteRawData(const char* detectors,
1004 const char* fileName,
1005 Bool_t deleteIntermediateFiles)
1007 // convert the digits to raw data
1008 // First DDL raw data files for the given detectors are created.
1009 // If a file name is given, the DDL files are then converted to a DATE file.
1010 // If deleteIntermediateFiles is true, the DDL raw files are deleted
1012 // If the file name has the extension ".root", the DATE file is converted
1014 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1016 TStopwatch stopwatch;
1019 if (!WriteRawFiles(detectors)) {
1020 if (fStopOnError) return kFALSE;
1023 TString dateFileName(fileName);
1024 if (!dateFileName.IsNull()) {
1025 Bool_t rootOutput = dateFileName.EndsWith(".root");
1026 if (rootOutput) dateFileName += ".date";
1027 if (!ConvertRawFilesToDate(dateFileName)) {
1028 if (fStopOnError) return kFALSE;
1030 if (deleteIntermediateFiles) {
1031 AliRunLoader* runLoader = LoadRun("READ");
1032 if (runLoader) for (Int_t iEvent = 0;
1033 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1035 sprintf(command, "rm -r raw%d", iEvent);
1036 gSystem->Exec(command);
1041 if (!ConvertDateToRoot(dateFileName, fileName)) {
1042 if (fStopOnError) return kFALSE;
1044 if (deleteIntermediateFiles) {
1045 gSystem->Unlink(dateFileName);
1050 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1051 stopwatch.RealTime(),stopwatch.CpuTime()));
1056 //_____________________________________________________________________________
1057 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1059 // convert the digits to raw data DDL files
1061 AliRunLoader* runLoader = LoadRun("READ");
1062 if (!runLoader) return kFALSE;
1064 // write raw data to DDL files
1065 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1066 AliInfo(Form("processing event %d", iEvent));
1067 runLoader->GetEvent(iEvent);
1068 TString baseDir = gSystem->WorkingDirectory();
1070 sprintf(dirName, "raw%d", iEvent);
1071 gSystem->MakeDirectory(dirName);
1072 if (!gSystem->ChangeDirectory(dirName)) {
1073 AliError(Form("couldn't change to directory %s", dirName));
1074 if (fStopOnError) return kFALSE; else continue;
1077 TString detStr = detectors;
1078 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1079 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1080 AliModule* det = (AliModule*) detArray->At(iDet);
1081 if (!det || !det->IsActive()) continue;
1082 if (IsSelected(det->GetName(), detStr)) {
1083 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1088 if (!WriteTriggerRawData())
1089 if (fStopOnError) return kFALSE;
1091 gSystem->ChangeDirectory(baseDir);
1092 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1093 AliError(Form("the following detectors were not found: %s",
1095 if (fStopOnError) return kFALSE;
1103 //_____________________________________________________________________________
1104 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1106 // convert raw data DDL files to a DATE file with the program "dateStream"
1108 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1110 AliError("the program dateStream was not found");
1111 if (fStopOnError) return kFALSE;
1116 AliRunLoader* runLoader = LoadRun("READ");
1117 if (!runLoader) return kFALSE;
1119 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1121 sprintf(command, "dateStream -D -o %s -# %d -C",
1122 dateFileName, runLoader->GetNumberOfEvents());
1123 FILE* pipe = gSystem->OpenPipe(command, "w");
1125 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1126 fprintf(pipe, "GDC\n");
1130 // loop over detectors and DDLs
1131 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1132 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1134 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1135 Int_t ldcID = Int_t(ldc + 0.0001);
1136 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1138 char rawFileName[256];
1139 sprintf(rawFileName, "raw%d/%s",
1140 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1142 // check existence and size of raw data file
1143 FILE* file = fopen(rawFileName, "rb");
1144 if (!file) continue;
1145 fseek(file, 0, SEEK_END);
1146 unsigned long size = ftell(file);
1148 if (!size) continue;
1150 if (ldcID != prevLDC) {
1151 fprintf(pipe, " LDC Id %d\n", ldcID);
1154 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1159 Int_t result = gSystem->ClosePipe(pipe);
1162 return (result == 0);
1165 //_____________________________________________________________________________
1166 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1167 const char* rootFileName)
1169 // convert a DATE file to a root file with the program "alimdc"
1172 const Int_t kDBSize = 1000000000;
1173 const Int_t kTagDBSize = 1000000000;
1174 const Bool_t kFilter = kFALSE;
1175 const Int_t kCompression = 1;
1177 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1179 AliError("the program alimdc was not found");
1180 if (fStopOnError) return kFALSE;
1185 AliInfo(Form("converting DATE file %s to root file %s",
1186 dateFileName, rootFileName));
1188 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1189 const char* tagDBFS = "/tmp/mdc1/tags";
1190 const char* runDBFS = "/tmp/mdc1/meta";
1192 // User defined file system locations
1193 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1194 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1195 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1196 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1197 if (gSystem->Getenv("ALIMDC_TAGDB"))
1198 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1199 if (gSystem->Getenv("ALIMDC_RUNDB"))
1200 runDBFS = gSystem->Getenv("ALIMDC_RUNDB");
1202 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1203 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1204 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1205 gSystem->Exec(Form("rm -rf %s",runDBFS));
1207 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1208 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1209 gSystem->Exec(Form("mkdir %s",tagDBFS));
1210 gSystem->Exec(Form("mkdir %s",runDBFS));
1212 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1213 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1214 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1216 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1217 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1218 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1219 gSystem->Exec(Form("rm -rf %s",runDBFS));
1221 return (result == 0);
1225 //_____________________________________________________________________________
1226 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1228 // delete existing run loaders, open a new one and load gAlice
1230 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1231 AliRunLoader* runLoader =
1232 AliRunLoader::Open(fGAliceFileName.Data(),
1233 AliConfig::GetDefaultEventFolderName(), mode);
1235 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1238 runLoader->LoadgAlice();
1239 gAlice = runLoader->GetAliRun();
1241 AliError(Form("no gAlice object found in file %s",
1242 fGAliceFileName.Data()));
1248 //_____________________________________________________________________________
1249 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1251 // get or calculate the number of signal events per background event
1253 if (!fBkgrdFileNames) return 1;
1254 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1255 if (nBkgrdFiles == 0) return 1;
1257 // get the number of signal events
1259 AliRunLoader* runLoader =
1260 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1261 if (!runLoader) return 1;
1262 nEvents = runLoader->GetNumberOfEvents();
1267 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1268 // get the number of background events
1269 const char* fileName = ((TObjString*)
1270 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1271 AliRunLoader* runLoader =
1272 AliRunLoader::Open(fileName, "BKGRD");
1273 if (!runLoader) continue;
1274 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1277 // get or calculate the number of signal per background events
1278 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1279 if (nSignalPerBkgrd <= 0) {
1280 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1281 } else if (result && (result != nSignalPerBkgrd)) {
1282 AliInfo(Form("the number of signal events per background event "
1283 "will be changed from %d to %d for stream %d",
1284 nSignalPerBkgrd, result, iBkgrdFile+1));
1285 nSignalPerBkgrd = result;
1288 if (!result) result = nSignalPerBkgrd;
1289 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1290 AliWarning(Form("not enough background events (%d) for %d signal events "
1291 "using %d signal per background events for stream %d",
1292 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1299 //_____________________________________________________________________________
1300 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1302 // check whether detName is contained in detectors
1303 // if yes, it is removed from detectors
1305 // check if all detectors are selected
1306 if ((detectors.CompareTo("ALL") == 0) ||
1307 detectors.BeginsWith("ALL ") ||
1308 detectors.EndsWith(" ALL") ||
1309 detectors.Contains(" ALL ")) {
1314 // search for the given detector
1315 Bool_t result = kFALSE;
1316 if ((detectors.CompareTo(detName) == 0) ||
1317 detectors.BeginsWith(detName+" ") ||
1318 detectors.EndsWith(" "+detName) ||
1319 detectors.Contains(" "+detName+" ")) {
1320 detectors.ReplaceAll(detName, "");
1324 // clean up the detectors string
1325 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1326 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1327 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);