CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HPDNoiseReader.cc
Go to the documentation of this file.
1 // --------------------------------------------------------
2 // Engine to read HPD noise events from the library
3 // Project: HPD noise library
4 // Author: F.Ratnikov UMd, Jan. 15, 2008
5 // $Id: HPDNoiseReader.cc,v 1.4 2008/08/04 22:07:08 fedor Exp $
6 // --------------------------------------------------------
7 
9 
10 #include "TFile.h"
11 #include "TTree.h"
12 #include "TBranch.h"
13 
14 #include <iostream>
15 
18 
20  mFile = new TFile (fFileName.c_str(), "READ");
21  mBuffer = new HPDNoiseData ();
23  mFile->GetObject (HPDNoiseDataCatalog::objectName(), catalog);
24  if (catalog) {
25  // initiate trees
26  const std::vector<std::string> names = catalog->allNames();
27  for (size_t i = 0; i < names.size(); ++i) {
28  TTree* tree = (TTree*) mFile->Get (names[i].c_str());
29  if (tree) {
30  mTrees.push_back (tree);
31  mDischargeRate.push_back (catalog->getDischargeRate (i));
34  mElectronEmissionRate.push_back(catalog->getElectronEmissionRate(i));
35  mCurrentEntry.push_back (0);
36  }
37  else {
38  std::cerr << "HPDNoiseReader::HPDNoiseReader-> Can not open tree " << names[i] << " in file " << fFileName << std::endl;
39  }
40  }
41  }
42  else {
43  std::cerr << "HPDNoiseReader::HPDNoiseReader-> Can not open catalog infile " << fFileName << std::endl;
44  }
45 }
46 
48  for (size_t i = 0; i < mTrees.size(); ++i) {
49  delete mTrees[i];
50  }
51  delete mFile;
52  delete mBuffer;
53 }
54 
55 
56 std::vector<std::string> HPDNoiseReader::allNames () const {
57  std::vector<std::string> result;
58  for (size_t i = 0; i < mTrees.size(); ++i) result.push_back (mTrees[i]->GetName ());
59  return result;
60 }
61 
63  for (size_t i = 0; i < mTrees.size(); ++i) {
64  if (std::string (mTrees[i]->GetName ()) == fName) return i;
65  }
66  return -1;
67 }
68 
69 float HPDNoiseReader::dischargeRate (Handle fHandle) const {
70  if (!valid (fHandle)) return 0;
71  return mDischargeRate[fHandle];
72 }
74  if (!valid (fHandle)) return 0;
75  return mIonFeedbackFirstPeakRate[fHandle];
76 }
78  if (!valid (fHandle)) return 0;
79  return mIonFeedbackSecondPeakRate[fHandle];
80 }
81 float HPDNoiseReader::emissionRate (Handle fHandle) const {
82  if (!valid (fHandle)) return 0;
83  return mElectronEmissionRate[fHandle];
84 }
85 
86 unsigned long HPDNoiseReader::totalEntries (Handle fHandle) const {
87  if (!valid (fHandle)) return 0;
88  return mTrees[fHandle]->GetEntries ();
89 }
90 
91 void HPDNoiseReader::grabEntry (Handle fHandle, unsigned long fEntry) {
92  if (!valid (fHandle)) return;
93  TBranch* branch = mTrees[fHandle]->GetBranch (HPDNoiseData::branchName());
94  branch->SetAddress (&mBuffer);
95  branch->GetEntry (fEntry);
96  mCurrentEntry [fHandle] = fEntry;
97 }
98 
99 void HPDNoiseReader::getEntry (Handle fHandle, HPDNoiseData** fData) {
100  if (!valid (fHandle)) return;
101  unsigned int entry = mCurrentEntry [fHandle];
102  if (++entry >= totalEntries (fHandle)) entry = 0; // roll over
103  grabEntry (fHandle, entry);
104  *fData = mBuffer;
105 }
106 
107 void HPDNoiseReader::getEntry (Handle fHandle, unsigned long fEntry, HPDNoiseData** fData) {
108  grabEntry (fHandle, fEntry);
109  *fData = mBuffer;
110 }
111 
int i
Definition: DBlmapReader.cc:9
static const HistoName names[]
float getElectronEmissionRate(size_t i) const
get thermal electron emission noise rate for the HPD instance
float dischargeRate(Handle fHandle) const
discharge rate for the instance
Handle getHandle(const std::string &fName)
get handle to access data for one HPD instance
std::vector< float > mElectronEmissionRate
unsigned long totalEntries(Handle fHandle) const
total number of noise events for the HPD instance in the library
static const char * objectName()
object name
float ionFeedbackSecondPeakRate(Handle fHandle) const
bool valid(Handle fHandle) const
test if handle is valid
float emissionRate(Handle fHandle) const
ithermal/field emission rate for the instance
std::vector< float > mIonFeedbackSecondPeakRate
HPDNoiseReader(const std::string &fFileName)
std::pair< std::string, MonitorElement * > entry
Definition: ME_MAP.h:8
tuple result
Definition: query.py:137
float getIonFeedbackSecondPeakRate(size_t i) const
float getIonFeedbackFirstPeakRate(size_t i) const
get ion feedback noise rate for the HPD instance
HPDNoiseData * mBuffer
const std::vector< std::string > & allNames() const
all HPD instance names
static const char * branchName()
branch name
Definition: HPDNoiseData.h:37
std::vector< float > mIonFeedbackFirstPeakRate
std::vector< std::string > allNames() const
all HPD instances in the library
std::pair< double, double > fEntry
Definition: CrossSection.h:40
float getDischargeRate(size_t i) const
get noise rate for the HPD instance
float ionFeedbackFirstPeakRate(Handle fHandle) const
ionfeedback rate for the instance
void getEntry(Handle fHandle, HPDNoiseData **fData)
retrive one entry from the sequentially
std::vector< TTree * > mTrees
std::vector< size_t > mCurrentEntry
void grabEntry(Handle fHandle, unsigned long fEntry)
std::vector< float > mDischargeRate