CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripPlotGain.cc
Go to the documentation of this file.
2 
7 
8 
9 
11  cacheID(0xFFFFFFFF)
12 {
13  //now do what ever initialization is needed
14  if(!edm::Service<SiStripDetInfoFileReader>().isAvailable()){
15  edm::LogError("TkLayerMap") <<
16  "\n------------------------------------------"
17  "\nUnAvailable Service SiStripDetInfoFileReader: please insert in the configuration file an instance like"
18  "\n\tprocess.SiStripDetInfoFileReader = cms.Service(\"SiStripDetInfoFileReader\")"
19  "\n------------------------------------------";
20  }
21 
23  file = new TFile("correlTest.root","RECREATE");
24  tkmap = new TrackerMap();
25 }
26 
27 
29 {}
30 
31 //
32 
33 void
35 
36  if(getCache(es)==cacheID )
37  return;
38  cacheID=getCache(es);
39 
40  edm::LogInfo("") << "[SiStripPlotGain::beginRun] cacheID " << cacheID << std::endl;
41 
42  es.get<SiStripApvGainRcd>().get(Handle_);
44 
45 
46 }
47 
48 void
50 
51  edm::LogInfo("") << "[Doanalysis]";
52 
53  std::vector<TH1F *>histos;
54 
57  iter=p.detid_begin;
58  iterE=p.detid_end;
59 
60  float value;
61 
62  //Divide result by d
63  for(;iter!=iterE;++iter){
64  getHistos(*iter,histos);
66 
67  edm::LogInfo("") << "[Doanalysis] detid " << *iter << " range " << range.second-range.first;
68  size_t apv=0, apvE= (range.second-range.first);
69  for (;apv<apvE;apv+=2){
70  value=gain.getApvGain(apv,range);
71  tkmap->fill(*iter,value);
72  for(size_t i=0;i<histos.size();++i)
73  histos[i]->Fill(value);
74  }
75 
76  }
77 }
78 
79 
80 void
81 SiStripPlotGain::getHistos(const uint32_t& detid,std::vector<TH1F*>& histos){
82 
83  histos.clear();
84 
85  int subdet=-999; int component=-999;
86  SiStripDetId a(detid);
87  if ( a.subdetId() == 3 ){
88  subdet=0;
89  component=TIBDetId(detid).layer();
90  } else if ( a.subdetId() == 4 ) {
91  subdet=1;
92  component=TIDDetId(detid).side()==2?TIDDetId(detid).wheel():TIDDetId(detid).wheel()+3;
93  } else if ( a.subdetId() == 5 ) {
94  subdet=2;
95  component=TOBDetId(detid).layer();
96  } else if ( a.subdetId() == 6 ) {
97  subdet=3;
98  component=TECDetId(detid).side()==2?TECDetId(detid).wheel():TECDetId(detid).wheel()+9;
99  }
100 
101  int index=100+subdet*100+component;
102 
103 
104  histos.push_back(getHisto(subdet+1));
105  histos.push_back(getHisto(index));
106 
107 }
108 
109 TH1F*
110 SiStripPlotGain::getHisto(const long unsigned int& index){
111  if(vTH1.size()<index+1)
112  vTH1.resize(index+1,0);
113 
114  if(vTH1[index]==0){
115  char name[128];
116  sprintf(name,"%lu",index);
117  edm::LogInfo("")<<"[getHisto] creating index " << index << std::endl;
118  vTH1[index]=new TH1F(name,name,150,0.,5.);
119  }
120 
121  return vTH1[index];
122 }
123 
124 void
126  for(size_t i=0;i<vTH1.size();i++)
127  if(vTH1[i]!=0)
128  vTH1[i]->Write();
129 
130  file->Write();
131  file->Close();
132 
133  tkmap->save(false,0,0,"testTkMap.png");
134 
135 }
136 
ContainerIterator getFirstElement(RegistryConstIterator &idet)
unsigned long long getCache(const edm::EventSetup &eSetup)
int i
Definition: DBlmapReader.cc:9
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
Registry::const_iterator RegistryConstIterator
RegistryConstIterator detid_end
RegistryConstIterator detid_begin
unsigned long long cacheID
unsigned int side() const
positive or negative id
Definition: TECDetId.h:47
virtual void endJob()
void getHistos(const uint32_t &detid, std::vector< TH1F * > &histos)
void save(bool print_total=true, float minval=0., float maxval=0., std::string s="svgmap.svg", int width=1500, int height=800)
Definition: TrackerMap.cc:632
std::pair< ContainerIterator, ContainerIterator > Range
void DoAnalysis(const SiStripApvGain &)
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
std::vector< TH1F * > vTH1
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
SiStripDetInfoFileReader * fr
dictionary histos
Definition: combine.py:3
ContainerIterator getLastElement(RegistryConstIterator &idet)
unsigned int side() const
positive or negative id
Definition: TIDDetId.h:45
const T & get() const
Definition: EventSetup.h:55
virtual void beginRun(const edm::Run &run, const edm::EventSetup &es)
T const * product() const
Definition: ESHandle.h:62
TH1F * getHisto(const long unsigned int &index)
TrackerMap * tkmap
unsigned int wheel() const
wheel id
Definition: TECDetId.h:52
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
double a
Definition: hdecay.h:121
SiStripPlotGain(const edm::ParameterSet &)
edm::ESHandle< SiStripApvGain > Handle_
RegistryPointers getRegistryPointers() const
float getApvGain(const uint16_t &apv, const Range &range) const
void fill(int layer, int ring, int nmod, float x)
Definition: TrackerMap.cc:2548
Definition: Run.h:31
unsigned int wheel() const
wheel id
Definition: TIDDetId.h:50