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"
131 #include "AliRawReaderFile.h"
133 #include "AliHeader.h"
134 #include "AliGenEventHeader.h"
136 ClassImp(AliSimulation)
139 //_____________________________________________________________________________
140 AliSimulation::AliSimulation(const char* configFileName, const char* cdbUri,
141 const char* name, const char* title) :
144 fRunGeneration(kTRUE),
145 fRunSimulation(kTRUE),
146 fLoadAlignFromCDB(kTRUE),
147 fLoadAlignData("ALL"),
151 fMakeDigitsFromHits(""),
153 fRawDataFileName(""),
154 fDeleteIntermediateFiles(kFALSE),
155 fStopOnError(kFALSE),
158 fConfigFileName(configFileName),
159 fGAliceFileName("galice.root"),
161 fBkgrdFileNames(NULL),
162 fAlignObjArray(NULL),
163 fUseBkgrdVertex(kTRUE),
164 fRegionOfInterest(kFALSE),
167 fEmbeddingFlag(kFALSE)
169 // create simulation object with default parameters
171 SetGAliceFile("galice.root");
174 //_____________________________________________________________________________
175 AliSimulation::AliSimulation(const AliSimulation& sim) :
178 fRunGeneration(sim.fRunGeneration),
179 fRunSimulation(sim.fRunSimulation),
180 fLoadAlignFromCDB(sim.fLoadAlignFromCDB),
181 fLoadAlignData(sim.fLoadAlignData),
182 fMakeSDigits(sim.fMakeSDigits),
183 fMakeDigits(sim.fMakeDigits),
184 fMakeTrigger(sim.fMakeTrigger),
185 fMakeDigitsFromHits(sim.fMakeDigitsFromHits),
186 fWriteRawData(sim.fWriteRawData),
187 fRawDataFileName(""),
188 fDeleteIntermediateFiles(kFALSE),
189 fStopOnError(sim.fStopOnError),
191 fNEvents(sim.fNEvents),
192 fConfigFileName(sim.fConfigFileName),
193 fGAliceFileName(sim.fGAliceFileName),
195 fBkgrdFileNames(NULL),
196 fAlignObjArray(NULL),
197 fUseBkgrdVertex(sim.fUseBkgrdVertex),
198 fRegionOfInterest(sim.fRegionOfInterest),
199 fCDBUri(sim.fCDBUri),
201 fEmbeddingFlag(sim.fEmbeddingFlag)
205 for (Int_t i = 0; i < sim.fEventsPerFile.GetEntriesFast(); i++) {
206 if (!sim.fEventsPerFile[i]) continue;
207 fEventsPerFile.Add(sim.fEventsPerFile[i]->Clone());
210 fBkgrdFileNames = new TObjArray;
211 for (Int_t i = 0; i < sim.fBkgrdFileNames->GetEntriesFast(); i++) {
212 if (!sim.fBkgrdFileNames->At(i)) continue;
213 fBkgrdFileNames->Add(sim.fBkgrdFileNames->At(i)->Clone());
216 for (Int_t i = 0; i < sim.fSpecCDBUri.GetEntriesFast(); i++) {
217 if (sim.fSpecCDBUri[i]) fSpecCDBUri.Add(sim.fSpecCDBUri[i]->Clone());
222 //_____________________________________________________________________________
223 AliSimulation& AliSimulation::operator = (const AliSimulation& sim)
225 // assignment operator
227 this->~AliSimulation();
228 new(this) AliSimulation(sim);
232 //_____________________________________________________________________________
233 AliSimulation::~AliSimulation()
237 fEventsPerFile.Delete();
238 // if(fAlignObjArray) fAlignObjArray->Delete(); // fAlignObjArray->RemoveAll() ???
239 // delete fAlignObjArray; fAlignObjArray=0;
241 if (fBkgrdFileNames) {
242 fBkgrdFileNames->Delete();
243 delete fBkgrdFileNames;
246 fSpecCDBUri.Delete();
250 //_____________________________________________________________________________
251 void AliSimulation::SetNumberOfEvents(Int_t nEvents)
253 // set the number of events for one run
258 //_____________________________________________________________________________
259 void AliSimulation::InitCDBStorage()
261 // activate a default CDB storage
262 // First check if we have any CDB storage set, because it is used
263 // to retrieve the calibration and alignment constants
265 AliCDBManager* man = AliCDBManager::Instance();
266 if (man->IsDefaultStorageSet())
268 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
269 AliWarning("Default CDB storage has been already set !");
270 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
271 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
275 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
276 AliDebug(2, Form("Default CDB storage is set to: %s",fCDBUri.Data()));
277 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
278 man->SetDefaultStorage(fCDBUri);
281 // Now activate the detector specific CDB storage locations
282 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
283 TObject* obj = fSpecCDBUri[i];
285 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
286 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
287 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
288 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
293 //_____________________________________________________________________________
294 void AliSimulation::SetDefaultStorage(const char* uri) {
295 // Store the desired default CDB storage location
296 // Activate it later within the Run() method
302 //_____________________________________________________________________________
303 void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
304 // Store a detector-specific CDB storage location
305 // Activate it later within the Run() method
307 AliCDBPath aPath(calibType);
308 if(!aPath.IsValid()){
309 AliError(Form("Not a valid path: %s", calibType));
313 TObject* obj = fSpecCDBUri.FindObject(calibType);
314 if (obj) fSpecCDBUri.Remove(obj);
315 fSpecCDBUri.Add(new TNamed(calibType, uri));
319 //_____________________________________________________________________________
320 void AliSimulation::SetConfigFile(const char* fileName)
322 // set the name of the config file
324 fConfigFileName = fileName;
327 //_____________________________________________________________________________
328 void AliSimulation::SetGAliceFile(const char* fileName)
330 // set the name of the galice file
331 // the path is converted to an absolute one if it is relative
333 fGAliceFileName = fileName;
334 if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
335 char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
337 fGAliceFileName = absFileName;
338 delete[] absFileName;
341 AliDebug(2, Form("galice file name set to %s", fileName));
344 //_____________________________________________________________________________
345 void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
348 // set the number of events per file for the given detector and data type
349 // ("Hits", "Summable Digits", "Digits", "Reconstructed Points" or "Tracks")
351 TNamed* obj = new TNamed(detector, type);
352 obj->SetUniqueID(nEvents);
353 fEventsPerFile.Add(obj);
356 //_____________________________________________________________________________
357 Bool_t AliSimulation::ApplyAlignObjsToGeom(TObjArray* alObjArray)
359 // Read collection of alignment objects (AliAlignObj derived) saved
360 // in the TClonesArray ClArrayName and apply them to the geometry
361 // manager singleton.
364 Int_t nvols = alObjArray->GetEntriesFast();
368 for(Int_t j=0; j<nvols; j++)
370 AliAlignObj* alobj = (AliAlignObj*) alObjArray->UncheckedAt(j);
371 if (alobj->ApplyToGeometry() == kFALSE) flag = kFALSE;
374 if (AliDebugLevelClass() >= 1) {
375 gGeoManager->GetTopNode()->CheckOverlaps(1);
376 TObjArray* ovexlist = gGeoManager->GetListOfOverlaps();
377 if(ovexlist->GetEntriesFast()){
378 AliError("The application of alignment objects to the geometry caused huge overlaps/extrusions!");
386 //_____________________________________________________________________________
387 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* fileName, const char* clArrayName)
389 // read collection of alignment objects (AliAlignObj derived) saved
390 // in the TClonesArray ClArrayName in the file fileName and apply
391 // them to the TGeo geometry passed as argument
394 TFile* inFile = TFile::Open(fileName,"READ");
395 if (!inFile || !inFile->IsOpen()) {
396 AliErrorClass(Form("Could not open file %s !",fileName));
400 TClonesArray* alObjArray = ((TClonesArray*) inFile->Get(clArrayName));
403 AliErrorClass(Form("Could not get array (%s) from file (%s) !",clArrayName,fileName));
407 return ApplyAlignObjsToGeom(alObjArray);
411 //_____________________________________________________________________________
412 Bool_t AliSimulation::ApplyAlignObjsToGeom(AliCDBParam* param, AliCDBId& Id)
414 // read collection of alignment objects (AliAlignObj derived) saved
415 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
416 // param (to get the AliCDBStorage) and Id; apply the alignment objects
417 // to the TGeo geometry passed as argument
420 AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage(param);
421 AliCDBEntry* entry = storage->Get(Id);
422 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
424 return ApplyAlignObjsToGeom(AlObjArray);
428 //_____________________________________________________________________________
429 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* uri, const char* path, Int_t runnum, Int_t version, Int_t sversion)
431 // read collection of alignment objects (AliAlignObj derived) saved
432 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
433 // param (to get the AliCDBStorage) and Id; apply the alignment objects
434 // to the TGeo geometry passed as argument
437 AliCDBParam* param = AliCDBManager::Instance()->CreateParameter(uri);
438 AliCDBId id(path, runnum, runnum, version, sversion);
440 return ApplyAlignObjsToGeom(param, id);
444 //_____________________________________________________________________________
445 Bool_t AliSimulation::ApplyAlignObjsToGeom(const char* detName, Int_t runnum, Int_t version, Int_t sversion)
447 // read collection of alignment objects (AliAlignObj derived) saved
448 // in the TClonesArray ClArrayName in the AliCDBEntry identified by
449 // param (to get the AliCDBStorage) and Id; apply the alignment objects
450 // to the TGeo geometry passed as argument
453 AliCDBPath path(detName,"Align","Data");
454 AliCDBEntry* entry = AliCDBManager::Instance()->Get(path.GetPath(),runnum,version,sversion);
456 if(!entry) return kFALSE;
457 TClonesArray* AlObjArray = ((TClonesArray*) entry->GetObject());
459 return ApplyAlignObjsToGeom(AlObjArray);
462 //_____________________________________________________________________________
463 Bool_t AliSimulation::SetAlignObjArraySingleDet(const char* detName)
465 // Fills array of single detector's alignable objects from CDB
467 AliDebug(2, Form("Loading alignment data for detector: %s",detName));
471 AliCDBPath path(detName,"Align","Data");
473 entry=AliCDBManager::Instance()->Get(path.GetPath());
475 AliDebug(2,Form("Couldn't load alignment data for detector %s",detName));
479 TClonesArray *alignArray = (TClonesArray*) entry->GetObject();
480 alignArray->SetOwner(0);
481 AliDebug(2,Form("Found %d alignment objects for %s",
482 alignArray->GetEntries(),detName));
484 AliAlignObj *alignObj=0;
485 TIter iter(alignArray);
487 // loop over align objects in detector
488 while( ( alignObj=(AliAlignObj *) iter.Next() ) ){
489 fAlignObjArray->Add(alignObj);
491 // delete entry --- Don't delete, it is cached!
493 AliDebug(2, Form("fAlignObjArray entries: %d",fAlignObjArray->GetEntries() ));
498 //_____________________________________________________________________________
499 Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
501 // Read the alignment objects from CDB.
502 // Each detector is supposed to have the
503 // alignment objects in DET/Align/Data CDB path.
504 // All the detector objects are then collected,
505 // sorted by geometry level (starting from ALIC) and
506 // then applied to the TGeo geometry.
507 // Finally an overlaps check is performed.
509 Bool_t delRunLoader = kFALSE;
511 runLoader = LoadRun("READ");
512 if (!runLoader) return kFALSE;
513 delRunLoader = kTRUE;
516 // Load alignment data from CDB and fill fAlignObjArray
517 if(fLoadAlignFromCDB){
518 if(!fAlignObjArray) fAlignObjArray = new TObjArray();
520 //fAlignObjArray->RemoveAll();
521 fAlignObjArray->Clear();
522 fAlignObjArray->SetOwner(0);
524 TString detStr = fLoadAlignData;
525 TString dataNotLoaded="";
526 TString dataLoaded="";
528 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
529 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
530 AliModule* det = (AliModule*) detArray->At(iDet);
531 if (!det || !det->IsActive()) continue;
532 if (IsSelected(det->GetName(), detStr)) {
533 if(!SetAlignObjArraySingleDet(det->GetName())){
534 dataNotLoaded += det->GetName();
535 dataNotLoaded += " ";
537 dataLoaded += det->GetName();
541 } // end loop over detectors
543 if ((detStr.CompareTo("ALL") == 0)) detStr = "";
544 dataNotLoaded += detStr;
545 if(!dataLoaded.IsNull()) AliInfo(Form("Alignment data loaded for: %s",
547 if(!dataNotLoaded.IsNull()) AliInfo(Form("Didn't/couldn't load alignment data for: %s",
548 dataNotLoaded.Data()));
549 } // fLoadAlignFromCDB flag
551 // Check if the array with alignment objects was
552 // provided by the user. If yes, apply the objects
553 // to the present TGeo geometry
554 if (fAlignObjArray) {
555 if (gGeoManager && gGeoManager->IsClosed()) {
556 if (ApplyAlignObjsToGeom(fAlignObjArray) == kFALSE) {
557 AliError("The misalignment of one or more volumes failed!"
558 "Compare the list of simulated detectors and the list of detector alignment data!");
559 if (delRunLoader) delete runLoader;
564 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
565 if (delRunLoader) delete runLoader;
571 // Update the internal geometry of modules (ITS needs it)
572 TString detStr = fLoadAlignData;
573 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
574 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
576 AliModule* det = (AliModule*) detArray->At(iDet);
577 if (!det || !det->IsActive()) continue;
578 if (IsSelected(det->GetName(), detStr)) {
579 det->UpdateInternalGeometry();
581 } // end loop over detectors
584 if (delRunLoader) delete runLoader;
590 //_____________________________________________________________________________
591 Bool_t AliSimulation::SetRunNumber()
593 // Set the CDB manager run number
594 // The run number is retrieved from gAlice
596 if(AliCDBManager::Instance()->GetRun() < 0) {
597 AliRunLoader* runLoader = LoadRun("READ");
598 if (!runLoader) return kFALSE;
600 AliCDBManager::Instance()->SetRun(runLoader->GetAliRun()->GetRunNumber());
601 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
608 //_____________________________________________________________________________
609 void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
611 // add a file with background events for merging
613 TObjString* fileNameStr = new TObjString(fileName);
614 fileNameStr->SetUniqueID(nSignalPerBkgrd);
615 if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
616 fBkgrdFileNames->Add(fileNameStr);
619 void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
621 // add a file with background events for embeddin
622 MergeWith(fileName, nSignalPerBkgrd);
623 fEmbeddingFlag = kTRUE;
626 //_____________________________________________________________________________
627 Bool_t AliSimulation::Run(Int_t nEvents)
629 // run the generation, simulation and digitization
633 if (nEvents > 0) fNEvents = nEvents;
635 // generation and simulation -> hits
636 if (fRunGeneration) {
637 if (!RunSimulation()) if (fStopOnError) return kFALSE;
640 // Set run number in CDBManager (if it is not already set in RunSimulation)
641 if (!SetRunNumber()) if (fStopOnError) return kFALSE;
643 // Load and misalign the geometry
645 TGeoManager::Import("geometry.root");
646 if (!gGeoManager) if (fStopOnError) return kFALSE;
647 if (!MisalignGeometry()) if (fStopOnError) return kFALSE;
650 // hits -> summable digits
651 if (!fMakeSDigits.IsNull()) {
652 if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
655 // summable digits -> digits
656 if (!fMakeDigits.IsNull()) {
657 if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
658 if (fStopOnError) return kFALSE;
663 if (!fMakeDigitsFromHits.IsNull()) {
664 if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
665 AliWarning(Form("Merging and direct creation of digits from hits "
666 "was selected for some detectors. "
667 "No merging will be done for the following detectors: %s",
668 fMakeDigitsFromHits.Data()));
670 if (!RunHitsDigitization(fMakeDigitsFromHits)) {
671 if (fStopOnError) return kFALSE;
676 if (!RunTrigger(fMakeTrigger)) {
677 if (fStopOnError) return kFALSE;
680 // digits -> raw data
681 if (!fWriteRawData.IsNull()) {
682 if (!WriteRawData(fWriteRawData, fRawDataFileName,
683 fDeleteIntermediateFiles)) {
684 if (fStopOnError) return kFALSE;
691 //_____________________________________________________________________________
692 Bool_t AliSimulation::RunTrigger(const char* descriptors)
696 TStopwatch stopwatch;
699 AliRunLoader* runLoader = LoadRun("READ");
700 if (!runLoader) return kFALSE;
701 TString des = descriptors;
704 if (gAlice->GetTriggerDescriptor() != "") {
705 des = gAlice->GetTriggerDescriptor();
708 AliWarning("No trigger descriptor is specified. Skipping the trigger simulation...");
713 runLoader->MakeTree( "GG" );
714 AliCentralTrigger* aCTP = runLoader->GetTrigger();
716 aCTP->LoadDescriptor( des );
719 if( !aCTP->RunTrigger( runLoader ) ) {
726 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
727 stopwatch.RealTime(),stopwatch.CpuTime()));
734 //_____________________________________________________________________________
735 Bool_t AliSimulation::WriteTriggerRawData()
737 // Writes the CTP (trigger) DDL raw data
738 // Details of the format are given in the
739 // trigger TDR - pages 134 and 135.
740 AliCTPRawData writer;
746 //_____________________________________________________________________________
747 Bool_t AliSimulation::RunSimulation(Int_t nEvents)
749 // run the generation and simulation
751 TStopwatch stopwatch;
755 AliError("no gAlice object. Restart aliroot and try again.");
758 if (gAlice->Modules()->GetEntries() > 0) {
759 AliError("gAlice was already run. Restart aliroot and try again.");
763 AliInfo(Form("initializing gAlice with config file %s",
764 fConfigFileName.Data()));
765 StdoutToAliInfo(StderrToAliError(
766 gAlice->Init(fConfigFileName.Data());
769 // Get the trigger descriptor string
770 // Either from AliSimulation or from
772 if (fMakeTrigger.IsNull()) {
773 if (gAlice->GetTriggerDescriptor() != "")
774 fMakeTrigger = gAlice->GetTriggerDescriptor();
777 gAlice->SetTriggerDescriptor(fMakeTrigger.Data());
779 // Set run number in CDBManager
780 AliCDBManager::Instance()->SetRun(gAlice->GetRunNumber());
781 AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
783 AliRunLoader* runLoader = gAlice->GetRunLoader();
785 AliError(Form("gAlice has no run loader object. "
786 "Check your config file: %s", fConfigFileName.Data()));
789 SetGAliceFile(runLoader->GetFileName());
791 // Export ideal geometry
792 if (gGeoManager) gGeoManager->Export("geometry.root");
795 // if (!MisalignGeometry(runLoader)) {
799 MisalignGeometry(runLoader);
801 // Export (mis)aligned geometry
802 if (gGeoManager) gGeoManager->Export("misaligned_geometry.root");
804 // AliRunLoader* runLoader = gAlice->GetRunLoader();
806 // AliError(Form("gAlice has no run loader object. "
807 // "Check your config file: %s", fConfigFileName.Data()));
810 // SetGAliceFile(runLoader->GetFileName());
812 if (!gAlice->Generator()) {
813 AliError(Form("gAlice has no generator object. "
814 "Check your config file: %s", fConfigFileName.Data()));
817 if (nEvents <= 0) nEvents = fNEvents;
819 // get vertex from background file in case of merging
820 if (fUseBkgrdVertex &&
821 fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
822 Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
823 const char* fileName = ((TObjString*)
824 (fBkgrdFileNames->At(0)))->GetName();
825 AliInfo(Form("The vertex will be taken from the background "
826 "file %s with nSignalPerBackground = %d",
827 fileName, signalPerBkgrd));
828 AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
829 gAlice->Generator()->SetVertexGenerator(vtxGen);
832 if (!fRunSimulation) {
833 gAlice->Generator()->SetTrackingFlag(0);
836 // set the number of events per file for given detectors and data types
837 for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
838 if (!fEventsPerFile[i]) continue;
839 const char* detName = fEventsPerFile[i]->GetName();
840 const char* typeName = fEventsPerFile[i]->GetTitle();
841 TString loaderName(detName);
842 loaderName += "Loader";
843 AliLoader* loader = runLoader->GetLoader(loaderName);
845 AliError(Form("RunSimulation", "no loader for %s found\n"
846 "Number of events per file not set for %s %s",
847 detName, typeName, detName));
850 AliDataLoader* dataLoader =
851 loader->GetDataLoader(typeName);
853 AliError(Form("no data loader for %s found\n"
854 "Number of events per file not set for %s %s",
855 typeName, detName, typeName));
858 dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
859 AliDebug(1, Form("number of events per file set to %d for %s %s",
860 fEventsPerFile[i]->GetUniqueID(), detName, typeName));
863 AliInfo("running gAlice");
864 StdoutToAliInfo(StderrToAliError(
865 gAlice->Run(nEvents);
870 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
871 stopwatch.RealTime(),stopwatch.CpuTime()));
876 //_____________________________________________________________________________
877 Bool_t AliSimulation::RunSDigitization(const char* detectors)
879 // run the digitization and produce summable digits
881 TStopwatch stopwatch;
884 AliRunLoader* runLoader = LoadRun();
885 if (!runLoader) return kFALSE;
887 TString detStr = detectors;
888 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
889 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
890 AliModule* det = (AliModule*) detArray->At(iDet);
891 if (!det || !det->IsActive()) continue;
892 if (IsSelected(det->GetName(), detStr)) {
893 AliInfo(Form("creating summable digits for %s", det->GetName()));
894 TStopwatch stopwatchDet;
895 stopwatchDet.Start();
897 AliInfo(Form("Execution time for %s: R:%.2fs C:%.2fs",
898 det->GetName(),stopwatchDet.RealTime(),stopwatchDet.CpuTime()));
902 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
903 AliError(Form("the following detectors were not found: %s",
905 if (fStopOnError) return kFALSE;
910 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
911 stopwatch.RealTime(),stopwatch.CpuTime()));
917 //_____________________________________________________________________________
918 Bool_t AliSimulation::RunDigitization(const char* detectors,
919 const char* excludeDetectors)
921 // run the digitization and produce digits from sdigits
923 TStopwatch stopwatch;
926 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
927 if (gAlice) delete gAlice;
931 if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
932 Int_t signalPerBkgrd = GetNSignalPerBkgrd();
933 AliRunDigitizer* manager = new AliRunDigitizer(nStreams, signalPerBkgrd);
934 // manager->SetEmbeddingFlag(fEmbeddingFlag);
935 manager->SetInputStream(0, fGAliceFileName.Data());
936 for (Int_t iStream = 1; iStream < nStreams; iStream++) {
937 const char* fileName = ((TObjString*)
938 (fBkgrdFileNames->At(iStream-1)))->GetName();
939 manager->SetInputStream(iStream, fileName);
942 TString detStr = detectors;
943 TString detExcl = excludeDetectors;
944 manager->GetInputStream(0)->ImportgAlice();
945 AliRunLoader* runLoader =
946 AliRunLoader::GetRunLoader(manager->GetInputStream(0)->GetFolderName());
947 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
948 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
949 AliModule* det = (AliModule*) detArray->At(iDet);
950 if (!det || !det->IsActive()) continue;
951 if (IsSelected(det->GetName(), detStr) &&
952 !IsSelected(det->GetName(), detExcl)) {
953 AliDigitizer* digitizer = det->CreateDigitizer(manager);
956 AliError(Form("no digitizer for %s", det->GetName()));
957 if (fStopOnError) return kFALSE;
959 digitizer->SetRegionOfInterest(fRegionOfInterest);
964 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
965 AliError(Form("the following detectors were not found: %s",
967 if (fStopOnError) return kFALSE;
970 if (!manager->GetListOfTasks()->IsEmpty()) {
971 AliInfo("executing digitization");
977 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
978 stopwatch.RealTime(),stopwatch.CpuTime()));
983 //_____________________________________________________________________________
984 Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
986 // run the digitization and produce digits from hits
988 TStopwatch stopwatch;
991 AliRunLoader* runLoader = LoadRun("READ");
992 if (!runLoader) return kFALSE;
994 TString detStr = detectors;
995 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
996 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
997 AliModule* det = (AliModule*) detArray->At(iDet);
998 if (!det || !det->IsActive()) continue;
999 if (IsSelected(det->GetName(), detStr)) {
1000 AliInfo(Form("creating digits from hits for %s", det->GetName()));
1005 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1006 AliError(Form("the following detectors were not found: %s",
1008 if (fStopOnError) return kFALSE;
1012 //PH Temporary fix to avoid interference with the PHOS loder/getter
1013 //PH The problem has to be solved in more general way 09/06/05
1015 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1016 stopwatch.RealTime(),stopwatch.CpuTime()));
1021 //_____________________________________________________________________________
1022 Bool_t AliSimulation::WriteRawData(const char* detectors,
1023 const char* fileName,
1024 Bool_t deleteIntermediateFiles)
1026 // convert the digits to raw data
1027 // First DDL raw data files for the given detectors are created.
1028 // If a file name is given, the DDL files are then converted to a DATE file.
1029 // If deleteIntermediateFiles is true, the DDL raw files are deleted
1031 // If the file name has the extension ".root", the DATE file is converted
1033 // If deleteIntermediateFiles is true, the DATE file is deleted afterwards.
1035 TStopwatch stopwatch;
1038 if (!WriteRawFiles(detectors)) {
1039 if (fStopOnError) return kFALSE;
1042 TString dateFileName(fileName);
1043 if (!dateFileName.IsNull()) {
1044 Bool_t rootOutput = dateFileName.EndsWith(".root");
1045 if (rootOutput) dateFileName += ".date";
1046 if (!ConvertRawFilesToDate(dateFileName)) {
1047 if (fStopOnError) return kFALSE;
1049 if (deleteIntermediateFiles) {
1050 AliRunLoader* runLoader = LoadRun("READ");
1051 if (runLoader) for (Int_t iEvent = 0;
1052 iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1054 sprintf(command, "rm -r raw%d", iEvent);
1055 gSystem->Exec(command);
1060 if (!ConvertDateToRoot(dateFileName, fileName)) {
1061 if (fStopOnError) return kFALSE;
1063 if (deleteIntermediateFiles) {
1064 gSystem->Unlink(dateFileName);
1069 AliInfo(Form("Execution time: R:%.2fs C:%.2fs",
1070 stopwatch.RealTime(),stopwatch.CpuTime()));
1075 //_____________________________________________________________________________
1076 Bool_t AliSimulation::WriteRawFiles(const char* detectors)
1078 // convert the digits to raw data DDL files
1080 AliRunLoader* runLoader = LoadRun("READ");
1081 if (!runLoader) return kFALSE;
1083 // write raw data to DDL files
1084 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1085 AliInfo(Form("processing event %d", iEvent));
1086 runLoader->GetEvent(iEvent);
1087 TString baseDir = gSystem->WorkingDirectory();
1089 sprintf(dirName, "raw%d", iEvent);
1090 gSystem->MakeDirectory(dirName);
1091 if (!gSystem->ChangeDirectory(dirName)) {
1092 AliError(Form("couldn't change to directory %s", dirName));
1093 if (fStopOnError) return kFALSE; else continue;
1096 TString detStr = detectors;
1097 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1098 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1099 AliModule* det = (AliModule*) detArray->At(iDet);
1100 if (!det || !det->IsActive()) continue;
1101 if (IsSelected(det->GetName(), detStr)) {
1102 AliInfo(Form("creating raw data from digits for %s", det->GetName()));
1107 if (!WriteTriggerRawData())
1108 if (fStopOnError) return kFALSE;
1110 gSystem->ChangeDirectory(baseDir);
1111 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1112 AliError(Form("the following detectors were not found: %s",
1114 if (fStopOnError) return kFALSE;
1122 //_____________________________________________________________________________
1123 Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName)
1125 // convert raw data DDL files to a DATE file with the program "dateStream"
1127 char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
1129 AliError("the program dateStream was not found");
1130 if (fStopOnError) return kFALSE;
1135 AliRunLoader* runLoader = LoadRun("READ");
1136 if (!runLoader) return kFALSE;
1138 AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
1140 // Note the option -s. It is used in order to avoid
1141 // the generation of SOR/EOR events.
1142 sprintf(command, "dateStream -s -D -o %s -# %d -C",
1143 dateFileName, runLoader->GetNumberOfEvents());
1144 FILE* pipe = gSystem->OpenPipe(command, "w");
1146 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1147 fprintf(pipe, "GDC\n");
1151 // loop over detectors and DDLs
1152 for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
1153 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
1155 Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
1156 Int_t ldcID = Int_t(ldc + 0.0001);
1157 ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
1159 char rawFileName[256];
1160 sprintf(rawFileName, "raw%d/%s",
1161 iEvent, AliDAQ::DdlFileName(iDet,iDDL));
1163 // check existence and size of raw data file
1164 FILE* file = fopen(rawFileName, "rb");
1165 if (!file) continue;
1166 fseek(file, 0, SEEK_END);
1167 unsigned long size = ftell(file);
1169 if (!size) continue;
1171 if (ldcID != prevLDC) {
1172 fprintf(pipe, " LDC Id %d\n", ldcID);
1175 fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
1180 Int_t result = gSystem->ClosePipe(pipe);
1183 return (result == 0);
1186 //_____________________________________________________________________________
1187 Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
1188 const char* rootFileName)
1190 // convert a DATE file to a root file with the program "alimdc"
1193 const Int_t kDBSize = 2000000000;
1194 const Int_t kTagDBSize = 1000000000;
1195 const Bool_t kFilter = kFALSE;
1196 const Int_t kCompression = 0;
1198 char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
1200 AliError("the program alimdc was not found");
1201 if (fStopOnError) return kFALSE;
1206 AliInfo(Form("converting DATE file %s to root file %s",
1207 dateFileName, rootFileName));
1209 const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
1210 const char* tagDBFS = "/tmp/mdc1/tags";
1211 const char* runDBFS = "/tmp/mdc1/meta";
1213 // User defined file system locations
1214 if (gSystem->Getenv("ALIMDC_RAWDB1"))
1215 rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
1216 if (gSystem->Getenv("ALIMDC_RAWDB2"))
1217 rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
1218 if (gSystem->Getenv("ALIMDC_TAGDB"))
1219 tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
1220 if (gSystem->Getenv("ALIMDC_RUNDB"))
1221 runDBFS = gSystem->Getenv("ALIMDC_RUNDB");
1223 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1224 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1225 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1226 gSystem->Exec(Form("rm -rf %s",runDBFS));
1228 gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
1229 gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
1230 gSystem->Exec(Form("mkdir %s",tagDBFS));
1231 gSystem->Exec(Form("mkdir %s",runDBFS));
1233 Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
1234 kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
1235 gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
1237 gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
1238 gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
1239 gSystem->Exec(Form("rm -rf %s",tagDBFS));
1240 gSystem->Exec(Form("rm -rf %s",runDBFS));
1242 return (result == 0);
1246 //_____________________________________________________________________________
1247 AliRunLoader* AliSimulation::LoadRun(const char* mode) const
1249 // delete existing run loaders, open a new one and load gAlice
1251 while (AliRunLoader::GetRunLoader()) delete AliRunLoader::GetRunLoader();
1252 AliRunLoader* runLoader =
1253 AliRunLoader::Open(fGAliceFileName.Data(),
1254 AliConfig::GetDefaultEventFolderName(), mode);
1256 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1259 runLoader->LoadgAlice();
1260 gAlice = runLoader->GetAliRun();
1262 AliError(Form("no gAlice object found in file %s",
1263 fGAliceFileName.Data()));
1269 //_____________________________________________________________________________
1270 Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
1272 // get or calculate the number of signal events per background event
1274 if (!fBkgrdFileNames) return 1;
1275 Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
1276 if (nBkgrdFiles == 0) return 1;
1278 // get the number of signal events
1280 AliRunLoader* runLoader =
1281 AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
1282 if (!runLoader) return 1;
1284 nEvents = runLoader->GetNumberOfEvents();
1289 for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
1290 // get the number of background events
1291 const char* fileName = ((TObjString*)
1292 (fBkgrdFileNames->At(iBkgrdFile)))->GetName();
1293 AliRunLoader* runLoader =
1294 AliRunLoader::Open(fileName, "BKGRD");
1295 if (!runLoader) continue;
1296 Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
1299 // get or calculate the number of signal per background events
1300 Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
1301 if (nSignalPerBkgrd <= 0) {
1302 nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
1303 } else if (result && (result != nSignalPerBkgrd)) {
1304 AliInfo(Form("the number of signal events per background event "
1305 "will be changed from %d to %d for stream %d",
1306 nSignalPerBkgrd, result, iBkgrdFile+1));
1307 nSignalPerBkgrd = result;
1310 if (!result) result = nSignalPerBkgrd;
1311 if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
1312 AliWarning(Form("not enough background events (%d) for %d signal events "
1313 "using %d signal per background events for stream %d",
1314 nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
1321 //_____________________________________________________________________________
1322 Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
1324 // check whether detName is contained in detectors
1325 // if yes, it is removed from detectors
1327 // check if all detectors are selected
1328 if ((detectors.CompareTo("ALL") == 0) ||
1329 detectors.BeginsWith("ALL ") ||
1330 detectors.EndsWith(" ALL") ||
1331 detectors.Contains(" ALL ")) {
1336 // search for the given detector
1337 Bool_t result = kFALSE;
1338 if ((detectors.CompareTo(detName) == 0) ||
1339 detectors.BeginsWith(detName+" ") ||
1340 detectors.EndsWith(" "+detName) ||
1341 detectors.Contains(" "+detName+" ")) {
1342 detectors.ReplaceAll(detName, "");
1346 // clean up the detectors string
1347 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1348 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1349 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1354 Bool_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName)
1357 // Steering routine to convert raw data in directory rawDirectory/ to fake SDigits.
1358 // These can be used for embedding of MC tracks into RAW data using the standard
1359 // merging procedure.
1361 // If an ESD file is given the reconstructed vertex is taken from it and stored in the event header.
1364 AliError("no gAlice object. Restart aliroot and try again.");
1367 if (gAlice->Modules()->GetEntries() > 0) {
1368 AliError("gAlice was already run. Restart aliroot and try again.");
1372 AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
1373 StdoutToAliInfo(StderrToAliError(gAlice->Init(fConfigFileName.Data());););
1377 AliCDBManager* man = AliCDBManager::Instance();
1378 man->SetRun(0); // Should this come from rawdata header ?
1382 // Get the runloader
1383 AliRunLoader* runLoader = gAlice->GetRunLoader();
1385 // Open esd file if available
1386 TFile* esdFile = TFile::Open(esdFileName);
1387 Bool_t esdOK = (esdFile != 0);
1388 AliESD* esd = new AliESD;
1391 treeESD = (TTree*) esdFile->Get("esdTree");
1393 AliWarning("No ESD tree found");
1396 treeESD->SetBranchAddress("ESD", &esd);
1400 // Create the RawReader
1401 AliRawReaderFile* rawReader = new AliRawReaderFile(rawDirectory);
1403 // Get list of detectors
1404 TObjArray* detArray = runLoader->GetAliRun()->Detectors();
1407 AliHeader* header = runLoader->GetHeader();
1412 if (!(rawReader->NextEvent())) break;
1415 for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
1416 AliModule* det = (AliModule*) detArray->At(iDet);
1417 AliInfo(Form("Calling Raw2SDigits for %s\n", det->GetName()));
1418 det->Raw2SDigits(rawReader);
1423 // If ESD information available obtain reconstructed vertex and store in header.
1425 treeESD->GetEvent(nev);
1426 const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
1427 Double_t position[3];
1428 esdVertex->GetXYZ(position);
1429 AliGenEventHeader* mcHeader = new AliGenEventHeader("ESD");
1432 for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
1433 mcHeader->SetPrimaryVertex(mcV);
1434 header->Reset(0,nev);
1435 header->SetGenEventHeader(mcHeader);
1436 printf("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]);
1441 runLoader->TreeE()->Fill();
1442 runLoader->SetNextEvent();
1448 runLoader->CdGAFile();
1449 runLoader->WriteHeader("OVERWRITE");
1450 runLoader->WriteRunLoader();