CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
OccupancyPlots.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: OccupancyPlots
4 // Class: OccupancyPlots
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Mon Oct 27 17:37:53 CET 2008
16 // $Id: OccupancyPlots.cc,v 1.3 2013/02/27 19:49:47 wmtan Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
25 
26 #include <vector>
27 #include <map>
28 #include <limits>
29 
30 #include "TProfile.h"
31 
42 
48 
53 
62 
63 //
64 // class decleration
65 //
66 
68  public:
69  explicit OccupancyPlots(const edm::ParameterSet&);
71 
72 
73 private:
74  virtual void beginJob() override ;
75  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
76  virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
77  virtual void endRun(const edm::Run&, const edm::EventSetup&) override;
78  virtual void endJob() override ;
79 
80  // ----------member data ---------------------------
81 
82  std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > > m_multiplicityMapTokens;
83  std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > > m_occupancyMapTokens;
85 
87  std::map<unsigned int,DetIdSelector> m_wantedsubdets;
88 
89  TProfile** m_avemultiplicity;
90  TProfile** m_aveoccupancy;
91 
94 
95  TProfile** m_averadius;
96  TProfile** m_avez;
97  TProfile** m_avex;
98  TProfile** m_avey;
99  TProfile** m_zavedr;
100  TProfile** m_zavedz;
101  TProfile** m_zavedrphi;
102  TProfile** m_yavedr;
103  TProfile** m_yavedz;
104  TProfile** m_yavedrphi;
105  TProfile** m_xavedr;
106  TProfile** m_xavedz;
107  TProfile** m_xavedrphi;
108 
109 
110 };
111 
112 //
113 // constants, enums and typedefs
114 //
115 
116 //
117 // static data member definitions
118 //
119 
120 //
121 // constructors and destructor
122 //
124  m_multiplicityMapTokens(edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag> >("multiplicityMaps"), [this](edm::InputTag const & tag){return consumes<std::map<unsigned int, int> >(tag);})),
125  m_occupancyMapTokens(edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag> >("occupancyMaps"), [this](edm::InputTag const & tag){return consumes<std::map<unsigned int, int> >(tag);})),
126  m_fp(iConfig.getUntrackedParameter<edm::FileInPath>("file",edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt"))),
127  m_rhm(consumesCollector()), m_wantedsubdets()
128 {
129  //now do what ever initialization is needed
130 
131  m_avemultiplicity = m_rhm.makeTProfile("avemult","Average Multiplicty",6000,0.5,6000.5);
132  m_aveoccupancy = m_rhm.makeTProfile("aveoccu","Average Occupancy",6000,0.5,6000.5);
133 
134  m_nchannels_ideal = m_rhm.makeTH1F("nchannels_ideal","Number of channels (ideal)",6000,0.5,6000.5);
135  m_nchannels_real = m_rhm.makeTH1F("nchannels_real","Number of channels (real)",6000,0.5,6000.5);
136 
137  m_averadius = m_rhm.makeTProfile("averadius","Average Module Radius",6000,0.5,6000.5);
138  m_avez = m_rhm.makeTProfile("avez","Average Module z coordinate",6000,0.5,6000.5);
139  m_avex = m_rhm.makeTProfile("avex","Average Module x coordinate",6000,0.5,6000.5);
140  m_avey = m_rhm.makeTProfile("avey","Average Module y coordinate",6000,0.5,6000.5);
141 
142  m_zavedr = m_rhm.makeTProfile("zavedr","Average z unit vector dr",6000,0.5,6000.5);
143  m_zavedz = m_rhm.makeTProfile("zavedz","Average z unit vector dz",6000,0.5,6000.5);
144  m_zavedrphi = m_rhm.makeTProfile("zavedrphi","Average z unit vector drphi",6000,0.5,6000.5);
145  m_xavedr = m_rhm.makeTProfile("xavedr","Average x unit vector dr",6000,0.5,6000.5);
146  m_xavedz = m_rhm.makeTProfile("xavedz","Average x unit vctor dz",6000,0.5,6000.5);
147  m_xavedrphi = m_rhm.makeTProfile("xavedrphi","Average Module x unit vector drphi",6000,0.5,6000.5);
148  m_yavedr = m_rhm.makeTProfile("yavedr","Average y unit vector dr",6000,0.5,6000.5);
149  m_yavedz = m_rhm.makeTProfile("yavedz","Average y unit vector dz",6000,0.5,6000.5);
150  m_yavedrphi = m_rhm.makeTProfile("yavedrphi","Average y unit vector drphi",6000,0.5,6000.5);
151 
152  std::vector<edm::ParameterSet> wantedsubdets_ps = iConfig.getParameter<std::vector<edm::ParameterSet> >("wantedSubDets");
153 
154  for(std::vector<edm::ParameterSet>::const_iterator wsdps = wantedsubdets_ps.begin();wsdps!=wantedsubdets_ps.end();++wsdps) {
155 
156  unsigned int detsel = wsdps->getParameter<unsigned int>("detSelection");
157  std::vector<std::string> selstr = wsdps->getUntrackedParameter<std::vector<std::string> >("selection");
158  m_wantedsubdets[detsel]=DetIdSelector(selstr);
159 
160  }
161 
162 
163 }
164 
165 
167 {
168 
169  // do anything here that needs to be done at desctruction time
170  // (e.g. close files, deallocate resources etc.)
171 
172 }
173 
174 
175 //
176 // member functions
177 //
178 
179 // ------------ method called to for each event ------------
180 void
182 {
183  using namespace edm;
184 
185 
186  for(std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > >::const_iterator mapToken = m_multiplicityMapTokens.begin();mapToken!=m_multiplicityMapTokens.end();++mapToken) {
187 
189  iEvent.getByToken(*mapToken,mults);
190 
191  for(std::map<unsigned int,int>::const_iterator mult=mults->begin();mult!=mults->end();mult++) {
192  if(m_avemultiplicity && *m_avemultiplicity) (*m_avemultiplicity)->Fill(mult->first,mult->second);
193  }
194  }
195 
196 
197 
198  for(std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > >::const_iterator mapToken = m_occupancyMapTokens.begin();mapToken!=m_occupancyMapTokens.end();++mapToken) {
199 
201  iEvent.getByToken(*mapToken,occus);
202 
203  for(std::map<unsigned int,int>::const_iterator occu=occus->begin();occu!=occus->end();occu++) {
204  if(m_aveoccupancy && *m_aveoccupancy) (*m_aveoccupancy)->Fill(occu->first,occu->second);
205  }
206  }
207 
208 }
209 
210 
211 // ------------ method called once each job just before starting event loop ------------
212 void
214 {
215 
216 }
217 
218 void
219 OccupancyPlots::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
220 
221  m_rhm.beginRun(iRun);
222 
223 }
224 
225 void
226 OccupancyPlots::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
227 
228 
229  // edm::ESHandle<GlobalTrackingGeometry> trkgeo;
230  // iSetup.get<GlobalTrackingGeometryRecord>().get("",trkgeo);
232  iSetup.get<TrackerDigiGeometryRecord>().get("",trkgeo);
233 
234  const Local2DPoint center(0.,0.);
235  const Local3DPoint locz(0.,0.,1.);
236  const Local3DPoint locx(1.,0.,0.);
237  const Local3DPoint locy(0.,1.,0.);
238  const GlobalPoint origin(0.,0.,0.);
239 
240  TrackingGeometry::DetIdContainer detunits = trkgeo->detUnitIds();
241 
242  for(TrackingGeometry::DetIdContainer::const_iterator det = detunits.begin(); det!=detunits.end(); ++det) {
243 
244  if(det->det()!=DetId::Tracker) continue;
245 
246  edm::LogInfo("DetIdFromGeometry") << det->rawId();
247 
248  GlobalPoint position = trkgeo->idToDet(*det)->toGlobal(center);
249  GlobalPoint zpos = trkgeo->idToDet(*det)->toGlobal(locz);
250  GlobalPoint xpos = trkgeo->idToDet(*det)->toGlobal(locx);
251  GlobalPoint ypos = trkgeo->idToDet(*det)->toGlobal(locy);
252  GlobalVector posvect = position - origin;
253  GlobalVector dz = zpos - position;
254  GlobalVector dx = xpos - position;
255  GlobalVector dy = ypos - position;
256 
257  double dzdr = posvect.perp()>0 ? (dz.x()*posvect.x()+dz.y()*posvect.y())/posvect.perp() : 0. ;
258  double dxdr = posvect.perp()>0 ? (dx.x()*posvect.x()+dx.y()*posvect.y())/posvect.perp() : 0. ;
259  double dydr = posvect.perp()>0 ? (dy.x()*posvect.x()+dy.y()*posvect.y())/posvect.perp() : 0. ;
260 
261  double dzdrphi = posvect.perp()>0 ? (dz.y()*posvect.x()-dz.x()*posvect.y())/posvect.perp() : 0. ;
262  double dxdrphi = posvect.perp()>0 ? (dx.y()*posvect.x()-dx.x()*posvect.y())/posvect.perp() : 0. ;
263  double dydrphi = posvect.perp()>0 ? (dy.y()*posvect.x()-dy.x()*posvect.y())/posvect.perp() : 0. ;
264 
265  for(std::map<unsigned int,DetIdSelector>::const_iterator sel=m_wantedsubdets.begin();sel!=m_wantedsubdets.end();++sel) {
266 
267  if(sel->second.isSelected(*det)) {
268  edm::LogInfo("SelectedDetId") << sel->first;
269  // average positions
270  if(m_averadius && *m_averadius) (*m_averadius)->Fill(sel->first,position.perp());
271  if(m_avez && *m_avez) (*m_avez)->Fill(sel->first,position.z());
272  if(m_avex && *m_avex) (*m_avex)->Fill(sel->first,position.x());
273  if(m_avey && *m_avey) (*m_avey)->Fill(sel->first,position.y());
274  if(m_zavedr && *m_zavedr) (*m_zavedr)->Fill(sel->first,dzdr);
275  if(m_zavedz && *m_zavedz) (*m_zavedz)->Fill(sel->first,dz.z());
276  if(m_zavedrphi && *m_zavedrphi) (*m_zavedrphi)->Fill(sel->first,dzdrphi);
277  if(m_xavedr && *m_xavedr) (*m_xavedr)->Fill(sel->first,dxdr);
278  if(m_xavedz && *m_xavedz) (*m_xavedz)->Fill(sel->first,dx.z());
279  if(m_xavedrphi && *m_xavedrphi) (*m_xavedrphi)->Fill(sel->first,dxdrphi);
280  if(m_yavedr && *m_yavedr) (*m_yavedr)->Fill(sel->first,dydr);
281  if(m_yavedz && *m_yavedz) (*m_yavedz)->Fill(sel->first,dy.z());
282  if(m_yavedrphi && *m_yavedrphi) (*m_yavedrphi)->Fill(sel->first,dydrphi);
283  }
284  }
285  }
286 
287 
289  iSetup.get<SiStripQualityRcd>().get("",quality);
290 
291 
293 
294  const std::vector<uint32_t>& detids = reader->getAllDetIds();
295 
296  for(std::vector<uint32_t>::const_iterator detid=detids.begin();detid!=detids.end();++detid) {
297 
298  int nchannideal = reader->getNumberOfApvsAndStripLength(*detid).first*128;
299  // int nchannreal = reader->getNumberOfApvsAndStripLength(*detid).first*128;
300  int nchannreal = 0;
301  for(int strip = 0; strip < nchannideal; ++strip) {
302  if(!quality->IsStripBad(*detid,strip)) ++nchannreal;
303  }
304 
305 
306  for(std::map<unsigned int,DetIdSelector>::const_iterator sel=m_wantedsubdets.begin();sel!=m_wantedsubdets.end();++sel) {
307 
308  if(sel->second.isSelected(*detid)) {
309  if(m_nchannels_ideal && *m_nchannels_ideal) (*m_nchannels_ideal)->Fill(sel->first,nchannideal);
310  if(m_nchannels_real && *m_nchannels_real) (*m_nchannels_real)->Fill(sel->first,nchannreal);
311  }
312  }
313 
314  }
315 
316 
318  iSetup.get<SiPixelQualityRcd>().get("",pxlquality);
319 
320 
322 
323  const std::vector<uint32_t>& pxldetids = pxlreader.getAllDetIds();
324 
325  for(std::vector<uint32_t>::const_iterator detid=pxldetids.begin();detid!=pxldetids.end();++detid) {
326 
327  int nchannideal = pxlreader.getDetUnitDimensions(*detid).first*pxlreader.getDetUnitDimensions(*detid).second;
328  int nchannreal = 0;
329  if(!pxlquality->IsModuleBad(*detid)) {
330  nchannreal = pxlreader.getDetUnitDimensions(*detid).first*pxlreader.getDetUnitDimensions(*detid).second;
331  }
332  /*
333  int nchannreal = 0;
334  for(int strip = 0; strip < nchannideal; ++strip) {
335  if(!quality->IsStripBad(*detid,strip)) ++nchannreal;
336  }
337  */
338 
339  for(std::map<unsigned int,DetIdSelector>::const_iterator sel=m_wantedsubdets.begin();sel!=m_wantedsubdets.end();++sel) {
340 
341  if(sel->second.isSelected(*detid)) {
342  if(m_nchannels_ideal && *m_nchannels_ideal) (*m_nchannels_ideal)->Fill(sel->first,nchannideal);
343  if(m_nchannels_real && *m_nchannels_real) (*m_nchannels_real)->Fill(sel->first,nchannreal);
344  }
345  }
346 
347  }
348 
349 }
350 // ------------ method called once each job just after ending the event loop ------------
351 void
353 }
354 //define this as a plug-in
OccupancyPlots(const edm::ParameterSet &)
TProfile ** m_zavedrphi
TProfile ** m_xavedz
T perp() const
Definition: PV3DBase.h:72
const std::vector< uint32_t > & getAllDetIds() const
TH1F ** m_nchannels_ideal
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void beginRun(const edm::Run &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
T y() const
Definition: PV3DBase.h:63
TProfile ** m_zavedr
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
virtual void beginJob() override
RunHistogramManager m_rhm
TH1F ** m_nchannels_real
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_multiplicityMapTokens
std::map< unsigned int, DetIdSelector > m_wantedsubdets
m_fp(iConfig.getUntrackedParameter< edm::FileInPath >("file", edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt")))
int iEvent
Definition: GenABIO.cc:230
TProfile ** m_avex
const std::vector< uint32_t > & getAllDetIds() const
TProfile ** m_aveoccupancy
TProfile ** m_averadius
m_rhm(consumesCollector())
T z() const
Definition: PV3DBase.h:64
TProfile ** m_avez
TProfile ** m_xavedr
TProfile ** m_avey
TProfile ** m_yavedr
m_occupancyMapTokens(edm::vector_transform(iConfig.getParameter< std::vector< edm::InputTag > >("occupancyMaps"), [this](edm::InputTag const &tag){return consumes< std::map< unsigned int, int > >(tag);}))
TProfile ** m_yavedz
TProfile ** m_avemultiplicity
m_wantedsubdets()
TProfile ** m_zavedz
const T & get() const
Definition: EventSetup.h:55
void beginRun(const edm::Run &iRun)
string const
Definition: compareJSON.py:14
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_occupancyMapTokens
TProfile ** m_yavedrphi
edm::FileInPath m_fp
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::string fullPath() const
Definition: FileInPath.cc:165
T x() const
Definition: PV3DBase.h:62
virtual void endJob() override
virtual void endRun(const edm::Run &, const edm::EventSetup &) override
Definition: Run.h:41
TProfile ** m_xavedrphi
std::vector< DetId > DetIdContainer
virtual void analyze(const edm::Event &, const edm::EventSetup &) override