CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTtTrigDBValidation.cc
Go to the documentation of this file.
2 
3 // Framework
8 
12 
13 // Geometry
17 
18 // DataFormats
20 
21 // tTrig
24 
25 #include <stdio.h>
26 #include <sstream>
27 #include <math.h>
28 #include "TFile.h"
29 
30 using namespace edm;
31 using namespace std;
32 
34  metname_("TTrigDBValidation"),
35  labelDBRef_(pset.getParameter<string>("labelDBRef")),
36  labelDB_(pset.getParameter<string>("labelDB")),
37  lowerLimit_(pset.getUntrackedParameter<int>("lowerLimit",1)),
38  higherLimit_(pset.getUntrackedParameter<int>("higherLimit",3))
39  {
40 
41  LogVerbatim(metname_) << "[DTtTrigDBValidation] Constructor called!";
42 }
43 
45 
47  edm::Run const &, edm::EventSetup const &setup) {
48 
49  LogVerbatim(metname_) << "[DTtTrigDBValidation] Parameters initialization";
50  iBooker.setCurrentFolder("DT/DtCalib/TTrigDBValidation");
51 
52  ESHandle<DTTtrig> tTrig_Ref;
53  setup.get<DTTtrigRcd>().get(labelDBRef_, tTrig_Ref);
54  const DTTtrig* DTTtrigRefMap = &*tTrig_Ref;
55  LogVerbatim(metname_) << "[DTtTrigDBValidation] reference Ttrig version: " << tTrig_Ref->version();
56 
57  ESHandle<DTTtrig> tTrig;
58  setup.get<DTTtrigRcd>().get(labelDB_, tTrig);
59  const DTTtrig* DTTtrigMap = &*tTrig;
60  LogVerbatim(metname_) << "[DTtTrigDBValidation] Ttrig to validate version: " << tTrig->version();
61 
62  //book&reset the summary histos
63  for(int wheel=-2; wheel<=2; wheel++){
64  bookHistos(iBooker, wheel);
65  tTrigDiffWheel_[wheel]->Reset();
66  }
67 
68  // Get the geometry
69  setup.get<MuonGeometryRecord>().get(dtGeom_);
70 
71  // Loop over Ref DB entries
72  for(DTTtrig::const_iterator it = DTTtrigRefMap->begin();
73  it != DTTtrigRefMap->end(); ++it) {
74  DTSuperLayerId slId((*it).first.wheelId,
75  (*it).first.stationId,
76  (*it).first.sectorId,
77  (*it).first.slId);
78  float tTrigMean;
79  float tTrigRms;
80  float kFactor;
81  DTTtrigRefMap->get(slId, tTrigMean, tTrigRms, kFactor, DTTimeUnits::ns);
82  float tTrigCorr = tTrigMean + kFactor*tTrigRms;
83  LogTrace(metname_)<< "Ref Superlayer: " << slId << "\n"
84  << " Ttrig mean (ns): " << tTrigMean
85  << " Ttrig rms (ns): " << tTrigRms
86  << " Ttrig k-Factor: " << kFactor
87  << " Ttrig value (ns): " << tTrigCorr;
88 
89  //tTrigRefMap[slId] = std::pair<float,float>(tTrigmean,tTrigrms);
90  tTrigRefMap_[slId] = pair<float,float>(tTrigCorr,tTrigRms);
91  }
92 
93  // Loop over Ref DB entries
94  for(DTTtrig::const_iterator it = DTTtrigMap->begin();
95  it != DTTtrigMap->end(); ++it) {
96  DTSuperLayerId slId((*it).first.wheelId,
97  (*it).first.stationId,
98  (*it).first.sectorId,
99  (*it).first.slId);
100  float tTrigMean;
101  float tTrigRms;
102  float kFactor;
103  DTTtrigMap->get(slId, tTrigMean, tTrigRms, kFactor, DTTimeUnits::ns);
104  float tTrigCorr = tTrigMean + kFactor*tTrigRms;
105  LogTrace(metname_)<< "Superlayer: " << slId << "\n"
106  << " Ttrig mean (ns): " << tTrigMean
107  << " Ttrig rms (ns): " << tTrigRms
108  << " Ttrig k-Factor: " << kFactor
109  << " Ttrig value (ns): " << tTrigCorr;
110 
111  //tTrigMap[slId] = std::pair<float,float>(tTrigmean,tTrigrms);
112  tTrigMap_[slId] = pair<float,float>(tTrigCorr,tTrigRms);
113  }
114 
115  for(map<DTSuperLayerId, pair<float,float> >::const_iterator it = tTrigRefMap_.begin();
116  it != tTrigRefMap_.end(); ++it) {
117  if(tTrigMap_.find((*it).first) == tTrigMap_.end()) continue;
118 
119  // compute the difference
120  float difference = tTrigMap_[(*it).first].first - (*it).second.first;
121 
122  //book histo
123  int wheel = (*it).first.chamberId().wheel();
124  int sector = (*it).first.chamberId().sector();
125  if(tTrigDiffHistos_.find(make_pair(wheel,sector)) == tTrigDiffHistos_.end()) bookHistos(iBooker, wheel, sector);
126 
127  LogTrace(metname_) << "Filling histos for super-layer: " << (*it).first << " difference: " << difference;
128 
129  // Fill the test histos
130  int entry = -1;
131  int station = (*it).first.chamberId().station();
132  if(station == 1) entry=0;
133  if(station == 2) entry=3;
134  if(station == 3) entry=6;
135  if(station == 4) entry=9;
136 
137  int slBin = entry + (*it).first.superLayer();
138  if(slBin == 12) slBin=11;
139 
140  tTrigDiffHistos_[make_pair(wheel,sector)]->setBinContent(slBin, difference);
141  if(abs(difference) < lowerLimit_){
142  tTrigDiffWheel_[wheel]->setBinContent(slBin,sector,1);
143  }else if(abs(difference) < higherLimit_){
144  tTrigDiffWheel_[wheel]->setBinContent(slBin,sector,2);
145  }else{
146  tTrigDiffWheel_[wheel]->setBinContent(slBin,sector,3);
147  }
148 
149 
150  } // Loop over the tTrig map reference
151 
152 }
153 
154 void DTtTrigDBValidation::bookHistos(DQMStore::IBooker &iBooker, int wheel, int sector) {
155 
156  LogTrace(metname_) << " Booking histos for Wheel, Sector: " << wheel << ", " << sector;
157 
158  // Compose the chamber name
159  stringstream str_wheel; str_wheel << wheel;
160  stringstream str_sector; str_sector << sector;
161 
162  string lHistoName = "_W" + str_wheel.str() + "_Sec" + str_sector.str();
163 
164  iBooker.setCurrentFolder("DT/DtCalib/TTrigDBValidation/Wheel" + str_wheel.str());
165 
166  // Create the monitor elements
167  MonitorElement * hDifference;
168  hDifference = iBooker.book1D("TTrigDifference"+lHistoName, "difference between the two tTrig values",11,0,11);
169 
170  pair<int,int> mypair(wheel,sector);
171  tTrigDiffHistos_[mypair] = hDifference;
172 
173  (tTrigDiffHistos_[mypair])->setBinLabel(1,"MB1_SL1",1);
174  (tTrigDiffHistos_[mypair])->setBinLabel(2,"MB1_SL2",1);
175  (tTrigDiffHistos_[mypair])->setBinLabel(3,"MB1_SL3",1);
176  (tTrigDiffHistos_[mypair])->setBinLabel(4,"MB2_SL1",1);
177  (tTrigDiffHistos_[mypair])->setBinLabel(5,"MB2_SL2",1);
178  (tTrigDiffHistos_[mypair])->setBinLabel(6,"MB2_SL3",1);
179  (tTrigDiffHistos_[mypair])->setBinLabel(7,"MB3_SL1",1);
180  (tTrigDiffHistos_[mypair])->setBinLabel(8,"MB3_SL2",1);
181  (tTrigDiffHistos_[mypair])->setBinLabel(9,"MB3_SL3",1);
182  (tTrigDiffHistos_[mypair])->setBinLabel(10,"MB4_SL1",1);
183  (tTrigDiffHistos_[mypair])->setBinLabel(11,"MB4_SL3",1);
184 
185 }
186 
187 // Book the summary histos
189 
190  stringstream wh; wh << wheel;
191 
192  iBooker.setCurrentFolder("DT/DtCalib/TTrigDBValidation");
193  tTrigDiffWheel_[wheel] = iBooker.book2D("TTrigDifference_W"+wh.str(), "W"+wh.str()+": summary of tTrig differences",11,1,12,14,1,15);
194  tTrigDiffWheel_[wheel]->setBinLabel(1,"MB1_SL1",1);
195  tTrigDiffWheel_[wheel]->setBinLabel(2,"MB1_SL2",1);
196  tTrigDiffWheel_[wheel]->setBinLabel(3,"MB1_SL3",1);
197  tTrigDiffWheel_[wheel]->setBinLabel(4,"MB2_SL1",1);
198  tTrigDiffWheel_[wheel]->setBinLabel(5,"MB2_SL2",1);
199  tTrigDiffWheel_[wheel]->setBinLabel(6,"MB2_SL3",1);
200  tTrigDiffWheel_[wheel]->setBinLabel(7,"MB3_SL1",1);
201  tTrigDiffWheel_[wheel]->setBinLabel(8,"MB3_SL2",1);
202  tTrigDiffWheel_[wheel]->setBinLabel(9,"MB3_SL3",1);
203  tTrigDiffWheel_[wheel]->setBinLabel(10,"MB4_SL1",1);
204  tTrigDiffWheel_[wheel]->setBinLabel(11,"MB4_SL3",1);
205 
206 }
207 
208 
210  return (int) (bin /3.1)+1;
211 }
212 
214  int ret = bin%3;
215  if(ret == 0 || bin == 11) ret = 3;
216 
217  return ret;
218 }
219 
220 
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Operations.
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
std::map< int, MonitorElement * > tTrigDiffWheel_
int slFromBin(int bin) const
std::map< std::pair< int, int >, MonitorElement * > tTrigDiffHistos_
std::vector< std::pair< DTTtrigId, DTTtrigData > >::const_iterator const_iterator
Access methods to data.
Definition: DTTtrig.h:183
edm::ESHandle< DTGeometry > dtGeom_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int stationFromBin(int bin) const
#define LogTrace(id)
virtual ~DTtTrigDBValidation()
Destructor.
std::map< DTSuperLayerId, std::pair< float, float > > tTrigRefMap_
DTtTrigDBValidation(const edm::ParameterSet &pset)
Constructor.
std::map< DTSuperLayerId, std::pair< float, float > > tTrigMap_
int get(int wheelId, int stationId, int sectorId, int slId, float &tTrig, float &tTrms, float &kFact, DTTimeUnits::type unit) const
get content
Definition: DTTtrig.cc:85
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
void bookHistos(DQMStore::IBooker &, int, int)
const T & get() const
Definition: EventSetup.h:55
const_iterator begin() const
Definition: DTTtrig.cc:354
const_iterator end() const
Definition: DTTtrig.cc:359
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
Definition: Run.h:41