//*-- Author : Manuel Sanchez //*-- Modified : 08/10/2003 by Romain Holzmann //*-- Modified : 24/02/2003 by Ilse Koenig //*-- Modified : 23/05/2002 by Manuel Sanchez //*-- Modified : 01/03/2000 by Manuel Sanchez //*-- Modified : 5/05/98 by Manuel Sanchez //*-- Copyright : GENP (Univ. Santiago de Compostela) //_HADES_CLASS_DESCRIPTION //////////////////////////////////////////////////////////////////////////// // HGeantMergeSource // // This data source can read data from a ROOT file generated by // HGEANT. N root files can be read in parallel and the content // of HGeantKine + geant detector hit categories will be merged // into one event. If only 1 input file is used the source behaves // like a standard HRootSource. // // Methods : // I. addFile() // + n file addAdditionalInput() // II. addMultFiles(TString comma_separated_List) // // The first input file will create the event struture // and the content of the other inputs will be merged into // this event. Take the file with largest content first! //////////////////////////////////////////////////////////////////////////// #include "hgeantmergesource.h" #include "hades.h" #include "hrecevent.h" #include "heventheader.h" #include "hgeantheader.h" #include "hpartialevent.h" #include "hcategory.h" #include "hlinearcategory.h" #include "hmatrixcategory.h" #include "hgeantmaxtrk.h" #include "hgeantdef.h" #include #include #include #include #include "TKey.h" using namespace std; ClassImp(HGeantMergeSource) HGeantMergeSource::HGeantMergeSource(Bool_t fPers, Bool_t fMerg) { // for fPers==kFALSE, input cats are suppressed in output // for fMerg==kTRUE, partial events in input tree are merged into current // event, if it exists already at init(). fInput = 0; fEventInFile = 0; fCursor = 0; fCurrentRunId = INT_MIN; fCurrentRefId = -1; fGlobalRefId = -1; fDirectory = "./"; fPersistency = fPers; fMerge = fMerg; fEventList = 0; fLastRunId = -1; overwriteVersion = kFALSE; replaceVersion = 0; } HGeantMergeSource::~HGeantMergeSource(void) { //Class destructor, it clears the data source. Clear(); } Bool_t HGeantMergeSource::createGeantEvent(HRecEvent* fCurrentEvent){ HPartialEvent* geantEvt = fCurrentEvent->getPartialEvent(catSimul); if(!geantEvt){ fCurrentEvent->addPartialEvent(catSimul,"Simul","Simulated event"); HPartialEvent* pKine = fCurrentEvent->getPartialEvent(catSimul); pKine->setSubHeader(new HGeantHeader); pKine->addCategory(catGeantKine,new HLinearCategory("HGeantKine",MAXTRKKINE) ); Int_t iniMdc[4]={6,4,7,MAXTRKMDC}; pKine->addCategory(catMdcGeantRaw,new HMatrixCategory("HGeantMdc",4,iniMdc,0.1) ); Int_t iniRich[2]={6,MAXCKOVRICH}; pKine->addCategory(catRichGeantRaw,new HMatrixCategory("HGeantRichPhoton",2,iniRich,0.1) ); Int_t ini1Rich[2]={6,MAXPARTRICH}; pKine->addCategory(catRichGeantRaw+1,new HMatrixCategory("HGeantRichDirect",2,ini1Rich,0.1) ); Int_t ini2Rich[2]={6,MAXMIRRRICH}; pKine->addCategory(catRichGeantRaw+2,new HMatrixCategory("HGeantRichMirror",2,ini2Rich,0.1) ); Int_t iniRpc[2]={6,MAXTRKRPC}; pKine->addCategory(catRpcGeantRaw,new HMatrixCategory("HGeantRpc",2,iniRpc,0.5) ); Int_t iniTof[2]={6,MAXTRKTOF}; pKine->addCategory(catTofGeantRaw,new HMatrixCategory("HGeantTof",2,iniTof,0.1) ); Int_t iniShower[2]={6,MAXTRKSHOW}; pKine->addCategory(catShowerGeantRaw,new HMatrixCategory("HGeantShower",2,iniShower,0.1) ); pKine->addCategory(catWallGeantRaw,new HLinearCategory("HGeantWall",MAXTRKWALL) ); Int_t iniStart=MAXTRKSTART*0.1; pKine->addCategory(catStartGeantRaw,new HLinearCategory("HGeantStart",iniStart) ); Int_t iniEmc[2]={6,MAXTRKEMC}; pKine->addCategory(catEmcGeantRaw,new HMatrixCategory("HGeantEmc",2,iniEmc,0.1) ); pKine->addCategory(catStsGeantRaw,new HLinearCategory("HGeantSts",MAXTRKSTS) ); pKine->addCategory(catFRpcGeantRaw,new HLinearCategory("HGeantFRpc",MAXTRKFRPC) ); return kTRUE; } return kFALSE; } Bool_t HGeantMergeSource::init(void) { // Initializes the data source connecting the parts of the input tree in // the input ROOT file to the event structure used for the analysis. If // no event structure has been set by the user, HGeantMergeSource::init() will // set one read from the input file. If it exists, it will be used as such // (fMerge==kFALSE) or merged with the one in the input file (fMerge==kTRUE). Bool_t r = kTRUE; if (fEventInFile != 0) { //If there is an active event structure HRecEvent* fCurrentEvent = (HRecEvent*)gHades->getCurrentEvent(); if (fCurrentEvent != 0) { if (fMerge) { // merge the 2 event structures Warning("init","Merging with predefined event structure"); ((HRecEvent*)fEventInFile)->merge(fCurrentEvent); } else { Warning("init","Using predefined event structure"); } } else { gHades->setEvent(fEventInFile); } if ( fInput!=0) { Char_t sl=*(strchr(fInput->GetTitle(),'.')+1); switch (sl) { case '0' : fSplitLevel = 0; break; case '1' : fSplitLevel = 1; break; case '2' : fSplitLevel = 2; break; default : fSplitLevel = 0; } fInput->SetBranchStatus("*",kFALSE); gHades->activateTree(fInput); if (fEventList) { fEntries = fEventList->GetN(); } else { fEntries=fInput->GetEntries(); } if (fCursor >= fEntries) { Error("init","Entry not existing %i",fCursor); return kFALSE; //Entry not existing } if (fEventList) { fInput->GetEvent(fEventList->GetEntry(fCursor)); } else { fInput->GetEvent(fCursor); } fCurrentRunId = gHades->getCurrentEvent()->getHeader()->getEventRunNumber(); fLastRunId = fCurrentRunId; if(overwriteVersion) { // fix not set version in simulation gHades->getCurrentEvent()->getHeader()->setVersion(replaceVersion); } if (fRefIds.find(fCurrentRunId) != fRefIds.end()) { fCurrentRefId = fRefIds[fCurrentRunId]; } else { fCurrentRefId = fGlobalRefId; } if (fPersistency == kFALSE) { // set all input categories non-persistent for(Int_t i = 0; i < 16; i ++) { // loop over partial events HPartialEvent* fPar = ((HRecEvent*)fEventInFile)->getPartialEvent(i<setPersistency(kFALSE); } } } setCursorToPreviousEvent(); r = kTRUE; } else { Warning("init","Not input"); Clear(); r = kFALSE; } } else { r = kFALSE; } return r; } EDsState HGeantMergeSource::getNextEvent(Bool_t doUnpack) { //Retrieves next event in the input file. Int_t bread = 0; if (fInput) { if (fSplitLevel < 2) (*fEventAddr)->clearAll(fSplitLevel); if (fCursor < fEntries) { Int_t entry = fCursor; if (fEventList) { entry = fEventList->GetEntry(fCursor) ; bread = fInput->GetEvent(entry); } else { bread = fInput->GetEvent(fCursor); } if(bread > 0) { //------------------------------------------------ // read other inputs to merge GEANT events Bool_t status = kTRUE; for(UInt_t i=0;igetEntry(entry) < 0 ) { status = kFALSE; break; } evt->mergeGeantEvent((HRecEvent*)gHades->getCurrentEvent()); } if(status == kFALSE) bread = 0; } if (bread == 0) { return kDsEndData; } fCursor ++; fCurrentRunId = gHades->getCurrentEvent()->getHeader()->getEventRunNumber(); if(overwriteVersion) { // fix not set version in simulation gHades->getCurrentEvent()->getHeader()->setVersion(replaceVersion); } if (fCurrentRunId != fLastRunId) { if (fRefIds.find(fCurrentRunId) != fRefIds.end()) { fCurrentRefId = fRefIds[fCurrentRunId]; } else { fCurrentRefId = fGlobalRefId; } fLastRunId = fCurrentRunId; return kDsEndFile; } } else { return kDsEndData; } } else { return kDsError; } return kDsOk; } Bool_t HGeantMergeSource::getEvent(Int_t eventN) { //Retrieves event in position eventN in the input file, copying the //information to the event structure. if (fInput) { Bool_t status = kFALSE; Int_t entry = eventN ; if (fEventList) { entry = fEventList->GetEntry(eventN) ; } if(fInput->GetEvent(entry) > 0) { status = kTRUE; for(UInt_t i=0;igetEntry(entry) < 0 ) { status = kFALSE; break; } evt->mergeGeantEvent((HRecEvent*)gHades->getCurrentEvent()); } } return status; } return kFALSE; } void HGeantMergeSource::Clear(void) { //Closes the input file. if (fInput) { delete fInput; fInput = 0; for(UInt_t i=0;icloseInputFile(); delete evt; } } } Bool_t HGeantMergeSource::addFile(const Text_t *file,Bool_t print) { Text_t treeName[] = "T"; TString fname = getFileName(file); if (fileExists(fname)) { if (!fInput) { //If chain doesn't already exist TFile *fileTemp; TKey *key = 0; TString title; if(print)Info("addFile()","Add file %s",fname.Data()); fileTemp = getFile(fname.Data()); //Obtain tree title //new TFile(fname.Data()); key = fileTemp->GetKey(treeName); //No need to delete this pointer fEventInFile = (HEvent *)fileTemp->Get("Event"); fCurrentRunId = fEventInFile->getHeader()->getEventRunNumber(); if (fRefIds.find(fCurrentRunId) != fRefIds.end()) { fCurrentRefId = fRefIds[fCurrentRunId]; } else { fCurrentRefId = fGlobalRefId; } title = key->GetTitle(); fileTemp->Close(); delete fileTemp; fInput = new TChain(treeName,title.Data()); } fInput->Add(fname.Data()); } else { Warning("addFile","File %s not found",fname.Data()); return kFALSE; } return kTRUE; } Bool_t HGeantMergeSource::addMultFiles(TString commaSeparatedList) { // set multiple paralellel input files format : // "file1_with_path,file2_with_path,file3_with_path" // the first file will be set a standard input. the content // of the other files will copied into the event of file1. // for performance reasons the file with the biggest event // should be set as first file. commaSeparatedList.ReplaceAll(" ",""); if(commaSeparatedList==""){ Error("addMultFiles()","No input file specified!"); return kFALSE; } TObjArray* ar = commaSeparatedList.Tokenize(","); if(ar){ Bool_t ok = kFALSE; TString infile; if(ar->GetEntries()>0) { infile = ((TObjString*)ar->At(0))->GetString(); Info("addMultFiles()","Add file %s",infile.Data()); ok = addFile((const Text_t*) infile.Data(),kFALSE); } if(ok){ for(Int_t i = 1; i < ar->GetEntries();i++){ infile = ((TObjString*)ar->At(i))->GetString(); Info("addMultFiles()","Add file %s",infile.Data()); ok = addAdditionalInput(infile,kFALSE); if(!ok) break; } } delete ar; return ok; } return kFALSE; } Bool_t HGeantMergeSource::addAdditionalInput(TString filename,Bool_t print) { HParallelEvent* evt = new HParallelEvent(filename,filename); Bool_t ok = evt->setInputFile(filename,print); if(!ok) { delete evt; return kFALSE; } fAdditionalInputs.push_back(evt); return kTRUE; }