CMS 3D CMS Logo

DTResidualCalibration.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  */
6 
8 
9 // Framework
15 
16 //Geometry
19 
20 //RecHit
23 
27 
28 #include "TFile.h"
29 #include "TH1F.h"
30 #include "TH2F.h"
31 
32 #include <algorithm>
33 
35  : histRange_(pset.getParameter<double>("histogramRange")),
36  segment4DToken_(consumes<DTRecSegment4DCollection>(pset.getParameter<edm::InputTag>("segment4DLabel"))),
37  rootBaseDir_(pset.getUntrackedParameter<std::string>("rootBaseDir", "DT/Residuals")),
38  detailedAnalysis_(pset.getUntrackedParameter<bool>("detailedAnalysis", false)),
39  dtGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
41  select_ = new DTSegmentSelector(pset, collector);
42 
43  LogDebug("Calibration") << "[DTResidualCalibration] Constructor called.";
44  std::string rootFileName = pset.getUntrackedParameter<std::string>("rootFileName", "residuals.root");
45  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
46  rootFile_->cd();
47 
48  segmok = 0;
49  segmbad = 0;
50  nevent = 0;
51 }
52 
54  delete select_;
55  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Destructor called.";
56  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Analyzed events: " << nevent;
57  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Good segments: " << segmok;
58  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Bad segments: " << segmbad;
59 }
60 
61 void DTResidualCalibration::beginJob() { TH1::SetDefaultSumw2(true); }
62 
64  // get the geometry
66  dtGeomH = setup.getHandle(dtGeomToken_);
67  dtGeom_ = dtGeomH.product();
68 
69  // Loop over all the chambers
70  if (histoMapTH1F_.empty()) {
71  for (auto ch_it : dtGeom_->chambers()) {
72  // Loop over the SLs
73  for (auto sl_it : ch_it->superLayers()) {
74  DTSuperLayerId slId = (sl_it)->id();
75  bookHistos(slId);
76  if (detailedAnalysis_) {
77  for (auto layer_it : (sl_it)->layers()) {
78  DTLayerId layerId = (layer_it)->id();
79  bookHistos(layerId);
80  }
81  }
82  }
83  }
84  }
85 }
86 
88  rootFile_->cd();
89  ++nevent;
90 
91  // Get the 4D rechits from the event
92  const edm::Handle<DTRecSegment4DCollection>& segments4D = event.getHandle(segment4DToken_);
93 
94  // Loop over segments by chamber
96  for (chamberIdIt = segments4D->id_begin(); chamberIdIt != segments4D->id_end(); ++chamberIdIt) {
97  const DTChamber* chamber = dtGeom_->chamber(*chamberIdIt);
98 
99  // Get the range for the corresponding ChamberId
100  DTRecSegment4DCollection::range range = segments4D->get((*chamberIdIt));
101  // Loop over the rechits of this DetUnit
102  for (DTRecSegment4DCollection::const_iterator segment = range.first; segment != range.second; ++segment) {
103  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
104  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
105 
106  if (!(*select_)(*segment, event, setup)) {
107  segmbad++;
108  continue;
109  }
110  segmok++;
111 
112  // Get all 1D RecHits at step 3 within the 4D segment
113  std::vector<DTRecHit1D> recHits1D_S3;
114 
115  if ((*segment).hasPhi()) {
116  const DTChamberRecSegment2D* phiSeg = (*segment).phiSegment();
117  const std::vector<DTRecHit1D>& phiRecHits = phiSeg->specificRecHits();
118  std::copy(phiRecHits.begin(), phiRecHits.end(), back_inserter(recHits1D_S3));
119  }
120 
121  if ((*segment).hasZed()) {
122  const DTSLRecSegment2D* zSeg = (*segment).zSegment();
123  const std::vector<DTRecHit1D>& zRecHits = zSeg->specificRecHits();
124  std::copy(zRecHits.begin(), zRecHits.end(), back_inserter(recHits1D_S3));
125  }
126 
127  // Loop over 1D RecHit inside 4D segment
128  for (std::vector<DTRecHit1D>::const_iterator recHit1D = recHits1D_S3.begin(); recHit1D != recHits1D_S3.end();
129  ++recHit1D) {
130  const DTWireId wireId = recHit1D->wireId();
131 
132  float segmDistance = segmentToWireDistance(*recHit1D, *segment);
133  if (segmDistance > 2.1)
134  LogTrace("Calibration") << "WARNING: segment-wire distance: " << segmDistance;
135  else
136  LogTrace("Calibration") << "segment-wire distance: " << segmDistance;
137 
138  float residualOnDistance = DTRecHitSegmentResidual().compute(dtGeom_, *recHit1D, *segment);
139  LogTrace("Calibration") << "Wire Id " << wireId << " residual on distance: " << residualOnDistance;
140 
141  fillHistos(wireId.superlayerId(), segmDistance, residualOnDistance);
142  if (detailedAnalysis_)
143  fillHistos(wireId.layerId(), segmDistance, residualOnDistance);
144  }
145  }
146  }
147 }
148 
150  // Get the layer and the wire position
151  const DTWireId wireId = recHit1D.wireId();
152  const DTLayer* layer = dtGeom_->layer(wireId);
153  float wireX = layer->specificTopology().wirePosition(wireId.wire());
154 
155  // Extrapolate the segment to the z of the wire
156  // Get wire position in chamber RF
157  // (y and z must be those of the hit to be coherent in the transf. of RF in case of rotations of the layer alignment)
158  LocalPoint wirePosInLay(wireX, recHit1D.localPosition().y(), recHit1D.localPosition().z());
159  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
160  const DTChamber* chamber = dtGeom_->chamber(wireId.layerId().chamberId());
161  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
162 
163  // Segment position at Wire z in chamber local frame
164  LocalPoint segPosAtZWire =
165  segment.localPosition() + segment.localDirection() * wirePosInChamber.z() / cos(segment.localDirection().theta());
166 
167  // Compute the distance of the segment from the wire
168  int sl = wireId.superlayer();
169  float segmDistance = -1;
170  if (sl == 1 || sl == 3)
171  segmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
172  else if (sl == 2)
173  segmDistance = fabs(segPosAtZWire.y() - wirePosInChamber.y());
174 
175  return segmDistance;
176 }
177 
179  LogDebug("Calibration") << "[DTResidualCalibration] Writing histos to file.";
180  rootFile_->cd();
181  rootFile_->Write();
182  rootFile_->Close();
183 
184  /*std::map<DTSuperLayerId, TH1F* >::const_iterator itSlHistos = histoMapTH1F_.begin();
185  std::map<DTSuperLayerId, TH1F* >::const_iterator itSlHistos_end = histoMapTH1F_.end();
186  for(; itSlHistos != itSlHistos_end; ++itSlHistos){
187  std::vector<TH1F*>::const_iterator itHistTH1F = (*itSlHistos).second.begin();
188  std::vector<TH1F*>::const_iterator itHistTH1F_end = (*itSlHistos).second.end();
189  for(; itHistTH1F != itHistTH1F_end; ++itHistTH1F) (*itHistTH1F)->Write();
190 
191  std::vector<TH2F*>::const_iterator itHistTH2F = histoMapTH2F_[(*itSlHistos).first].begin();
192  std::vector<TH2F*>::const_iterator itHistTH2F_end = histoMapTH2F_[(*itSlHistos).first].end();
193  for(; itHistTH2F != itHistTH2F_end; ++itHistTH2F) (*itHistTH2F)->Write();
194  }*/
195 }
196 
198  TH1AddDirectorySentry addDir;
199  rootFile_->cd();
200 
201  LogDebug("Calibration") << "[DTResidualCalibration] Booking histos for SL: " << slId;
202 
203  // Compose the chamber name
204  // Define the step
205  int step = 3;
206 
207  std::string wheelStr = std::to_string(slId.wheel());
208  std::string stationStr = std::to_string(slId.station());
209  std::string sectorStr = std::to_string(slId.sector());
210 
211  std::string slHistoName = "_STEP" + std::to_string(step) + "_W" + wheelStr + "_St" + stationStr + "_Sec" + sectorStr +
212  "_SL" + std::to_string(slId.superlayer());
213 
214  LogDebug("Calibration") << "Accessing " << rootBaseDir_;
215  TDirectory* baseDir = rootFile_->GetDirectory(rootBaseDir_.c_str());
216  if (!baseDir)
217  baseDir = rootFile_->mkdir(rootBaseDir_.c_str());
218  LogDebug("Calibration") << "Accessing " << ("Wheel" + wheelStr);
219  TDirectory* wheelDir = baseDir->GetDirectory(("Wheel" + wheelStr).c_str());
220  if (!wheelDir)
221  wheelDir = baseDir->mkdir(("Wheel" + wheelStr).c_str());
222  LogDebug("Calibration") << "Accessing " << ("Station" + stationStr);
223  TDirectory* stationDir = wheelDir->GetDirectory(("Station" + stationStr).c_str());
224  if (!stationDir)
225  stationDir = wheelDir->mkdir(("Station" + stationStr).c_str());
226  LogDebug("Calibration") << "Accessing " << ("Sector" + sectorStr);
227  TDirectory* sectorDir = stationDir->GetDirectory(("Sector" + sectorStr).c_str());
228  if (!sectorDir)
229  sectorDir = stationDir->mkdir(("Sector" + sectorStr).c_str());
230 
231  sectorDir->cd();
232 
233  // Create the monitor elements
234  TH1F* histosTH1F = new TH1F(("hResDist" + slHistoName).c_str(),
235  "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
236  200,
237  -histRange_,
238  histRange_);
239  TH2F* histosTH2F = new TH2F(("hResDistVsDist" + slHistoName).c_str(),
240  "Residuals on the dist. (cm) from wire (rec_hit - segm_extr) vs dist. (cm)",
241  100,
242  0,
243  2.5,
244  200,
245  -histRange_,
246  histRange_);
247  histoMapTH1F_[slId] = histosTH1F;
248  histoMapTH2F_[slId] = histosTH2F;
249 }
250 
252  TH1AddDirectorySentry addDir;
253  rootFile_->cd();
254 
255  LogDebug("Calibration") << "[DTResidualCalibration] Booking histos for layer: " << layerId;
256 
257  // Compose the chamber name
258  std::string wheelStr = std::to_string(layerId.wheel());
259  std::string stationStr = std::to_string(layerId.station());
260  std::string sectorStr = std::to_string(layerId.sector());
261  std::string superLayerStr = std::to_string(layerId.superlayer());
262  std::string layerStr = std::to_string(layerId.layer());
263  // Define the step
264  int step = 3;
265 
266  std::string layerHistoName = "_STEP" + std::to_string(step) + "_W" + wheelStr + "_St" + stationStr + "_Sec" +
267  sectorStr + "_SL" + superLayerStr + "_Layer" + layerStr;
268 
269  LogDebug("Calibration") << "Accessing " << rootBaseDir_;
270  TDirectory* baseDir = rootFile_->GetDirectory(rootBaseDir_.c_str());
271  if (!baseDir)
272  baseDir = rootFile_->mkdir(rootBaseDir_.c_str());
273  LogDebug("Calibration") << "Accessing " << ("Wheel" + wheelStr);
274  TDirectory* wheelDir = baseDir->GetDirectory(("Wheel" + wheelStr).c_str());
275  if (!wheelDir)
276  wheelDir = baseDir->mkdir(("Wheel" + wheelStr).c_str());
277  LogDebug("Calibration") << "Accessing " << ("Station" + stationStr);
278  TDirectory* stationDir = wheelDir->GetDirectory(("Station" + stationStr).c_str());
279  if (!stationDir)
280  stationDir = wheelDir->mkdir(("Station" + stationStr).c_str());
281  LogDebug("Calibration") << "Accessing " << ("Sector" + sectorStr);
282  TDirectory* sectorDir = stationDir->GetDirectory(("Sector" + sectorStr).c_str());
283  if (!sectorDir)
284  sectorDir = stationDir->mkdir(("Sector" + sectorStr).c_str());
285  LogDebug("Calibration") << "Accessing " << ("SL" + superLayerStr);
286  TDirectory* superLayerDir = sectorDir->GetDirectory(("SL" + superLayerStr).c_str());
287  if (!superLayerDir)
288  superLayerDir = sectorDir->mkdir(("SL" + superLayerStr).c_str());
289 
290  superLayerDir->cd();
291  // Create histograms
292  TH1F* histosTH1F = new TH1F(("hResDist" + layerHistoName).c_str(),
293  "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
294  200,
295  -histRange_,
296  histRange_);
297  TH2F* histosTH2F = new TH2F(("hResDistVsDist" + layerHistoName).c_str(),
298  "Residuals on the dist. (cm) from wire (rec_hit - segm_extr) vs dist. (cm)",
299  100,
300  0,
301  2.5,
302  200,
303  -histRange_,
304  histRange_);
305  histoMapPerLayerTH1F_[layerId] = histosTH1F;
306  histoMapPerLayerTH2F_[layerId] = histosTH2F;
307 }
308 
309 // Fill a set of histograms for a given SL
310 void DTResidualCalibration::fillHistos(DTSuperLayerId slId, float distance, float residualOnDistance) {
311  histoMapTH1F_[slId]->Fill(residualOnDistance);
312  histoMapTH2F_[slId]->Fill(distance, residualOnDistance);
313 }
314 
315 // Fill a set of histograms for a given layer
316 void DTResidualCalibration::fillHistos(DTLayerId layerId, float distance, float residualOnDistance) {
317  histoMapPerLayerTH1F_[layerId]->Fill(residualOnDistance);
318  histoMapPerLayerTH2F_[layerId]->Fill(distance, residualOnDistance);
319 }
Log< level::Info, true > LogVerbatim
int station() const
Return the station number.
Definition: DTChamberId.h:42
std::map< DTSuperLayerId, TH1F * > histoMapTH1F_
std::map< DTLayerId, TH2F * > histoMapPerLayerTH2F_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
DTSegmentSelector * select_
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
int wire() const
Return the wire number.
Definition: DTWireId.h:42
void fillHistos(DTSuperLayerId slId, float distance, float residualOnDistance)
void bookHistos(DTSuperLayerId slId)
DTResidualCalibration(const edm::ParameterSet &pset)
Constructor.
T z() const
Definition: PV3DBase.h:61
LocalVector localDirection() const override
Local direction in Chamber frame.
identifier iterator
Definition: RangeMap.h:130
LocalPoint localPosition() const override
Local position in Chamber frame.
static std::string to_string(const XMLCh *ch)
#define LogTrace(id)
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T const * product() const
Definition: ESHandle.h:86
DTChamberId chamberId() const
Return the corresponding ChamberId.
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
void beginRun(const edm::Run &, const edm::EventSetup &) override
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
const std::string rootBaseDir_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Transition
Definition: Transition.h:12
~DTResidualCalibration() override
Destructor.
std::map< DTSuperLayerId, TH2F * > histoMapTH2F_
int superlayer() const
Return the superlayer number (deprecated method name)
void analyze(const edm::Event &event, const edm::EventSetup &setup) override
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
float compute(const DTGeometry *, const DTRecHit1D &, const DTRecSegment4D &)
int layer() const
Return the layer number.
Definition: DTLayerId.h:42
const DTGeometry * dtGeom_
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
HLT enums.
DTWireId wireId() const
Return the wireId.
Definition: DTRecHit1D.h:76
int sector() const
Definition: DTChamberId.h:49
const edm::EDGetTokenT< DTRecSegment4DCollection > segment4DToken_
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:45
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
step
Definition: StallMonitor.cc:83
float segmentToWireDistance(const DTRecHit1D &recHit1D, const DTRecSegment4D &segment)
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:45
std::map< DTLayerId, TH1F * > histoMapPerLayerTH1F_
LocalPoint localPosition() const override
Return the 3-dimensional local position.
Definition: DTRecHit1D.h:47
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90
Definition: event.py:1
Definition: Run.h:45
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
#define LogDebug(id)
const DTLayer * layer(const DTLayerId &id) const
Return a layer given its id.
Definition: DTGeometry.cc:96