CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MonitorTrackResiduals.cc
Go to the documentation of this file.
23 
25  : dqmStore_( edm::Service<DQMStore>().operator->() )
26  , conf_(iConfig), m_cacheID_(0)
27  , genTriggerEventFlag_(new GenericTriggerEventFlag(iConfig, consumesCollector(), *this))
28  , avalidator_(iConfig, consumesCollector()) {
29  ModOn = conf_.getParameter<bool>("Mod_On");
30 }
31 
34 }
35 
36 
38 }
39 
41 {
42  unsigned long long cacheID = iSetup.get<SiStripDetCablingRcd>().cacheIdentifier();
43  if (m_cacheID_ != cacheID) {
44  m_cacheID_ = cacheID;
45  this->createMEs( ibooker , iSetup);
46  }
47 }
48 
50 
51  // Initialize the GenericTriggerEventFlag
52  if ( genTriggerEventFlag_->on() ) genTriggerEventFlag_->initRun( run, iSetup );
53 }
54 
56 
57  //Retrieve tracker topology from geometry
59  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
60  const TrackerTopology* const tTopo = tTopoHandle.product();
61 
62  Parameters = conf_.getParameter<edm::ParameterSet>("TH1ResModules");
63  int32_t i_residuals_Nbins = Parameters.getParameter<int32_t>("Nbinx");
64  double d_residual_xmin = Parameters.getParameter<double>("xmin");
65  double d_residual_xmax = Parameters.getParameter<double>("xmax");
66  Parameters = conf_.getParameter<edm::ParameterSet>("TH1NormResModules");
67  int32_t i_normres_Nbins = Parameters.getParameter<int32_t>("Nbinx");
68  double d_normres_xmin = Parameters.getParameter<double>("xmin");
69  double d_normres_xmax = Parameters.getParameter<double>("xmax");
70 
71 
72  // use SistripHistoId for producing histogram id (and title)
73  SiStripHistoId hidmanager;
74  folder_organizer.setSiStripFolder(); // top SiStrip folder
75 
76  // take from eventSetup the SiStripDetCabling object
78  iSetup.get<SiStripDetCablingRcd>().get(tkmechstruct);
79 
80  // get list of active detectors from SiStripDetCabling
81  std::vector<uint32_t> activeDets;
82  activeDets.clear(); // just in case
83  tkmechstruct->addActiveDetectorsRawIds(activeDets);
84 
85  // use SiStripSubStructure for selecting certain regions
86  SiStripSubStructure substructure;
87  std::vector<uint32_t> DetIds = activeDets;
88 
89  // book histo per each detector module
90  for (std::vector<uint32_t>::const_iterator DetItr=activeDets.begin(),
91  DetItrEnd = activeDets.end(); DetItr!=DetItrEnd; ++DetItr)
92  {
93  uint ModuleID = (*DetItr);
94 
95  // is this a StripModule?
96  if( SiStripDetId(ModuleID).subDetector() != 0 ) {
97 
98  folder_organizer.setDetectorFolder(ModuleID, tTopo); // detid sets appropriate detector folder
99  // Book module histogramms?
100  if (ModOn) {
101  std::string hid = hidmanager.createHistoId("HitResiduals","det",ModuleID);
102  std::string normhid = hidmanager.createHistoId("NormalizedHitResiduals","det",ModuleID);
103  HitResidual[ModuleID] = ibooker.book1D(hid, hid,
104  i_residuals_Nbins,d_residual_xmin,d_residual_xmax);
105  HitResidual[ModuleID]->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
106  NormedHitResiduals[ModuleID] = ibooker.book1D(normhid, normhid,
107  i_normres_Nbins,d_normres_xmin,d_normres_xmax);
108  NormedHitResiduals[ModuleID]->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
109  }
110  // book layer level histogramms
111  std::pair<std::string,int32_t> subdetandlayer = folder_organizer.GetSubDetAndLayer(ModuleID, tTopo);
112  folder_organizer.setLayerFolder(ModuleID,tTopo,subdetandlayer.second);
113  if(! m_SubdetLayerResiduals[subdetandlayer ] ) {
114  // book histogramms on layer level, check for barrel for correct labeling
115  std::string histoname = Form(subdetandlayer.first.find("B") != std::string::npos ?
116  "HitResiduals_%s__Layer__%d" : "HitResiduals_%s__wheel__%d" ,
117  subdetandlayer.first.c_str(),std::abs(subdetandlayer.second));
118  std::string normhistoname =
119  Form(subdetandlayer.first.find("B") != std::string::npos ?
120  "NormalizedHitResidual_%s__Layer__%d" : "NormalizedHitResidual_%s__wheel__%d" ,
121  subdetandlayer.first.c_str(),std::abs(subdetandlayer.second));
122  m_SubdetLayerResiduals[subdetandlayer] =
123  ibooker.book1D(histoname.c_str(),histoname.c_str(),
124  i_residuals_Nbins,d_residual_xmin,d_residual_xmax);
125  m_SubdetLayerResiduals[subdetandlayer]->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
126  m_SubdetLayerNormedResiduals[subdetandlayer] =
127  ibooker.book1D(normhistoname.c_str(),normhistoname.c_str(),
128  i_normres_Nbins,d_normres_xmin,d_normres_xmax);
129  m_SubdetLayerNormedResiduals[subdetandlayer]->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
130  }
131  } // end 'is strip module'
132  } // end loop over activeDets
133 }
134 
135 
137 }
138 
139 
141 
142  //dqmStore_->showDirStructure();
143  bool outputMEsInRootFile = conf_.getParameter<bool>("OutputMEsInRootFile");
145  if(outputMEsInRootFile){
146  dqmStore_->save(outputFileName);
147  }
148 }
149 
150 
152 
153  // Filter out events if Trigger Filtering is requested
154  if (genTriggerEventFlag_->on()&& ! genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
155 
156  //Retrieve tracker topology from geometry
157  edm::ESHandle<TrackerTopology> tTopoHandle;
158  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
159  const TrackerTopology* const tTopo = tTopoHandle.product();
160 
161  std::vector<TrackerValidationVariables::AVHitStruct> v_hitstruct;
162  avalidator_.fillHitQuantities(iEvent,v_hitstruct);
163  for (std::vector<TrackerValidationVariables::AVHitStruct>::const_iterator it = v_hitstruct.begin(),
164  itEnd = v_hitstruct.end(); it != itEnd; ++it) {
165  uint RawId = it->rawDetId;
166 
167  // fill if hit belongs to StripDetector and its error is not zero
168  if( it->resXprimeErr != 0 && SiStripDetId(RawId).subDetector() != 0 ) {
169  if (ModOn && HitResidual[RawId]) {
170  HitResidual[RawId]->Fill(it->resXprime);
171  NormedHitResiduals[RawId]->Fill(it->resXprime/it->resXprimeErr);
172  }
173  std::pair<std::string, int32_t> subdetandlayer = folder_organizer.GetSubDetAndLayer(RawId, tTopo);
174  if(m_SubdetLayerResiduals[subdetandlayer]) {
175  m_SubdetLayerResiduals[subdetandlayer]->Fill(it->resXprime);
176  m_SubdetLayerNormedResiduals[subdetandlayer]->Fill(it->resXprime/it->resXprimeErr);
177  }
178  }
179  }
180 
181 }
182 
183 
184 
186 
T getParameter(std::string const &) const
virtual void analyze(const edm::Event &, const edm::EventSetup &)
void dqmBeginRun(const edm::Run &, const edm::EventSetup &)
unsigned long long m_cacheID_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::pair< std::string, int32_t > GetSubDetAndLayer(const uint32_t &detid, const TrackerTopology *tTopo, bool ring_flag=0)
void fillHitQuantities(const Trajectory *trajectory, std::vector< AVHitStruct > &v_avhitout)
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
void createMEs(DQMStore::IBooker &, const edm::EventSetup &)
SiStripFolderOrganizer folder_organizer
void setDetectorFolder(uint32_t rawdetid, const TrackerTopology *tTopo)
TrackerValidationVariables avalidator_
int iEvent
Definition: GenABIO.cc:230
vector< ParameterSet > Parameters
bool accept(const edm::Event &event, const edm::EventSetup &setup)
To be called from analyze/filter() methods.
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual void endRun(const edm::Run &, const edm::EventSetup &)
virtual void beginJob(void)
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
SubDetector subDetector() const
Definition: SiStripDetId.h:114
void setLayerFolder(uint32_t rawdetid, const TrackerTopology *tTopo, int32_t layer=0, bool ring_flag=0)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
std::map< std::pair< std::string, int32_t >, MonitorElement * > m_SubdetLayerNormedResiduals
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, const uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE", const bool resetMEsAfterWriting=false)
Definition: DQMStore.cc:2597
GenericTriggerEventFlag * genTriggerEventFlag_
std::map< std::pair< std::string, int32_t >, MonitorElement * > m_SubdetLayerResiduals
MonitorTrackResiduals(const edm::ParameterSet &)
void initRun(const edm::Run &run, const edm::EventSetup &setup)
To be called from beginRun() methods.
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
Definition: Run.h:41