//////////////////////////////////// // KRATTA DST reader // for the Asy-Eos experiment // TKratDstReader DEFINITION // Feb 2013 // revison 02/2013 // E.d.F ver 1.0 // sebastian.kupny@uj.edu.pl // Changes: //////////////////////////////////// #include "TKratDstReader.h" using std::cout; using std::endl; //_________________________TKratDstProgressBar__________________________ //______________________________________________________________________ TKratDstProgressBar::TKratDstProgressBar (int max, int len) : fLen(len), fMax(max), fMeter(new char[len+1]), fUpdateFrequency(10000) { Reset( max ); } //______________________________________________________________________ void TKratDstProgressBar::Reset (int max ) { fMax = max; memset (fMeter, ' ', fLen); fMeter[fLen] = 0; time(&fStart); } //______________________________________________________________________ bool TKratDstProgressBar::UpdateAndWrite(int iev){ iev++; if (iev % fUpdateFrequency == 0 || iev== fMax){ int proc = iev*100/fMax; double elapsed = difftime(time(0), fStart); time_t estimated = elapsed/iev*fMax - elapsed; struct tm * T = gmtime(&estimated); memset( fMeter,'#',proc*fLen/100); std::cout << "\rProgress: " << iev << "/" << fMax << " (" << proc << " %)" << " in " <tm_hour << ":" << std::setfill ('0') << std::setw (2) << T->tm_min << ":" << std::setfill ('0') << std::setw (2) << T->tm_sec << " " << " [" << fMeter << "]" << std::flush; if(iev==fMax){ std::cout << "\rCompleted"<< std::endl; return false; } } return true; } //___________________________TKratDstReader_____________________________ //______________________________________________________________________ TKratDstReader::TKratDstReader (const char* name, TString runlist, Int_t iVerbose ) : FairTask( name, iVerbose ), fSaveOutputToTree(0), fCurrentFileNumberOnTheList(0), fNumberOfEventsInAnalysedFiles(0), fEntryIndex(0) { fRunlist = runlist; fVerbose = iVerbose; fFileMainTreeName = "cbmsim"; fASYEventFromFileName = "KRATTA_EVENT_PARAM_CLONE"; fASYEventCopyName = "KRATTA_ASYEVENT_CLONE"; } //______________________________________________________________________ TKratDstReader::~TKratDstReader() { delete fProgressBar; } //______________________________________________________________________ void TKratDstReader::SetParContainers() { // Get run and runtime database FairRunAna* run = FairRunAna::Instance(); if ( ! run ) Fatal("SetParContainers", "No analysis run"); FairRuntimeDb* rtdb = run->GetRuntimeDb(); if ( ! rtdb ) Fatal("SetParContainers", "No runtime database"); } //______________________________________________________________________ Long64_t TKratDstReader::GetNumberOfEventsInInputFiles() { return fNumberOfEventsInAllFiles; } //______________________________________________________________________ Long64_t TKratDstReader::LoadFilesFromRunlist() { if(fRunlist.Length() == 0 ){ cout << "[TKratDstReader::Init():] Errror: File name with list of run not set or it is an empty string."<> fNFiles; Long64_t nNumberOfEventsInOneFile = 0; for( Int_t i = 0; i < fNFiles; i++) { fileMainTree = 0; fInputDstFileName = ""; nNumberOfEventsInOneFile = 0; inputFileWithRunList >> fInputDstFileName; if( fInputDstFileName.Length() > 0) { if ( CheckIfFileExist( fInputDstFileName.Data() )){ ///NOTE: If something is wrong one can remove bad file from the list. inputDstFile = new TFile( fInputDstFileName ); if ( inputDstFile->IsOpen() ){ fileMainTree = (TTree*) inputDstFile->Get( fFileMainTreeName.Data() ); if ( fileMainTree ){ nNumberOfEventsInOneFile = fileMainTree->GetEntries(); fListOfInputFiles.push_back( fInputDstFileName ); rRealNFiles++; fNumberOfEventsInAllFiles+=nNumberOfEventsInOneFile; if ( fVerbose > 1 ){ cout << "[TKratDstReader::Init():] " << fInputDstFileName << " with " << nNumberOfEventsInOneFile; cout << " entry/ies added to the list - [OK]" << endl; } }else{ if ( fVerbose > 1 ){ cout << "[TKratDstReader::Init():] Warning: in file " << fInputDstFileName; cout << " the tree " << fFileMainTreeName.Data() << " not found [BAD]" << endl; } } }else{ if ( fVerbose > 1 ) cout << "[TKratDstReader::Init():] Warning: " << fInputDstFileName << " couldn't open in root [BAD]" << endl; } inputDstFile-> Close(); delete inputDstFile; }else{ if ( fVerbose > 1 ) cout << "[TKratDstReader::Init():] Warning: " << fInputDstFileName << " not found [BAD]" << endl; } } } }else{ if ( fVerbose > 1 ) cout << "[TKratDstReader::Init():] Error: could not open the file: " << fRunlist << endl; } inputFileWithRunList.close(); if ( fNFiles != rRealNFiles ){ if ( fVerbose > 1 ){ cout << "[TKratDstReader::Init():] Warning: " << fRunlist << " contains only " << rRealNFiles; cout << " correct files what is less than " << fNFiles << " as you want! Changed " << endl; } fNFiles = rRealNFiles; } if ( fVerbose > 2 ){ cout << "[TKratDstReader::Init():] Files from run list " << fRunlist << " able to analyse:" << endl; for(Int_t i = 0; i< fListOfInputFiles.size(); i++) { cout << "[TKratDstReader::Init():] ["<< i << "] " << fListOfInputFiles[i] << endl; } } if ( fVerbose > 1 ) cout << "[TKratDstReader::Init():] Total number of events: " << fNumberOfEventsInAllFiles << endl; return fNumberOfEventsInAllFiles; } //______________________________________________________________________ InitStatus TKratDstReader::Init() { //fLogger->Info(MESSAGE_ORIGIN," TKratDstReader::Init()------------------Start "); //Farroot version: v-13.05 cout << "TKratDstReader::Init()------------------Start" << endl; /// Read input files if ( fVerbose > 0 ) cout <<"runlist= " << fRunlist << endl; LoadFilesFromRunlist(); if ( fNumberOfEventsInAllFiles < 1 ){ Fatal("Init", "No events to analyse"); } /// Configure Fairroot instances FairRootManager* ioman = FairRootManager::Instance(); if ( ! ioman ) Fatal("Init", "No FairRootManager"); /// Create an instance of progress bar (before function GetNextDstFile) fProgressBar = new TKratDstProgressBar( fNumberOfEventsInAllFiles ); ///Open first file to analyse if ( fListOfInputFiles.size() > 0 ){ GetNextDstFile(); } else{ if ( fVerbose > 0 )cout << "[TKratDstReader::Init:] Warning: No files found which can be analysd. Finish " << endl; } fASYEventCopy = new TClonesArray("ASYEvent"); fASYEventCopy = (TClonesArray*) fASYEventFromFile; ioman->Register( fASYEventCopyName, "Kratta ASYEvent Copy", fASYEventCopy, fSaveOutputToTree); cout << "TKratDstReader::Init()------------------End" << endl; return kSUCCESS; } //______________________________________________________________________ void TKratDstReader::Exec(Option_t* opt) { if ( fFileMainTree ){ if ( fEntryIndex < fNEventsInTree ){ fFileMainTree-> GetEntry( fEntryIndex ); if ( fVerbose > 0 ){ fProgressBar->UpdateAndWrite( fNumberOfEventsInAnalysedFiles + fEntryIndex ); } /// DOES IT ANYTHING IS NEEDED HERE, WHEN WE PUT THE fASYEvent STRUCTURE TO FAIRROOT MANAGER - NEED TO CHECK WITH NEXT TASK. /// Need to define copy constructor /// Debug level 6 or higher. First 1000 events and every 1000th event will be written if ( fVerbose > 5 && ((fNumberOfEventsInAnalysedFiles + fEntryIndex) % 1000 == 0 || (fNumberOfEventsInAnalysedFiles + fEntryIndex) < 1000 )) { cout << "Global Index=\t" << fNumberOfEventsInAnalysedFiles + fEntryIndex; cout << ", File Index=\t"<< fEntryIndex << "\t" << *fASYEventFromFile << endl; } //*fASYEvent = *fASYEventFromFile; /// IT IS ALSO POSSIBLE /* FOR DEBUGGING fKrattaAsyPeakClonesArray = fASYEventFromFile->peak; if ( fVerbose > 3 )fKrattaAsyPeakClonesArray->Print(); Int_t nPeaks = fKrattaAsyPeakClonesArray -> GetEntries(); if ( fVerbose > 3 )cout << fEntryIndex << "ORYG Run number = " << fASYEventFromFile->run << " " << fASYEventFromFile->evt << " peaks=" << nPeaks << endl; if (fASYEvent) delete fASYEvent; Int_t nPeaksNew = fKrattaAsyPeakClonesArray -> GetEntries(); if ( fVerbose > 3 )cout << fEntryIndex << "COPY Run number = " << fASYEvent->run << " " << fASYEvent->evt << " peaks=" << nPeaksNew << endl; if ( fVerbose > 3 ) for (Int_t pi = 0; pi < nPeaks; pi++ ) { ASYFadcPeak *krattaPt = (ASYFadcPeak*) (*fKrattaAsyPeakClonesArray)[pi]; cout << "\t" << pi << " " << krattaPt->s28 << " " << krattaPt->pds1 << " " << krattaPt->m11 << endl; } */ fEntryIndex++; }else{ ///cout << endl; fNumberOfEventsInAnalysedFiles += fEntryIndex; fEntryIndex = 0; if ( fVerbose > 1 ) cout << "[TKratDstReader::Exec:] No more entry in the tree. Analysed events: " << fNumberOfEventsInAnalysedFiles << ", reset index: " << fEntryIndex << endl; GetNextDstFile(); } }else{ if ( fVerbose > 1 ) cout << "[TKratDstReader::Exec:] No more entry in the tree " << endl; //Finish(); FairRunAna::Instance()->TerminateRun(); ///NOTE:Using that, the data will be wrong written to file with no possiblilty to open then. } } //______________________________________________________________________ void TKratDstReader::Reset() { // Clear the structure // cout << " -I- Digit Reset() called " << endl; // pc if (fLandDigi ) fLandDigi->Clear(); } //______________________________________________________________________ void TKratDstReader::Finish() { /// TODO: FREE ALL USED MEMORY, FREE POINTERS, ETC. // THIS PART OF CODE CAN BE LAUNCHED MANY TIMES /// OR BETTER PUT THAT TO DESTRUCTOR ?? ///cout << "FINISH" << endl; /// YOUR NEVER HERE RICH THIS POINT BECAUSE OF TerminateRun :) } //______________________________________________________________________ void TKratDstReader::GetNextDstFile() { ///Close the old one if ( fInputDstFile ){ if ( fInputDstFile -> IsOpen() ){ fInputDstFile-> Close(); delete fInputDstFile; fInputDstFile = 0; /// necessary fFileMainTree = 0; } } ///Open next file to analyse if( fCurrentFileNumberOnTheList < fListOfInputFiles.size() ) { fInputDstFileName = fListOfInputFiles[ fCurrentFileNumberOnTheList ]; fInputDstFile = new TFile( fInputDstFileName ); fCurrentFileNumberOnTheList++; fFileMainTree = (TTree*) fInputDstFile->Get( fFileMainTreeName.Data() ); fFileMainTree -> SetBranchAddress( fASYEventFromFileName.Data() , &fASYEventFromFile); fNEventsInTree = fFileMainTree->GetEntries(); if ( fVerbose > 1 ) cout << "[TKratDstReader::GetNextDstFile:] Number of entries in file." << fInputDstFileName << ": " << fNEventsInTree << endl; } else { if ( fVerbose > 1 ) cout << "[TKratDstReader::GetNextDstFile:] No more files." << endl; } } //______________________________________________________________________ bool TKratDstReader::CheckIfFileExist(const char *filename) { bool result = false; ifstream ifile(filename); if (ifile) { result = true; } ifile.close(); return result; } ClassImp( TKratDstProgressBar ) ClassImp( TKratDstReader )