#include <TFile.h>
#include <TGeoGlobalMagField.h>
#include <TGeoManager.h>
#include <TObjString.h>
#include <TROOT.h>
#include <TSystem.h>
#include <TVirtualMC.h>
#include <TVirtualMCApplication.h>
#include <TDatime.h>
#include <TInterpreter.h>
#include "AliAlignObj.h"
#include "AliCDBEntry.h"
#include "AliCDBManager.h"
#include "AliGRPManager.h"
#include "AliCDBStorage.h"
#include "AliCTPRawData.h"
#include "AliCentralTrigger.h"
#include "AliCentralTrigger.h"
#include "AliCodeTimer.h"
#include "AliDAQ.h"
#include "AliDigitizer.h"
#include "AliESDEvent.h"
#include "AliGRPObject.h"
#include "AliGenEventHeader.h"
#include "AliGenerator.h"
#include "AliGeomManager.h"
#include "AliHLTSimulation.h"
#include "AliHeader.h"
#include "AliLego.h"
#include "AliLegoGenerator.h"
#include "AliLog.h"
#include "AliMC.h"
#include "AliMagF.h"
#include "AliModule.h"
#include "AliPDG.h"
#include "AliRawReaderDate.h"
#include "AliRawReaderFile.h"
#include "AliRawReaderRoot.h"
#include "AliRun.h"
#include "AliDigitizationInput.h"
#include "AliRunLoader.h"
#include "AliStack.h"
#include "AliSimulation.h"
#include "AliSysInfo.h"
#include "AliVertexGenFile.h"
using std::ofstream;
ClassImp(AliSimulation)
AliSimulation *AliSimulation::fgInstance = 0;
const char* AliSimulation::fgkDetectorName[AliSimulation::fgkNDetectors] = {"ITS", "TPC", "TRD",
"TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE","AD",
"FIT","MFT","HLT"};
AliSimulation::AliSimulation(const char* configFileName,
const char* name, const char* title) :
TNamed(name, title),
fRunGeneratorOnly(kFALSE),
fRunGeneration(kTRUE),
fRunSimulation(kTRUE),
fLoadAlignFromCDB(kTRUE),
fLoadAlObjsListOfDets("ALL"),
fMakeSDigits("ALL"),
fMakeDigits("ALL"),
fTriggerConfig(""),
fMakeDigitsFromHits(""),
fWriteRawData(""),
fRawDataFileName(""),
fDeleteIntermediateFiles(kFALSE),
fWriteSelRawData(kFALSE),
fStopOnError(kFALSE),
fUseMonitoring(kFALSE),
fNEvents(1),
fConfigFileName(configFileName),
fGAliceFileName("galice.root"),
fEventsPerFile(),
fBkgrdFileNames(NULL),
fAlignObjArray(NULL),
fUseBkgrdVertex(kTRUE),
fRegionOfInterest(kFALSE),
fCDBUri(""),
fQARefUri(""),
fSpecCDBUri(),
fRun(-1),
fSeed(0),
fInitCDBCalled(kFALSE),
fInitRunNumberCalled(kFALSE),
fSetRunNumberFromDataCalled(kFALSE),
fEmbeddingFlag(kFALSE),
fLego(NULL),
fKey(0),
fUseVertexFromCDB(0),
fUseMagFieldFromGRP(0),
fGRPWriteLocation(Form("local://%s", gSystem->pwd())),
fUseTimeStampFromCDB(0),
fTimeStart(0),
fTimeEnd(0),
fQADetectors("ALL"),
fQATasks("ALL"),
fRunQA(kTRUE),
fEventSpecie(AliRecoParam::kDefault),
fWriteQAExpertData(kTRUE),
fGeometryFile(),
fRunHLT("default"),
fpHLT(NULL),
fWriteGRPEntry(kTRUE)
{
fgInstance = this;
SetGAliceFile("galice.root");
AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
qam->SetActiveDetectors(fQADetectors) ;
fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
qam->SetTasks(fQATasks) ;
}
AliSimulation::~AliSimulation()
{
fEventsPerFile.Delete();
if (fBkgrdFileNames) {
fBkgrdFileNames->Delete();
delete fBkgrdFileNames;
}
fSpecCDBUri.Delete();
if (fgInstance==this) fgInstance = 0;
AliQAManager::QAManager()->ShowQA() ;
AliQAManager::Destroy() ;
AliCodeTimer::Instance()->Print();
}
void AliSimulation::SetNumberOfEvents(Int_t nEvents)
{
fNEvents = nEvents;
}
void AliSimulation::InitQA()
{
if (fInitCDBCalled) return;
fInitCDBCalled = kTRUE;
AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kSIMMODE) ;
qam->SetActiveDetectors(fQADetectors) ;
fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
qam->SetTasks(fQATasks) ;
if (fWriteQAExpertData)
qam->SetWriteExpert() ;
if (qam->IsDefaultStorageSet()) {
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
AliWarning("Default QA reference storage has been already set !");
AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fQARefUri.Data()));
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
fQARefUri = qam->GetDefaultStorage()->GetURI();
} else {
if (fQARefUri.Length() > 0) {
AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
} else {
fQARefUri="local://$ALICE_ROOT/QARef";
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
AliWarning("Default QA reference storage not yet set !!!!");
AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
}
qam->SetDefaultStorage(fQARefUri);
}
}
void AliSimulation::InitCDB()
{
if (fInitCDBCalled) return;
fInitCDBCalled = kTRUE;
AliCDBManager* man = AliCDBManager::Instance();
if (man->IsDefaultStorageSet())
{
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
AliWarning("Default CDB storage has been already set !");
AliWarning(Form("Ignoring the default storage declared in AliSimulation: %s",fCDBUri.Data()));
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
fCDBUri = man->GetDefaultStorage()->GetURI();
}
else {
if (fCDBUri.Length() > 0)
{
AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
} else {
fCDBUri="local://$ALICE_ROOT/OCDB";
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
AliWarning("Default CDB storage not yet set !!!!");
AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
}
man->SetDefaultStorage(fCDBUri);
}
for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
TObject* obj = fSpecCDBUri[i];
if (!obj) continue;
AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
}
}
void AliSimulation::InitRunNumber(){
if (fInitRunNumberCalled) return;
fInitRunNumberCalled = kTRUE;
AliCDBManager* man = AliCDBManager::Instance();
if (man->GetRun() >= 0)
{
AliFatal(Form("Run number cannot be set in AliCDBManager before start of simulation: "
"Use external variable DC_RUN or AliSimulation::SetRun()!"));
}
if(fRun >= 0) {
AliDebug(2, Form("Setting CDB run number to: %d",fRun));
} else {
fRun=0;
AliWarning(Form("Run number not yet set !!!! Setting it now to: %d",
fRun));
}
man->SetRun(fRun);
man->Print();
}
void AliSimulation::SetCDBLock() {
ULong64_t key = AliCDBManager::Instance()->SetLock(1);
if (key) fKey = key;
}
void AliSimulation::SetDefaultStorage(const char* uri) {
fCDBUri = uri;
}
void AliSimulation::SetQARefDefaultStorage(const char* uri) {
fQARefUri = uri;
AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
}
void AliSimulation::SetSpecificStorage(const char* calibType, const char* uri) {
AliCDBPath aPath(calibType);
if(!aPath.IsValid()){
AliError(Form("Not a valid path: %s", calibType));
return;
}
TObject* obj = fSpecCDBUri.FindObject(calibType);
if (obj) fSpecCDBUri.Remove(obj);
fSpecCDBUri.Add(new TNamed(calibType, uri));
}
void AliSimulation::SetRunNumber(Int_t run)
{
fRun = run;
}
void AliSimulation::SetSeed(Int_t seed)
{
fSeed = seed;
}
Bool_t AliSimulation::SetRunNumberFromData()
{
if (fSetRunNumberFromDataCalled) return kTRUE;
fSetRunNumberFromDataCalled = kTRUE;
AliCDBManager* man = AliCDBManager::Instance();
Int_t runData = -1, runCDB = -1;
AliRunLoader* runLoader = LoadRun("READ");
if (!runLoader) return kFALSE;
else {
runData = runLoader->GetHeader()->GetRun();
delete runLoader;
}
runCDB = man->GetRun();
if(runCDB >= 0) {
if (runCDB != runData) {
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
AliWarning(Form("A run number was previously set in AliCDBManager: %d !", runCDB));
AliWarning(Form("It will be replaced with the run number got from run header: %d !", runData));
AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
}
}
man->SetRun(runData);
fRun = runData;
if(man->GetRun() < 0) {
AliError("Run number not properly initalized!");
return kFALSE;
}
man->Print();
return kTRUE;
}
void AliSimulation::SetConfigFile(const char* fileName)
{
fConfigFileName = fileName;
}
void AliSimulation::SetGAliceFile(const char* fileName)
{
fGAliceFileName = fileName;
if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) {
char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
fGAliceFileName);
fGAliceFileName = absFileName;
delete[] absFileName;
}
AliDebug(2, Form("galice file name set to %s", fileName));
}
void AliSimulation::SetEventsPerFile(const char* detector, const char* type,
Int_t nEvents)
{
TNamed* obj = new TNamed(detector, type);
obj->SetUniqueID(nEvents);
fEventsPerFile.Add(obj);
}
Bool_t AliSimulation::MisalignGeometry(AliRunLoader *runLoader)
{
if (!AliGeomManager::GetGeometry() || !AliGeomManager::GetGeometry()->IsClosed()) {
AliError("Can't apply the misalignment! Geometry is not loaded or it is still opened!");
return kFALSE;
}
InitCDB();
SetCDBLock();
Bool_t delRunLoader = kFALSE;
if (!runLoader) {
runLoader = LoadRun("READ");
if (!runLoader) return kFALSE;
delRunLoader = kTRUE;
}
if(!IsGeometryFromFile()) AliGeomManager::GetGeometry()->Export("geometry.root");
if(fLoadAlignFromCDB){
TString detStr = fLoadAlObjsListOfDets;
TString loadAlObjsListOfDets = "";
TObjArray* detArray = runLoader->GetAliRun()->Detectors();
for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
AliModule* det = (AliModule*) detArray->At(iDet);
if (!det || !det->IsActive()) continue;
if (IsSelected(det->GetName(), detStr)) {
loadAlObjsListOfDets += det->GetName();
loadAlObjsListOfDets += " ";
}
}
loadAlObjsListOfDets.Prepend("GRP ");
AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
}else{
if (fAlignObjArray) {
if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
AliError("The misalignment of one or more volumes failed!"
"Compare the list of simulated detectors and the list of detector alignment data!");
if (delRunLoader) delete runLoader;
return kFALSE;
}
}
}
TString detStr = fLoadAlObjsListOfDets;
TObjArray* detArray = runLoader->GetAliRun()->Detectors();
for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
AliModule* det = (AliModule*) detArray->At(iDet);
if (!det || !det->IsActive()) continue;
if (IsSelected(det->GetName(), detStr)) {
det->UpdateInternalGeometry();
}
}
if (delRunLoader) delete runLoader;
return kTRUE;
}
void AliSimulation::MergeWith(const char* fileName, Int_t nSignalPerBkgrd)
{
TObjString* fileNameStr = new TObjString(fileName);
fileNameStr->SetUniqueID(nSignalPerBkgrd);
if (!fBkgrdFileNames) fBkgrdFileNames = new TObjArray;
fBkgrdFileNames->Add(fileNameStr);
}
void AliSimulation::EmbedInto(const char* fileName, Int_t nSignalPerBkgrd)
{
MergeWith(fileName, nSignalPerBkgrd);
fEmbeddingFlag = kTRUE;
}
Bool_t AliSimulation::Run(Int_t nEvents)
{
AliCodeTimerAuto("",0)
AliSysInfo::AddStamp("Start_Run");
ProcessEnvironmentVars();
AliSysInfo::AddStamp("ProcessEnvironmentVars");
gRandom->SetSeed(fSeed);
if (nEvents > 0) fNEvents = nEvents;
if (fRunGeneratorOnly)
{
if(!RunGeneratorOnly())
{
if (fStopOnError) return kFALSE;
}
else
return kTRUE;
}
if (!fRunHLT.IsNull() && !CreateHLT()) {
if (fStopOnError) return kFALSE;
fRunHLT="";
}
if (fRunGeneration) {
if (!RunSimulation()) if (fStopOnError) return kFALSE;
}
AliSysInfo::AddStamp("RunSimulation");
InitCDB();
AliSysInfo::AddStamp("InitCDB");
if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
SetCDBLock();
if (!AliGeomManager::GetGeometry()) {
AliGeomManager::LoadGeometry("geometry.root");
AliSysInfo::AddStamp("GetGeometry");
if(!AliGeomManager::CheckSymNamesLUT("ALL"))
AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
if (!AliGeomManager::GetGeometry()) if (fStopOnError) return kFALSE;
if(!MisalignGeometry()) if (fStopOnError) return kFALSE;
}
AliSysInfo::AddStamp("MissalignGeometry");
AliSysInfo::AddStamp("Start_sdigitization");
if (!fMakeSDigits.IsNull()) {
if (!RunSDigitization(fMakeSDigits)) if (fStopOnError) return kFALSE;
}
AliSysInfo::AddStamp("Stop_sdigitization");
AliSysInfo::AddStamp("Start_digitization");
if (!fMakeDigits.IsNull()) {
if (!RunDigitization(fMakeDigits, fMakeDigitsFromHits)) {
if (fStopOnError) return kFALSE;
}
}
AliSysInfo::AddStamp("Stop_digitization");
if (!fMakeDigitsFromHits.IsNull()) {
if (fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
AliWarning(Form("Merging and direct creation of digits from hits "
"was selected for some detectors. "
"No merging will be done for the following detectors: %s",
fMakeDigitsFromHits.Data()));
}
if (!RunHitsDigitization(fMakeDigitsFromHits)) {
if (fStopOnError) return kFALSE;
}
}
AliSysInfo::AddStamp("Hits2Digits");
if (!fTriggerConfig.EqualTo("none",TString::kIgnoreCase) &&
!RunTrigger(fTriggerConfig,fMakeDigits)) {
if (fStopOnError) return kFALSE;
}
AliSysInfo::AddStamp("RunTrigger");
if (!fWriteRawData.IsNull()) {
if (!WriteRawData(fWriteRawData, fRawDataFileName,
fDeleteIntermediateFiles,fWriteSelRawData)) {
if (fStopOnError) return kFALSE;
}
}
AliSysInfo::AddStamp("WriteRaw");
if (!fRunHLT.IsNull() && fWriteRawData.IsNull()) {
if (!RunHLT()) {
if (fStopOnError) return kFALSE;
}
}
AliSysInfo::AddStamp("RunHLT");
if (fRunQA) {
Bool_t rv = RunQA() ;
if (!rv)
if (fStopOnError)
return kFALSE ;
}
AliSysInfo::AddStamp("RunQA");
StoreUsedCDBMaps();
TString snapshotFileOut("");
if(TString(gSystem->Getenv("OCDB_SNAPSHOT_CREATE")) == TString("kTRUE")){
AliInfo(" ******** Creating the snapshot! *********");
TString snapshotFile(gSystem->Getenv("OCDB_SNAPSHOT_FILENAME"));
if(!(snapshotFile.IsNull() || snapshotFile.IsWhitespace()))
snapshotFileOut = snapshotFile;
else
snapshotFileOut="OCDB.root";
AliCDBManager::Instance()->DumpToSnapshotFile(snapshotFileOut.Data(),kFALSE);
}
AliCDBManager::Instance()->ClearCache();
return kTRUE;
}
Bool_t AliSimulation::RunLego(const char *setup, Int_t nc1, Float_t c1min,
Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener, Int_t nev)
{
/*
<img src="picts/AliRunLego1.gif">
*/
//End_Html
/*
<img src="picts/AliRunLego2.gif">
*/
//End_Html
/*
<img src="picts/AliRunLego3.gif">
*/
//End_Html
AliCodeTimerAuto("",0)
InitCDB();
InitRunNumber();
SetCDBLock();
if (!gAlice) {
AliError("no gAlice object. Restart aliroot and try again.");
return kFALSE;
}
if (gAlice->Modules()->GetEntries() > 0) {
AliError("gAlice was already run. Restart aliroot and try again.");
return kFALSE;
}
AliInfo(Form("initializing gAlice with config file %s",
fConfigFileName.Data()));
if (nev == -1) nev = nc1 * nc2;
if (!gener) gener = new AliLegoGenerator();
gener->SetRadiusRange(rmin, rmax);
gener->SetZMax(zmax);
gener->SetCoor1Range(nc1, c1min, c1max);
gener->SetCoor2Range(nc2, c2min, c2max);
fLego = new AliLego("lego",gener);
gAlice->Announce();
if (fUseMagFieldFromGRP) {
AliGRPManager grpM;
grpM.ReadGRPEntry();
grpM.SetMagField();
AliInfo("Field is locked now. It cannot be changed in Config.C");
}
gROOT->LoadMacro(setup);
gInterpreter->ProcessLine(gAlice->GetConfigFunction());
if(AliCDBManager::Instance()->GetRun() >= 0) {
SetRunNumber(AliCDBManager::Instance()->GetRun());
} else {
AliWarning("Run number not initialized!!");
}
AliRunLoader::Instance()->CdGAFile();
AliPDG::AddParticlesToPdgDataBase();
TVirtualMC::GetMC()->SetMagField(TGeoGlobalMagField::Instance()->GetField());
gAlice->GetMCApp()->Init();
gAlice->InitLoaders();
AliRunLoader::Instance()->MakeTree("E");
AliRunLoader::Instance()->CdGAFile();
gAlice->Write();
AliGenerator *gen=gAlice->GetMCApp()->Generator();
gAlice->GetMCApp()->ResetGenerator(gener);
TVirtualMC::GetMC()->InitLego();
AliRunLoader::Instance()->SetNumberOfEventsPerFile(nev);
TVirtualMC::GetMC()->ProcessRun(nev);
FinishRun();
gAlice->GetMCApp()->ResetGenerator(gen);
delete fLego;
return kTRUE;
}
Bool_t AliSimulation::RunTrigger(const char* config, const char* detectors)
{
AliCodeTimerAuto("",0)
InitCDB();
if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
SetCDBLock();
AliRunLoader* runLoader = LoadRun("READ");
if (!runLoader) return kFALSE;
TString trconfiguration = config;
if (trconfiguration.IsNull()) {
if(!fTriggerConfig.IsNull()) {
trconfiguration = fTriggerConfig;
}
else
AliWarning("No trigger descriptor is specified. Loading the one that is in the CDB.");
}
runLoader->MakeTree( "GG" );
AliCentralTrigger* aCTP = runLoader->GetTrigger();
if (!aCTP->LoadConfiguration( trconfiguration ))
return kFALSE;
if( !aCTP->RunTrigger( runLoader , detectors ) ) {
if (fStopOnError) {
return kFALSE;
}
}
delete runLoader;
return kTRUE;
}
Bool_t AliSimulation::WriteTriggerRawData()
{
AliCTPRawData writer;
writer.RawDataRun2();
return kTRUE;
}
Bool_t AliSimulation::RunSimulation(Int_t nEvents)
{
AliCodeTimerAuto("",0)
AliSysInfo::AddStamp("RunSimulation_Begin");
InitCDB();
AliSysInfo::AddStamp("RunSimulation_InitCDB");
InitRunNumber();
SetCDBLock();
AliSysInfo::AddStamp("RunSimulation_SetCDBLock");
if (!gAlice) {
AliError("no gAlice object. Restart aliroot and try again.");
return kFALSE;
}
if (gAlice->Modules()->GetEntries() > 0) {
AliError("gAlice was already run. Restart aliroot and try again.");
return kFALSE;
}
gAlice->GetMCApp()->SetUseMonitoring(fUseMonitoring);
AliInfo(Form("initializing gAlice with config file %s",
fConfigFileName.Data()));
gAlice->Announce();
if (fUseMagFieldFromGRP) {
AliGRPManager grpM;
grpM.ReadGRPEntry();
grpM.SetMagField();
AliInfo("Field is locked now. It cannot be changed in Config.C");
}
TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
gROOT->LoadMacro(fConfigFileName.Data());
Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
}
AliSysInfo::AddStamp("RunSimulation_Config");
if (fUseVertexFromCDB) {
Double_t vtxPos[3] = {0., 0., 0.};
Double_t vtxSig[3] = {0., 0., 0.};
AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
if (entry) {
AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
if (vertex) {
if(vertex->GetXRes()>2.8) {
entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
if (entry) vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
}
}
if (vertex) {
vertex->GetXYZ(vtxPos);
vertex->GetSigmaXYZ(vtxSig);
AliInfo("Overwriting Config.C vertex settings !");
AliInfo(Form("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
vtxPos[0], vtxPos[1], vtxPos[2], vtxSig[2]));
AliGenerator *gen = gAlice->GetMCApp()->Generator();
gen->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]);
gen->SetSigmaZ(vtxSig[2]);
}
}
}
if (fUseTimeStampFromCDB) {
AliGRPManager grpM;
grpM.ReadGRPEntry();
const AliGRPObject *grpObj = grpM.GetGRPData();
if (!grpObj || (grpObj->GetTimeEnd() <= grpObj->GetTimeStart())) {
AliError("Missing GRP or bad SOR/EOR time-stamps! Switching off the time-stamp generation from GRP!");
fTimeStart = fTimeEnd = 0;
fUseTimeStampFromCDB = kFALSE;
}
else {
fTimeStart = grpObj->GetTimeStart();
fTimeEnd = grpObj->GetTimeEnd();
}
}
if(AliCDBManager::Instance()->GetRun() >= 0) {
AliRunLoader::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
AliRunLoader::Instance()->SetNumberOfEventsPerRun(fNEvents);
} else {
AliWarning("Run number not initialized!!");
}
AliRunLoader::Instance()->CdGAFile();
AliPDG::AddParticlesToPdgDataBase();
TVirtualMC::GetMC()->SetMagField(TGeoGlobalMagField::Instance()->GetField());
AliSysInfo::AddStamp("RunSimulation_GetField");
gAlice->GetMCApp()->Init();
AliSysInfo::AddStamp("RunSimulation_InitMCApp");
gAlice->InitLoaders();
AliRunLoader::Instance()->MakeTree("E");
AliRunLoader::Instance()->LoadKinematics("RECREATE");
AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
AliRunLoader::Instance()->LoadHits("all","RECREATE");
AliRunLoader::Instance()->CdGAFile();
gAlice->Write();
gAlice->SetEventNrInRun(-1);
AliSysInfo::AddStamp("RunSimulation_InitLoaders");
AliSysInfo::AddStamp("RunSimulation_TriggerDescriptor");
AliInfo(Form("Run number: %d",AliCDBManager::Instance()->GetRun()));
AliRunLoader* runLoader = AliRunLoader::Instance();
if (!runLoader) {
AliError(Form("gAlice has no run loader object. "
"Check your config file: %s", fConfigFileName.Data()));
return kFALSE;
}
SetGAliceFile(runLoader->GetFileName());
#if ROOT_VERSION_CODE < 331527
AliGeomManager::SetGeometry(gGeoManager);
TString detsToBeChecked = "";
TObjArray* detArray = runLoader->GetAliRun()->Detectors();
for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
AliModule* det = (AliModule*) detArray->At(iDet);
if (!det || !det->IsActive()) continue;
detsToBeChecked += det->GetName();
detsToBeChecked += " ";
}
if(!AliGeomManager::CheckSymNamesLUT(detsToBeChecked.Data()))
AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
MisalignGeometry(runLoader);
AliSysInfo::AddStamp("RunSimulation_MisalignGeometry");
#endif
if (!gAlice->GetMCApp()->Generator()) {
AliError(Form("gAlice has no generator object. "
"Check your config file: %s", fConfigFileName.Data()));
return kFALSE;
}
if (fWriteGRPEntry)
WriteGRPEntry();
AliSysInfo::AddStamp("RunSimulation_WriteGRP");
if (nEvents <= 0) nEvents = fNEvents;
if (fUseBkgrdVertex &&
fBkgrdFileNames && (fBkgrdFileNames->GetEntriesFast() > 0)) {
Int_t signalPerBkgrd = GetNSignalPerBkgrd(nEvents);
const char* fileName = ((TObjString*)
(fBkgrdFileNames->At(0)))->GetName();
AliInfo(Form("The vertex will be taken from the background "
"file %s with nSignalPerBackground = %d",
fileName, signalPerBkgrd));
AliVertexGenFile* vtxGen = new AliVertexGenFile(fileName, signalPerBkgrd);
gAlice->GetMCApp()->Generator()->SetVertexGenerator(vtxGen);
}
if (!fRunSimulation) {
gAlice->GetMCApp()->Generator()->SetTrackingFlag(0);
}
for (Int_t i = 0; i < fEventsPerFile.GetEntriesFast(); i++) {
if (!fEventsPerFile[i]) continue;
const char* detName = fEventsPerFile[i]->GetName();
const char* typeName = fEventsPerFile[i]->GetTitle();
TString loaderName(detName);
loaderName += "Loader";
AliLoader* loader = runLoader->GetLoader(loaderName);
if (!loader) {
AliError(Form("RunSimulation no loader for %s found\n Number of events per file not set for %s %s",
detName, typeName, detName));
continue;
}
AliDataLoader* dataLoader =
loader->GetDataLoader(typeName);
if (!dataLoader) {
AliError(Form("no data loader for %s found\n"
"Number of events per file not set for %s %s",
typeName, detName, typeName));
continue;
}
dataLoader->SetNumberOfEventsPerFile(fEventsPerFile[i]->GetUniqueID());
AliDebug(1, Form("number of events per file set to %d for %s %s",
fEventsPerFile[i]->GetUniqueID(), detName, typeName));
}
AliInfo("running gAlice");
AliSysInfo::AddStamp("Start_ProcessRun");
TVirtualMC::GetMC()->ProcessRun(nEvents);
if(nEvents>0) FinishRun();
AliSysInfo::AddStamp("Stop_ProcessRun");
delete runLoader;
return kTRUE;
}
Bool_t AliSimulation::RunGeneratorOnly()
{
TInterpreter::EErrorCode interpreterError=TInterpreter::kNoError;
gROOT->LoadMacro(fConfigFileName.Data());
Long_t interpreterResult=gInterpreter->ProcessLine(gAlice->GetConfigFunction(), &interpreterError);
if (interpreterResult!=0 || interpreterError!=TInterpreter::kNoError) {
AliFatal(Form("execution of config file \"%s\" failed with error %d", fConfigFileName.Data(), (int)interpreterError));
}
AliRunLoader* runLoader = AliRunLoader::Instance();
AliGenerator* generator = gAlice->GetMCApp()->Generator();
if (!runLoader) {
AliError(Form("gAlice has no run loader object. "
"Check your config file: %s", fConfigFileName.Data()));
return kFALSE;
}
if (!generator) {
AliError(Form("gAlice has no generator object. "
"Check your config file: %s", fConfigFileName.Data()));
return kFALSE;
}
runLoader->LoadKinematics("RECREATE");
runLoader->MakeTree("E");
runLoader->MakeStack();
AliStack* stack = runLoader->Stack();
AliHeader* header = runLoader->GetHeader();
generator->Init();
generator->SetStack(stack);
for (Int_t iev=0; iev<fNEvents; iev++)
{
header->Reset(0,iev);
runLoader->SetEventNumber(iev);
stack->Reset();
runLoader->MakeTree("K");
generator->Generate();
header->SetNprimary(stack->GetNprimary());
header->SetNtrack(stack->GetNtrack());
stack->FinishEvent();
header->SetStack(stack);
runLoader->TreeE()->Fill();
runLoader->WriteKinematics("OVERWRITE");
}
generator->FinishRun();
runLoader->WriteHeader("OVERWRITE");
generator->Write();
runLoader->Write();
return kTRUE;
}
Bool_t AliSimulation::RunSDigitization(const char* detectors)
{
static Int_t eventNr=0;
AliCodeTimerAuto("",0) ;
InitCDB();
if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
SetCDBLock();
AliRunLoader* runLoader = LoadRun();
if (!runLoader) return kFALSE;
TString detStr = detectors;
TObjArray* detArray = runLoader->GetAliRun()->Detectors();
for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
AliModule* det = (AliModule*) detArray->At(iDet);
if (!det || !det->IsActive()) continue;
if (IsSelected(det->GetName(), detStr)) {
AliInfo(Form("creating summable digits for %s", det->GetName()));
AliCodeTimerStart(Form("creating summable digits for %s", det->GetName()));
det->Hits2SDigits();
AliCodeTimerStop(Form("creating summable digits for %s", det->GetName()));
AliSysInfo::AddStamp(Form("SDigit_%s_%d",det->GetName(),eventNr), 0,1, eventNr);
}
}
if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
AliError(Form("the following detectors were not found: %s",
detStr.Data()));
if (fStopOnError) return kFALSE;
}
eventNr++;
delete runLoader;
return kTRUE;
}
Bool_t AliSimulation::RunDigitization(const char* detectors,
const char* excludeDetectors)
{
AliCodeTimerAuto("",0)
InitCDB();
if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
SetCDBLock();
delete AliRunLoader::Instance();
delete gAlice;
gAlice = NULL;
Int_t nStreams = 1;
if (fBkgrdFileNames) nStreams = fBkgrdFileNames->GetEntriesFast() + 1;
Int_t signalPerBkgrd = GetNSignalPerBkgrd();
AliDigitizationInput digInp(nStreams, signalPerBkgrd);
digInp.SetRegionOfInterest(fRegionOfInterest);
digInp.SetInputStream(0, fGAliceFileName.Data());
for (Int_t iStream = 1; iStream < nStreams; iStream++) {
const char* fileName = ((TObjString*)(fBkgrdFileNames->At(iStream-1)))->GetName();
digInp.SetInputStream(iStream, fileName);
}
TObjArray detArr;
detArr.SetOwner(kTRUE);
TString detStr = detectors;
TString detExcl = excludeDetectors;
if (!static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice()) {
AliError("Error occured while getting gAlice from Input 0");
return kFALSE;
}
AliRunLoader* runLoader = AliRunLoader::GetRunLoader(digInp.GetInputStream(0)->GetFolderName());
TObjArray* detArray = runLoader->GetAliRun()->Detectors();
for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
AliModule* det = (AliModule*) detArray->At(iDet);
if (!det || !det->IsActive()) continue;
if (!IsSelected(det->GetName(), detStr) || IsSelected(det->GetName(), detExcl)) continue;
AliDigitizer* digitizer = det->CreateDigitizer(&digInp);
if (!digitizer || !digitizer->Init()) {
AliError(Form("no digitizer for %s", det->GetName()));
if (fStopOnError) return kFALSE;
else continue;
}
detArr.AddLast(digitizer);
AliInfo(Form("Created digitizer from SDigits -> Digits for %s", det->GetName()));
}
if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
AliError(Form("the following detectors were not found: %s", detStr.Data()));
if (fStopOnError) return kFALSE;
}
Int_t ndigs = detArr.GetEntriesFast();
Int_t eventsCreated = 0;
AliRunLoader* outRl = digInp.GetOutRunLoader();
while ((eventsCreated++ < fNEvents) || (fNEvents < 0)) {
if (!digInp.ConnectInputTrees()) break;
digInp.InitEvent();
if (outRl) outRl->SetEventNumber(eventsCreated-1);
static_cast<AliStream*>(digInp.GetInputStream(0))->ImportgAlice();
for (int id=0;id<ndigs;id++) {
((AliDigitizer*)detArr[id])->Digitize("");
AliSysInfo::AddStamp(Form("Digit_%s_%d",detArr[id]->GetName(),eventsCreated), 0,2, eventsCreated);
}
digInp.FinishEvent();
};
digInp.FinishGlobal();
return kTRUE;
}
Bool_t AliSimulation::RunHitsDigitization(const char* detectors)
{
AliCodeTimerAuto("",0)
InitCDB();
if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
SetCDBLock();
AliRunLoader* runLoader = LoadRun("READ");
if (!runLoader) return kFALSE;
TString detStr = detectors;
TObjArray* detArray = runLoader->GetAliRun()->Detectors();
for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
AliModule* det = (AliModule*) detArray->At(iDet);
if (!det || !det->IsActive()) continue;
if (IsSelected(det->GetName(), detStr)) {
AliInfo(Form("creating digits from hits for %s", det->GetName()));
det->Hits2Digits();
}
}
if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
AliError(Form("the following detectors were not found: %s",
detStr.Data()));
if (fStopOnError) return kFALSE;
}
return kTRUE;
}
Bool_t AliSimulation::WriteRawData(const char* detectors,
const char* fileName,
Bool_t deleteIntermediateFiles,
Bool_t selrawdata)
{
AliCodeTimerAuto("",0)
AliSysInfo::AddStamp("WriteRawData_Start");
TString detStr = detectors;
if (!WriteRawFiles(detStr.Data())) {
if (fStopOnError) return kFALSE;
}
AliSysInfo::AddStamp("WriteRawFiles");
detStr.ReplaceAll("HLT", "");
if (!fRunHLT.IsNull()) {
if (!RunHLT()) {
if (fStopOnError) return kFALSE;
}
}
AliSysInfo::AddStamp("WriteRawData_RunHLT");
TString dateFileName(fileName);
if (!dateFileName.IsNull()) {
Bool_t rootOutput = dateFileName.EndsWith(".root");
if (rootOutput) dateFileName += ".date";
TString selDateFileName;
if (selrawdata) {
selDateFileName = "selected.";
selDateFileName+= dateFileName;
}
if (!ConvertRawFilesToDate(dateFileName,selDateFileName)) {
if (fStopOnError) return kFALSE;
}
AliSysInfo::AddStamp("ConvertRawFilesToDate");
if (deleteIntermediateFiles) {
AliRunLoader* runLoader = LoadRun("READ");
if (runLoader) for (Int_t iEvent = 0;
iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
char command[256];
snprintf(command, 256, "rm -r raw%d", iEvent);
gSystem->Exec(command);
}
delete runLoader;
}
if (rootOutput) {
if (!ConvertDateToRoot(dateFileName, fileName)) {
if (fStopOnError) return kFALSE;
}
AliSysInfo::AddStamp("ConvertDateToRoot");
if (deleteIntermediateFiles) {
gSystem->Unlink(dateFileName);
}
if (selrawdata) {
TString selFileName = "selected.";
selFileName += fileName;
if (!ConvertDateToRoot(selDateFileName, selFileName)) {
if (fStopOnError) return kFALSE;
}
if (deleteIntermediateFiles) {
gSystem->Unlink(selDateFileName);
}
}
}
}
return kTRUE;
}
Bool_t AliSimulation::WriteRawFiles(const char* detectors)
{
AliCodeTimerAuto("",0)
AliRunLoader* runLoader = LoadRun("READ");
if (!runLoader) return kFALSE;
for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
AliInfo(Form("processing event %d", iEvent));
runLoader->GetEvent(iEvent);
TString baseDir = gSystem->WorkingDirectory();
char dirName[256];
snprintf(dirName, 256, "raw%d", iEvent);
gSystem->MakeDirectory(dirName);
if (!gSystem->ChangeDirectory(dirName)) {
AliError(Form("couldn't change to directory %s", dirName));
if (fStopOnError) return kFALSE; else continue;
}
ofstream runNbFile(Form("run%u",runLoader->GetHeader()->GetRun()));
runNbFile.close();
TString detStr = detectors;
if (IsSelected("HLT", detStr)) {
}
TObjArray* detArray = runLoader->GetAliRun()->Detectors();
for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
AliModule* det = (AliModule*) detArray->At(iDet);
if (!det || !det->IsActive()) continue;
if (IsSelected(det->GetName(), detStr)) {
AliInfo(Form("creating raw data from digits for %s", det->GetName()));
det->Digits2Raw();
}
}
if (!WriteTriggerRawData())
if (fStopOnError) return kFALSE;
gSystem->ChangeDirectory(baseDir);
if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
AliError(Form("the following detectors were not found: %s",
detStr.Data()));
if (fStopOnError) return kFALSE;
}
}
delete runLoader;
return kTRUE;
}
Bool_t AliSimulation::ConvertRawFilesToDate(const char* dateFileName,
const char* selDateFileName)
{
AliCodeTimerAuto("",0)
char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
if (!path) {
AliError("the program dateStream was not found");
if (fStopOnError) return kFALSE;
} else {
delete[] path;
}
AliRunLoader* runLoader = LoadRun("READ");
if (!runLoader) return kFALSE;
AliInfo(Form("converting raw data DDL files to DATE file %s", dateFileName));
Bool_t selrawdata = kFALSE;
if (strcmp(selDateFileName,"") != 0) selrawdata = kTRUE;
char command[256];
snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d",
dateFileName, runLoader->GetNumberOfEvents(),runLoader->GetHeader()->GetRun());
FILE* pipe = gSystem->OpenPipe(command, "w");
if (!pipe) {
AliError(Form("Cannot execute command: %s",command));
return kFALSE;
}
Int_t selEvents = 0;
for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
UInt_t detectorPattern = 0;
runLoader->GetEvent(iEvent);
if (!runLoader->LoadTrigger()) {
AliCentralTrigger *aCTP = runLoader->GetTrigger();
detectorPattern = aCTP->GetClusterMask();
if (selrawdata) {
if (aCTP->GetClassMask()) selEvents++;
}
}
else {
AliWarning("No trigger can be loaded! Some fields in the event header will be empty !");
if (selrawdata) {
AliWarning("No trigger can be loaded! Writing of selected raw data is abandoned !");
selrawdata = kFALSE;
}
}
fprintf(pipe, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
Float_t ldc = 0;
Int_t prevLDC = -1;
for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
Int_t ldcID = Int_t(ldc + 0.0001);
ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
char rawFileName[256];
snprintf(rawFileName, 256, "raw%d/%s",
iEvent, AliDAQ::DdlFileName(iDet,iDDL));
FILE* file = fopen(rawFileName, "rb");
if (!file) continue;
fseek(file, 0, SEEK_END);
unsigned long size = ftell(file);
fclose(file);
if (!size) continue;
if (ldcID != prevLDC) {
fprintf(pipe, " LDC Id %d\n", ldcID);
prevLDC = ldcID;
}
fprintf(pipe, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
}
}
}
Int_t result = gSystem->ClosePipe(pipe);
if (!(selrawdata && selEvents > 0)) {
delete runLoader;
return (result == 0);
}
AliInfo(Form("converting selected by trigger cluster raw data DDL files to DATE file %s", selDateFileName));
snprintf(command, 256, "dateStream -c -s -D -o %s -# %d -C -run %d",
selDateFileName,selEvents,runLoader->GetHeader()->GetRun());
FILE* pipe2 = gSystem->OpenPipe(command, "w");
for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
UInt_t detectorPattern = 0;
TString detClust;
runLoader->GetEvent(iEvent);
if (!runLoader->LoadTrigger()) {
AliCentralTrigger *aCTP = runLoader->GetTrigger();
if (aCTP->GetClassMask() == 0) continue;
detectorPattern = aCTP->GetClusterMask();
detClust = AliDAQ::ListOfTriggeredDetectors(detectorPattern);
AliInfo(Form("List of detectors to be read out: %s",detClust.Data()));
}
fprintf(pipe2, "GDC DetectorPattern %u Timestamp %ld\n", detectorPattern, runLoader->GetHeader()->GetTimeStamp());
Float_t ldc = 0;
Int_t prevLDC = -1;
for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
if (!IsSelected(AliDAQ::DetectorName(iDet),detClust)) continue;
for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(iDet); iDDL++) {
Int_t ddlID = AliDAQ::DdlID(iDet,iDDL);
Int_t ldcID = Int_t(ldc + 0.0001);
ldc += AliDAQ::NumberOfLdcs(iDet) / AliDAQ::NumberOfDdls(iDet);
char rawFileName[256];
snprintf(rawFileName, 256, "raw%d/%s",
iEvent, AliDAQ::DdlFileName(iDet,iDDL));
FILE* file = fopen(rawFileName, "rb");
if (!file) continue;
fseek(file, 0, SEEK_END);
unsigned long size = ftell(file);
fclose(file);
if (!size) continue;
if (ldcID != prevLDC) {
fprintf(pipe2, " LDC Id %d\n", ldcID);
prevLDC = ldcID;
}
fprintf(pipe2, " Equipment Id %d Payload %s\n", ddlID, rawFileName);
}
}
}
Int_t result2 = gSystem->ClosePipe(pipe2);
delete runLoader;
return ((result == 0) && (result2 == 0));
}
Bool_t AliSimulation::ConvertDateToRoot(const char* dateFileName,
const char* rootFileName)
{
const Int_t kDBSize = 2000000000;
const Int_t kTagDBSize = 1000000000;
const Bool_t kFilter = kFALSE;
const Int_t kCompression = 1;
char* path = gSystem->Which(gSystem->Getenv("PATH"), "alimdc");
if (!path) {
AliError("the program alimdc was not found");
if (fStopOnError) return kFALSE;
} else {
delete[] path;
}
AliInfo(Form("converting DATE file %s to root file %s",
dateFileName, rootFileName));
const char* rawDBFS[2] = { "/tmp/mdc1", "/tmp/mdc2" };
const char* tagDBFS = "/tmp/mdc1/tags";
if (gSystem->Getenv("ALIMDC_RAWDB1"))
rawDBFS[0] = gSystem->Getenv("ALIMDC_RAWDB1");
if (gSystem->Getenv("ALIMDC_RAWDB2"))
rawDBFS[1] = gSystem->Getenv("ALIMDC_RAWDB2");
if (gSystem->Getenv("ALIMDC_TAGDB"))
tagDBFS = gSystem->Getenv("ALIMDC_TAGDB");
gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
gSystem->Exec(Form("rm -rf %s",tagDBFS));
gSystem->Exec(Form("mkdir %s",rawDBFS[0]));
gSystem->Exec(Form("mkdir %s",rawDBFS[1]));
gSystem->Exec(Form("mkdir %s",tagDBFS));
Int_t result = gSystem->Exec(Form("alimdc %d %d %d %d %s",
kDBSize, kTagDBSize, kFilter, kCompression, dateFileName));
gSystem->Exec(Form("mv %s/*.root %s", rawDBFS[0], rootFileName));
gSystem->Exec(Form("rm -rf %s",rawDBFS[0]));
gSystem->Exec(Form("rm -rf %s",rawDBFS[1]));
gSystem->Exec(Form("rm -rf %s",tagDBFS));
return (result == 0);
}
AliRunLoader* AliSimulation::LoadRun(const char* mode) const
{
delete AliRunLoader::Instance();
AliRunLoader* runLoader =
AliRunLoader::Open(fGAliceFileName.Data(),
AliConfig::GetDefaultEventFolderName(), mode);
if (!runLoader) {
AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
return NULL;
}
runLoader->LoadgAlice();
runLoader->LoadHeader();
gAlice = runLoader->GetAliRun();
if (!gAlice) {
AliError(Form("no gAlice object found in file %s",
fGAliceFileName.Data()));
return NULL;
}
return runLoader;
}
Int_t AliSimulation::GetNSignalPerBkgrd(Int_t nEvents) const
{
if (!fBkgrdFileNames) return 1;
Int_t nBkgrdFiles = fBkgrdFileNames->GetEntriesFast();
if (nBkgrdFiles == 0) return 1;
if (nEvents <= 0) {
AliRunLoader* runLoader =
AliRunLoader::Open(fGAliceFileName.Data(), "SIGNAL");
if (!runLoader) return 1;
nEvents = runLoader->GetNumberOfEvents();
delete runLoader;
}
Int_t result = 0;
for (Int_t iBkgrdFile = 0; iBkgrdFile < nBkgrdFiles; iBkgrdFile++) {
const char* fileName = ((TObjString*)
(fBkgrdFileNames->At(iBkgrdFile)))->GetName();
AliRunLoader* runLoader =
AliRunLoader::Open(fileName, "BKGRD");
if (!runLoader) continue;
Int_t nBkgrdEvents = runLoader->GetNumberOfEvents();
delete runLoader;
Int_t nSignalPerBkgrd = fBkgrdFileNames->At(iBkgrdFile)->GetUniqueID();
if (nSignalPerBkgrd <= 0) {
nSignalPerBkgrd = (nEvents-1) / nBkgrdEvents + 1;
} else if (result && (result != nSignalPerBkgrd)) {
AliInfo(Form("the number of signal events per background event "
"will be changed from %d to %d for stream %d",
nSignalPerBkgrd, result, iBkgrdFile+1));
nSignalPerBkgrd = result;
}
if (!result) result = nSignalPerBkgrd;
if (nSignalPerBkgrd * nBkgrdEvents < nEvents) {
AliWarning(Form("not enough background events (%d) for %d signal events "
"using %d signal per background events for stream %d",
nBkgrdEvents, nEvents, nSignalPerBkgrd, iBkgrdFile+1));
}
}
return result;
}
Bool_t AliSimulation::IsSelected(TString detName, TString& detectors) const
{
if ((detectors.CompareTo("ALL") == 0) ||
detectors.BeginsWith("ALL ") ||
detectors.EndsWith(" ALL") ||
detectors.Contains(" ALL ")) {
detectors = "ALL";
return kTRUE;
}
Bool_t result = kFALSE;
if ((detectors.CompareTo(detName) == 0) ||
detectors.BeginsWith(detName+" ") ||
detectors.EndsWith(" "+detName) ||
detectors.Contains(" "+detName+" ")) {
detectors.ReplaceAll(detName, "");
result = kTRUE;
}
while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
return result;
}
Int_t AliSimulation::ConvertRaw2SDigits(const char* rawDirectory, const char* esdFileName, Int_t N, Int_t nSkip)
{
if (!gAlice) {
AliError("no gAlice object. Restart aliroot and try again.");
return kFALSE;
}
if (gAlice->Modules()->GetEntries() > 0) {
AliError("gAlice was already run. Restart aliroot and try again.");
return kFALSE;
}
AliInfo(Form("initializing gAlice with config file %s",fConfigFileName.Data()));
gAlice->Announce();
gROOT->LoadMacro(fConfigFileName.Data());
gInterpreter->ProcessLine(gAlice->GetConfigFunction());
if(AliCDBManager::Instance()->GetRun() >= 0) {
SetRunNumber(AliCDBManager::Instance()->GetRun());
} else {
AliWarning("Run number not initialized!!");
}
AliRunLoader::Instance()->CdGAFile();
AliPDG::AddParticlesToPdgDataBase();
TVirtualMC::GetMC()->SetMagField(TGeoGlobalMagField::Instance()->GetField());
gAlice->GetMCApp()->Init();
gAlice->InitLoaders();
AliRunLoader::Instance()->MakeTree("E");
AliRunLoader::Instance()->LoadKinematics("RECREATE");
AliRunLoader::Instance()->LoadTrackRefs("RECREATE");
AliRunLoader::Instance()->LoadHits("all","RECREATE");
AliRunLoader::Instance()->CdGAFile();
gAlice->Write();
InitCDB();
Int_t iDet;
AliRunLoader* runLoader = AliRunLoader::Instance();
TFile* esdFile = 0;
TTree* treeESD = 0;
AliESDEvent* esd = 0;
if (esdFileName && (strlen(esdFileName)>0)) {
esdFile = TFile::Open(esdFileName);
if (esdFile) {
esd = new AliESDEvent();
esdFile->GetObject("esdTree", treeESD);
if (treeESD) {
esd->ReadFromTree(treeESD);
if (nSkip>0) {
AliInfo(Form("Asking to skip first %d ESDs events",nSkip));
} else {
nSkip=0;
}
}
}
}
TString fileName(rawDirectory);
AliRawReader* rawReader = AliRawReader::Create(fileName.Data());
if (!rawReader) return (kFALSE);
TObjArray* detArray = runLoader->GetAliRun()->Detectors();
AliHeader* header = runLoader->GetHeader();
Int_t nev = 0;
while(kTRUE) {
if (!(rawReader->NextEvent())) break;
runLoader->SetEventNumber(nev);
runLoader->GetHeader()->Reset(rawReader->GetRunNumber(),
nev, nev);
runLoader->GetEvent(nev);
AliInfo(Form("We are at event %d",nev));
TString detStr = fMakeSDigits;
for (iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
AliModule* det = (AliModule*) detArray->At(iDet);
if (!det || !det->IsActive()) continue;
if (IsSelected(det->GetName(), detStr)) {
AliInfo(Form("Calling Raw2SDigits for %s", det->GetName()));
det->Raw2SDigits(rawReader);
rawReader->Reset();
}
}
if (treeESD) {
AliInfo(Form("Selected event %d correspond to event %d in raw and to %d in esd",nev,rawReader->GetEventIndex(),nSkip+rawReader->GetEventIndex()));
treeESD->GetEvent(nSkip+rawReader->GetEventIndex());
const AliESDVertex* esdVertex = esd->GetPrimaryVertex();
Double_t position[3];
esdVertex->GetXYZ(position);
AliGenEventHeader* mcHeader = new AliGenEventHeader("ESD");
TArrayF mcV;
mcV.Set(3);
for (Int_t i = 0; i < 3; i++) mcV[i] = position[i];
mcHeader->SetPrimaryVertex(mcV);
header->Reset(0,nev);
header->SetGenEventHeader(mcHeader);
AliInfo(Form("***** Saved vertex %f %f %f \n", position[0], position[1], position[2]));
}
runLoader->TreeE()->Fill();
AliInfo(Form("Finished event %d",nev));
nev++;
if (N>0&&nev>=N)
break;
}
delete rawReader;
runLoader->CdGAFile();
runLoader->WriteHeader("OVERWRITE");
runLoader->WriteRunLoader();
return nev;
}
void AliSimulation::FinishRun()
{
if(IsLegoRun())
{
AliDebug(1, "Finish Lego");
AliRunLoader::Instance()->CdGAFile();
fLego->FinishRun();
}
TIter next(gAlice->Modules());
AliModule *detector;
while((detector = dynamic_cast<AliModule*>(next()))) {
AliDebug(2, Form("%s->FinishRun()", detector->GetName()));
detector->FinishRun();
}
AliDebug(1, "AliRunLoader::Instance()->WriteHeader(OVERWRITE)");
AliRunLoader::Instance()->WriteHeader("OVERWRITE");
AliRunLoader::Instance()->CdGAFile();
gAlice->Write(0,TObject::kOverwrite);
AliRunLoader::Instance()->Write(0,TObject::kOverwrite);
if(gAlice->GetMCApp()) gAlice->GetMCApp()->FinishRun();
AliRunLoader::Instance()->Synchronize();
}
Int_t AliSimulation::GetDetIndex(const char* detector)
{
Int_t index = -1 ;
for (index = 0; index < fgkNDetectors ; index++) {
if ( strcmp(detector, fgkDetectorName[index]) == 0 )
break ;
}
return index ;
}
Bool_t AliSimulation::CreateHLT()
{
gSystem->Load(ALIHLTSIMULATION_LIBRARY);
AliHLTSimulationGetLibraryVersion_t fctVersion=(AliHLTSimulationGetLibraryVersion_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_GET_LIBRARY_VERSION));
if (!fctVersion) {
AliError(Form("can not load library %s", ALIHLTSIMULATION_LIBRARY));
return kFALSE;
}
if (fctVersion()!= ALIHLTSIMULATION_LIBRARY_VERSION) {
AliWarning(Form("%s version does not match: compiled for version %d, loaded %d", ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_LIBRARY_VERSION, fctVersion()));
}
typedef void (*CompileInfo)( const char*& date, const char*& time);
CompileInfo fctInfo=(CompileInfo)gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, "CompileInfo");
if (fctInfo) {
const char* date="";
const char* time="";
(*fctInfo)(date, time);
if (!date) date="unknown";
if (!time) time="unknown";
AliInfo(Form("%s build on %s (%s)", ALIHLTSIMULATION_LIBRARY, date, time));
} else {
AliInfo(Form("no build info available for %s", ALIHLTSIMULATION_LIBRARY));
}
AliHLTSimulationCreateInstance_t fctCreate=(AliHLTSimulationCreateInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_CREATE_INSTANCE));
if (fctCreate==NULL || (fpHLT=(fctCreate()))==NULL) {
AliError(Form("can not create instance of HLT simulation (creator %p)", fctCreate));
return kFALSE;
}
TString specObjects;
for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
if (specObjects.Length()>0) specObjects+=" ";
specObjects+=fSpecCDBUri[i]->GetName();
}
AliHLTSimulationSetup_t fctSetup=(AliHLTSimulationSetup_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_SETUP));
if (fctSetup==NULL || fctSetup(fpHLT, this, specObjects.Data())<0) {
AliWarning(Form("failed to setup HLT simulation (function %p)", fctSetup));
}
return kTRUE;
}
Bool_t AliSimulation::RunHLT()
{
int iResult=0;
if (!fpHLT && !CreateHLT()) {
return kFALSE;
}
AliHLTSimulation* pHLT=fpHLT;
AliRunLoader* pRunLoader = LoadRun("READ");
if (!pRunLoader) return kFALSE;
InitCDB();
if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
SetCDBLock();
TString options;
if (fRunHLT.CompareTo("default")!=0) options=fRunHLT;
TString detStr = fWriteRawData;
if (!IsSelected("HLT", detStr)) {
options+=" writerawfiles=";
} else {
options+=" writerawfiles=HLT";
}
if (!detStr.IsNull() && !options.Contains("rawfile=")) {
options+=" rawfile=./";
}
AliHLTSimulationInit_t fctInit=(AliHLTSimulationInit_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_INIT));
if (fctInit==NULL || (iResult=(fctInit(pHLT, pRunLoader, options.Data())))<0) {
AliError(Form("can not init HLT simulation: error %d (init %p)", iResult, fctInit));
} else {
AliHLTSimulationRun_t fctRun=(AliHLTSimulationRun_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_RUN));
if (fctRun==NULL || (iResult=(fctRun(pHLT, pRunLoader)))<0) {
AliError(Form("can not run HLT simulation: error %d (run %p)", iResult, fctRun));
}
}
AliHLTSimulationDeleteInstance_t fctDelete=(AliHLTSimulationDeleteInstance_t)(gSystem->DynFindSymbol(ALIHLTSIMULATION_LIBRARY, ALIHLTSIMULATION_DELETE_INSTANCE));
if (fctDelete==NULL || fctDelete(pHLT)<0) {
AliError(Form("can not delete instance of HLT simulation (creator %p)", fctDelete));
}
pHLT=NULL;
return iResult>=0?kTRUE:kFALSE;
}
Bool_t AliSimulation::RunQA()
{
AliQAManager::QAManager()->SetRunLoader(AliRunLoader::Instance()) ;
TString detectorsw("") ;
Bool_t rv = kTRUE ;
AliQAManager::QAManager()->SetEventSpecie(fEventSpecie) ;
detectorsw = AliQAManager::QAManager()->Run(fQADetectors.Data()) ;
if ( detectorsw.IsNull() )
rv = kFALSE ;
return rv ;
}
Bool_t AliSimulation::SetRunQA(TString detAndAction)
{
if (!detAndAction.Contains(":")) {
AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
fRunQA = kFALSE ;
return kFALSE ;
}
Int_t colon = detAndAction.Index(":") ;
fQADetectors = detAndAction(0, colon) ;
if (fQADetectors.Contains("ALL") ){
TString tmp = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
Int_t minus = fQADetectors.Last('-') ;
TString toKeep = Form("%s %s", fMakeDigits.Data(), fMakeDigitsFromHits.Data()) ;
TString toRemove("") ;
while (minus >= 0) {
toRemove = fQADetectors(minus+1, fQADetectors.Length()) ;
toRemove = toRemove.Strip() ;
toKeep.ReplaceAll(toRemove, "") ;
fQADetectors.ReplaceAll(Form("-%s", toRemove.Data()), "") ;
minus = fQADetectors.Last('-') ;
}
fQADetectors = toKeep ;
}
fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
if (fQATasks.Contains("ALL") ) {
fQATasks = Form("%d %d %d", AliQAv1::kHITS, AliQAv1::kSDIGITS, AliQAv1::kDIGITS) ;
} else {
fQATasks.ToUpper() ;
TString tempo("") ;
if ( fQATasks.Contains("HIT") )
tempo = Form("%d ", AliQAv1::kHITS) ;
if ( fQATasks.Contains("SDIGIT") )
tempo += Form("%d ", AliQAv1::kSDIGITS) ;
if ( fQATasks.Contains("DIGIT") )
tempo += Form("%d ", AliQAv1::kDIGITS) ;
fQATasks = tempo ;
if (fQATasks.IsNull()) {
AliInfo("No QA requested\n") ;
fRunQA = kFALSE ;
return kTRUE ;
}
}
TString tempo(fQATasks) ;
tempo.ReplaceAll(Form("%d", AliQAv1::kHITS), AliQAv1::GetTaskName(AliQAv1::kHITS)) ;
tempo.ReplaceAll(Form("%d", AliQAv1::kSDIGITS), AliQAv1::GetTaskName(AliQAv1::kSDIGITS)) ;
tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITS), AliQAv1::GetTaskName(AliQAv1::kDIGITS)) ;
AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
fRunQA = kTRUE ;
AliQAManager::QAManager()->SetActiveDetectors(fQADetectors) ;
AliQAManager::QAManager()->SetTasks(fQATasks) ;
for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++)
AliQAManager::QAManager()->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
return kTRUE;
}
void AliSimulation::ProcessEnvironmentVars()
{
AliInfo("Processing environment variables");
if (fSeed == 0) {
if (gSystem->Getenv("CONFIG_SEED")) {
fSeed = atoi(gSystem->Getenv("CONFIG_SEED"));
}
} else {
if (gSystem->Getenv("CONFIG_SEED")) {
AliInfo(Form("Seed for random number generation already set (%d)"
": CONFIG_SEED variable ignored!", fSeed));
}
}
AliInfo(Form("Seed for random number generation = %d ", fSeed));
if(fRun < 0) {
if (gSystem->Getenv("DC_RUN")) {
fRun = atoi(gSystem->Getenv("DC_RUN"));
}
} else {
if (gSystem->Getenv("DC_RUN")) {
AliInfo(Form("Run number already set (%d): DC_RUN variable ignored!", fRun));
}
}
AliInfo(Form("Run number = %d", fRun));
}
void AliSimulation::WriteGRPEntry()
{
AliInfo("Writing global run parameters entry into the OCDB");
AliGRPObject* grpObj = new AliGRPObject();
grpObj->SetRunType("PHYSICS");
grpObj->SetTimeStart(fTimeStart);
grpObj->SetTimeEnd(fTimeEnd);
grpObj->SetBeamEnergyIsSqrtSHalfGeV();
const AliGenerator *gen = gAlice->GetMCApp()->Generator();
Int_t a = 0;
Int_t z = 0;
if (gen) {
TString projectile;
gen->GetProjectile(projectile,a,z);
TString target;
gen->GetTarget(target,a,z);
TString beamType = projectile + "-" + target;
beamType.ReplaceAll(" ","");
if (!beamType.CompareTo("-")) {
grpObj->SetBeamType("UNKNOWN");
grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2);
}
else {
grpObj->SetBeamType(beamType);
if (z != 0) {
grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 * a / z);
} else {
grpObj->SetBeamEnergy(gen->GetEnergyCMS()/2 );
}
fEventSpecie = AliRecoParam::kHighMult;
if ((strcmp(beamType,"p-p") == 0) ||
(strcmp(beamType,"p-") == 0) ||
(strcmp(beamType,"-p") == 0) ||
(strcmp(beamType,"P-P") == 0) ||
(strcmp(beamType,"P-") == 0) ||
(strcmp(beamType,"-P") == 0)) {
fEventSpecie = AliRecoParam::kLowMult;
}
}
} else {
AliWarning("Unknown beam type and energy! Setting energy to 0");
grpObj->SetBeamEnergy(0);
grpObj->SetBeamType("UNKNOWN");
}
UInt_t detectorPattern = 0;
Int_t nDets = 0;
TObjArray *detArray = gAlice->Detectors();
for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors-1; iDet++) {
if (detArray->FindObject(AliDAQ::OfflineModuleName(iDet))) {
AliDebug(1, Form("Detector #%d found: %s", iDet, AliDAQ::OfflineModuleName(iDet)));
detectorPattern |= (1 << iDet);
nDets++;
}
}
if (!fTriggerConfig.IsNull())
detectorPattern |= (1 << AliDAQ::DetectorID("TRG"));
if (!fRunHLT.IsNull())
detectorPattern |= (1 << AliDAQ::kHLTId);
grpObj->SetNumberOfDetectors((Char_t)nDets);
grpObj->SetDetectorMask((Int_t)detectorPattern);
grpObj->SetLHCPeriod("LHC08c");
grpObj->SetLHCState("STABLE_BEAMS");
AliMagF *field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
Float_t solenoidField = field ? TMath::Abs(field->SolenoidField()) : 0;
Float_t factorSol = field ? field->GetFactorSol() : 0;
Float_t currentSol = TMath::Abs(factorSol)>1E-6 ?
TMath::Nint(TMath::Abs(solenoidField/factorSol))/5.*30000.*TMath::Abs(factorSol) : 0;
Float_t factorDip = field ? field->GetFactorDip() : 0;
Float_t currentDip = 6000.*TMath::Abs(factorDip);
grpObj->SetL3Current(currentSol,(AliGRPObject::Stats)0);
grpObj->SetDipoleCurrent(currentDip,(AliGRPObject::Stats)0);
grpObj->SetL3Polarity(factorSol>0 ? 0:1);
grpObj->SetDipolePolarity(factorDip>0 ? 0:1);
if (field) grpObj->SetUniformBMap(field->IsUniform());
grpObj->SetPolarityConventionLHC();
grpObj->SetCavernTemperature(0,(AliGRPObject::Stats)0);
AliCDBManager* man = AliCDBManager::Instance();
man->SetLock(0, fKey);
AliCDBStorage* sto = man->GetStorage(fGRPWriteLocation.Data());
AliCDBId id("GRP/GRP/Data", man->GetRun(), man->GetRun(), 1, 1);
AliCDBMetaData *metadata= new AliCDBMetaData();
metadata->SetResponsible("alice-off@cern.ch");
metadata->SetComment("Automatically produced GRP entry for Monte Carlo");
sto->Put(grpObj,id,metadata);
man->SetLock(1, fKey);
}
time_t AliSimulation::GenerateTimeStamp() const
{
if (fUseTimeStampFromCDB)
return fTimeStart + gRandom->Integer(fTimeEnd-fTimeStart);
else
return 0;
}
void AliSimulation::StoreUsedCDBMaps() const
{
AliRunLoader* runLoader = LoadRun();
if (!runLoader) {
AliError("Failed to open gAlice.root in write mode");
return;
}
const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
cdbMapCopy->SetOwner(1);
TIter iter(cdbMap->GetTable());
TPair* pair = 0;
while((pair = dynamic_cast<TPair*> (iter.Next()))){
TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
if (keyStr && valStr)
cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
}
TList *cdbListCopy = new TList();
cdbListCopy->SetOwner(1);
TIter iter2(cdbList);
AliCDBId* id=0;
while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
cdbListCopy->Add(new TObjString(id->ToString().Data()));
}
AliRunLoader::Instance()->CdGAFile();
gDirectory->WriteObject(cdbMapCopy,"cdbMap","kSingleKey");
gDirectory->WriteObject(cdbListCopy,"cdbList","kSingleKey");
delete runLoader;
AliInfo(Form("Stored used OCDB entries as TMap %s and TList %s in %s",
cdbMapCopy->GetName(),
cdbListCopy->GetName(),
fGAliceFileName.Data()));
}