CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DigiVertexCorrHistogramMaker.cc
Go to the documentation of this file.
7 #include "TH2F.h"
8 #include "TProfile.h"
9 
12 
13 
15  m_fhm(),
16  m_runHisto(false),
17  m_hitname(), m_nbins(500), m_scalefact(), m_binmax(), m_labels(), m_nmultvsnvtx(), m_nmultvsnvtxprof(), m_nmultvsnvtxvsbxprofrun(), m_subdirs() { }
18 
20  m_fhm(),
21  m_runHisto(iConfig.getUntrackedParameter<bool>("runHisto",false)),
22  m_hitname(iConfig.getUntrackedParameter<std::string>("hitName","digi")),
23  m_nbins(iConfig.getUntrackedParameter<int>("numberOfBins",500)),
24  m_scalefact(iConfig.getUntrackedParameter<int>("scaleFactor",5)),
25  m_labels(), m_nmultvsnvtx(), m_nmultvsnvtxprof(), m_subdirs()
26 {
27 
28  std::vector<edm::ParameterSet>
29  wantedsubds(iConfig.getUntrackedParameter<std::vector<edm::ParameterSet> >("wantedSubDets",std::vector<edm::ParameterSet>()));
30 
31  for(std::vector<edm::ParameterSet>::iterator ps=wantedsubds.begin();ps!=wantedsubds.end();++ps) {
32  m_labels[ps->getParameter<unsigned int>("detSelection")] = ps->getParameter<std::string>("detLabel");
33  m_binmax[ps->getParameter<unsigned int>("detSelection")] = ps->getParameter<int>("binMax");
34  }
35 
36 
37 }
38 
39 
41 
42  for(std::map<unsigned int,std::string>::const_iterator lab=m_labels.begin();lab!=m_labels.end();lab++) {
43 
44  const unsigned int i = lab->first; const std::string slab = lab->second;
45 
46  delete m_subdirs[i];
47  delete m_fhm[i];
48  }
49 
50 }
51 
52 
53 
54 void DigiVertexCorrHistogramMaker::book(const std::string dirname, const std::map<unsigned int, std::string>& labels) {
55 
56  m_labels = labels;
57  book(dirname);
58 
59 }
60 
62 
64  TFileDirectory subev = tfserv->mkdir(dirname);
65 
66  SiStripTKNumbers trnumb;
67 
68  edm::LogInfo("NumberOfBins") << "Number of Bins: " << m_nbins;
69  edm::LogInfo("ScaleFactors") << "y-axis range scale factor: " << m_scalefact;
70  edm::LogInfo("BinMaxValue") << "Setting bin max values";
71 
72  for(std::map<unsigned int,std::string>::const_iterator lab=m_labels.begin();lab!=m_labels.end();lab++) {
73 
74  const unsigned int i = lab->first; const std::string slab = lab->second;
75 
76  if(m_binmax.find(i)==m_binmax.end()) {
77  edm::LogVerbatim("NotConfiguredBinMax") << "Bin max for " << lab->second
78  << " not configured: " << trnumb.nstrips(i) << " used";
79  m_binmax[i] = trnumb.nstrips(i);
80  }
81 
82  edm::LogVerbatim("BinMaxValue") << "Bin max for " << lab->second << " is " << m_binmax[i];
83 
84  }
85 
86  for(std::map<unsigned int,std::string>::const_iterator lab=m_labels.begin();lab!=m_labels.end();++lab) {
87 
88  const int i = lab->first; const std::string slab = lab->second;
89 
90  char name[200];
91  char title[500];
92 
93  m_subdirs[i] = new TFileDirectory(subev.mkdir(slab.c_str()));
94  m_fhm[i] = new RunHistogramManager(true);
95 
96  if(m_subdirs[i]) {
97  sprintf(name,"n%sdigivsnvtx",slab.c_str());
98  sprintf(title,"%s %s multiplicity vs Nvtx",slab.c_str(),m_hitname.c_str());
99  m_nmultvsnvtx[i] = m_subdirs[i]->make<TH2F>(name,title,60,-0.5,59.5,m_nbins,0.,m_binmax[i]/(m_scalefact*m_nbins)*m_nbins);
100  m_nmultvsnvtx[i]->GetXaxis()->SetTitle("Number of Vertices"); m_nmultvsnvtx[i]->GetYaxis()->SetTitle("Number of Hits");
101  sprintf(name,"n%sdigivsnvtxprof",slab.c_str());
102  m_nmultvsnvtxprof[i] = m_subdirs[i]->make<TProfile>(name,title,60,-0.5,59.5);
103  m_nmultvsnvtxprof[i]->GetXaxis()->SetTitle("Number of Vertices"); m_nmultvsnvtxprof[i]->GetYaxis()->SetTitle("Number of Hits");
104 
105  if(m_runHisto) {
106  edm::LogInfo("RunHistos") << "Pseudo-booking run histos " << slab.c_str();
107  sprintf(name,"n%sdigivsnvtxvsbxprofrun",slab.c_str());
108  sprintf(title,"%s %s multiplicity vs Nvtx vs BX",slab.c_str(),m_hitname.c_str());
109  m_nmultvsnvtxvsbxprofrun[i] = m_fhm[i]->makeTProfile2D(name,title,3564,0.5,3563.5,60,-0.5,59.5);
110  }
111 
112  }
113 
114  }
115 
116 
117 }
118 
120 
122 
123 
124  for(std::map<unsigned int,std::string>::const_iterator lab=m_labels.begin();lab!=m_labels.end();++lab) {
125  const int i = lab->first; const std::string slab = lab->second;
126  m_fhm[i]->beginRun(iRun,*m_subdirs[i]);
127  if(m_runHisto) {
128  (*m_nmultvsnvtxvsbxprofrun[i])->GetXaxis()->SetTitle("BX"); (*m_nmultvsnvtxvsbxprofrun[i])->GetYaxis()->SetTitle("Nvertices");
129  }
130  }
131 
132 }
133 
134 void DigiVertexCorrHistogramMaker::fill(const edm::Event& iEvent, const unsigned int nvtx, const std::map<unsigned int,int>& ndigi) {
135 
136 
138 
139 
140  for(std::map<unsigned int,int>::const_iterator digi=ndigi.begin();digi!=ndigi.end();digi++) {
141 
142  if(m_labels.find(digi->first) != m_labels.end()) {
143 
144  const unsigned int i=digi->first;
145  m_nmultvsnvtx[i]->Fill(nvtx,digi->second);
146  m_nmultvsnvtxprof[i]->Fill(nvtx,digi->second);
147 
148  if(m_nmultvsnvtxvsbxprofrun[i] && *m_nmultvsnvtxvsbxprofrun[i]) (*m_nmultvsnvtxvsbxprofrun[i])->Fill(iEvent.bunchCrossing(),nvtx,digi->second);
149 
150  }
151 
152  }
153 }
154 
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::map< unsigned int, TProfile2D ** > m_nmultvsnvtxvsbxprofrun
std::map< unsigned int, std::string > m_labels
int bunchCrossing() const
Definition: EventBase.h:62
std::map< unsigned int, int > m_binmax
int nstrips(const SiStripDetId &detid) const
int iEvent
Definition: GenABIO.cc:243
std::map< unsigned int, RunHistogramManager * > m_fhm
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void book(const std::string dirname, const std::map< unsigned int, std::string > &labels)
std::map< unsigned int, TFileDirectory * > m_subdirs
std::map< unsigned int, TProfile * > m_nmultvsnvtxprof
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
void fill(const edm::Event &iEvent, const unsigned int nvtx, const std::map< unsigned int, int > &ndigi)
std::map< unsigned int, TH2F * > m_nmultvsnvtx
Definition: Run.h:33