CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTResidualCalibration.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * $Date: $
6  * $Revision: $
7  */
8 
10 
11 // Framework
17 
18 //Geometry
21 
22 //RecHit
25 
29 
30 #include "TFile.h"
31 #include "TH1F.h"
32 #include "TH2F.h"
33 
34 #include <algorithm>
35 
37  select_(pset),
38  segment4DLabel_(pset.getParameter<edm::InputTag>("segment4DLabel")),
39  rootBaseDir_(pset.getUntrackedParameter<std::string>("rootBaseDir","DT/Residuals")) {
40 
41  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Constructor called.";
42 
43  std::string rootFileName = pset.getUntrackedParameter<std::string>("rootFileName","residuals.root");
44  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
45  rootFile_->cd();
46 }
47 
49  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Destructor called.";
50 }
51 
53  TH1::SetDefaultSumw2(true);
54 }
55 
57 
58  // get the geometry
60  setup.get<MuonGeometryRecord>().get(dtGeomH);
61  dtGeom_ = dtGeomH.product();
62 
63  // Loop over all the chambers
64  if(histoMapTH1F_.size() == 0){
65  std::vector<DTChamber*>::const_iterator ch_it = dtGeom_->chambers().begin();
66  std::vector<DTChamber*>::const_iterator ch_end = dtGeom_->chambers().end();
67  for (; ch_it != ch_end; ++ch_it) {
68  DTChamberId chId = (*ch_it)->id();
69  std::vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
70  std::vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
71  // Loop over the SLs
72  for(; sl_it != sl_end; ++sl_it) {
73  DTSuperLayerId slId = (*sl_it)->id();
74  bookHistos(slId);
75  }
76  }
77  }
78 }
79 
81  rootFile_->cd();
82 
83  // Get the 4D rechits from the event
85  event.getByLabel(segment4DLabel_, segment4Ds);
86 
87  // Loop over segments by chamber
89  for(chamberIdIt = segment4Ds->id_begin(); chamberIdIt != segment4Ds->id_end(); ++chamberIdIt){
90 
91  const DTChamber* chamber = dtGeom_->chamber(*chamberIdIt);
92 
93  // Get the range for the corresponding ChamberId
94  DTRecSegment4DCollection::range range = segment4Ds->get((*chamberIdIt));
95  // Loop over the rechits of this DetUnit
96  for(DTRecSegment4DCollection::const_iterator segment = range.first;
97  segment != range.second; ++segment){
98 
99  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
100  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
101 
102  if( !select_(*segment, event, setup) ) continue;
103 
104  // Get all 1D RecHits at step 3 within the 4D segment
105  std::vector<DTRecHit1D> recHits1D_S3;
106 
107  if( (*segment).hasPhi() ){
108  const DTChamberRecSegment2D* phiSeg = (*segment).phiSegment();
109  const std::vector<DTRecHit1D>& phiRecHits = phiSeg->specificRecHits();
110  std::copy(phiRecHits.begin(), phiRecHits.end(), back_inserter(recHits1D_S3));
111  }
112 
113  if( (*segment).hasZed() ){
114  const DTSLRecSegment2D* zSeg = (*segment).zSegment();
115  const std::vector<DTRecHit1D>& zRecHits = zSeg->specificRecHits();
116  std::copy(zRecHits.begin(), zRecHits.end(), back_inserter(recHits1D_S3));
117  }
118 
119  // Loop over 1D RecHit inside 4D segment
120  for(std::vector<DTRecHit1D>::const_iterator recHit1D = recHits1D_S3.begin();
121  recHit1D != recHits1D_S3.end(); ++recHit1D) {
122  const DTWireId wireId = recHit1D->wireId();
123 
124  float segmDistance = segmentToWireDistance(*recHit1D,*segment);
125  if(segmDistance > 2.1) LogTrace("Calibration") << "WARNING: segment-wire distance: " << segmDistance;
126  else LogTrace("Calibration") << "segment-wire distance: " << segmDistance;
127 
128  float residualOnDistance = DTRecHitSegmentResidual().compute(dtGeom_,*recHit1D,*segment);
129  LogTrace("Calibration") << "Wire Id " << wireId << " residual on distance: " << residualOnDistance;
130 
131  fillHistos(wireId.superlayerId(), segmDistance, residualOnDistance);
132  }
133  }
134  }
135 
136 }
137 
139 
140  // Get the layer and the wire position
141  const DTWireId wireId = recHit1D.wireId();
142  const DTLayer* layer = dtGeom_->layer(wireId);
143  float wireX = layer->specificTopology().wirePosition(wireId.wire());
144 
145  // Extrapolate the segment to the z of the wire
146  // Get wire position in chamber RF
147  // (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)
148  LocalPoint wirePosInLay(wireX,recHit1D.localPosition().y(),recHit1D.localPosition().z());
149  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
150  const DTChamber* chamber = dtGeom_->chamber(wireId.layerId().chamberId());
151  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
152 
153  // Segment position at Wire z in chamber local frame
154  LocalPoint segPosAtZWire = segment.localPosition() + segment.localDirection()*wirePosInChamber.z()/cos(segment.localDirection().theta());
155 
156  // Compute the distance of the segment from the wire
157  int sl = wireId.superlayer();
158  float segmDistance = -1;
159  if(sl == 1 || sl == 3) segmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
160  else if(sl == 2) segmDistance = fabs(segPosAtZWire.y() - wirePosInChamber.y());
161 
162  return segmDistance;
163 }
164 
166 
167  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Writing histos to file.";
168  rootFile_->cd();
169  rootFile_->Write();
170  rootFile_->Close();
171 
172  /*std::map<DTSuperLayerId, std::vector<TH1F*> >::const_iterator itSlHistos = histoMapTH1F_.begin();
173  std::map<DTSuperLayerId, std::vector<TH1F*> >::const_iterator itSlHistos_end = histoMapTH1F_.end();
174  for(; itSlHistos != itSlHistos_end; ++itSlHistos){
175  std::vector<TH1F*>::const_iterator itHistTH1F = (*itSlHistos).second.begin();
176  std::vector<TH1F*>::const_iterator itHistTH1F_end = (*itSlHistos).second.end();
177  for(; itHistTH1F != itHistTH1F_end; ++itHistTH1F) (*itHistTH1F)->Write();
178 
179  std::vector<TH2F*>::const_iterator itHistTH2F = histoMapTH2F_[(*itSlHistos).first].begin();
180  std::vector<TH2F*>::const_iterator itHistTH2F_end = histoMapTH2F_[(*itSlHistos).first].end();
181  for(; itHistTH2F != itHistTH2F_end; ++itHistTH2F) (*itHistTH2F)->Write();
182  }*/
183 
184 }
185 
187  TH1AddDirectorySentry addDir;
188  rootFile_->cd();
189 
190  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Booking histos for SL: " << slId;
191 
192  // Compose the chamber name
193  std::stringstream wheelStr; wheelStr << slId.wheel();
194  std::stringstream stationStr; stationStr << slId.station();
195  std::stringstream sectorStr; sectorStr << slId.sector();
196  std::stringstream superLayerStr; superLayerStr << slId.superlayer();
197  // Define the step
198  int step = 3;
199  std::stringstream stepStr; stepStr << step;
200 
201  std::string slHistoName =
202  "_STEP" + stepStr.str() +
203  "_W" + wheelStr.str() +
204  "_St" + stationStr.str() +
205  "_Sec" + sectorStr.str() +
206  "_SL" + superLayerStr.str();
207 
208  edm::LogVerbatim("Calibration") << "Accessing " << rootBaseDir_;
209  TDirectory* baseDir = rootFile_->GetDirectory(rootBaseDir_.c_str());
210  if(!baseDir) baseDir = rootFile_->mkdir(rootBaseDir_.c_str());
211  edm::LogVerbatim("Calibration") << "Accessing " << ("Wheel" + wheelStr.str());
212  TDirectory* wheelDir = baseDir->GetDirectory(("Wheel" + wheelStr.str()).c_str());
213  if(!wheelDir) wheelDir = baseDir->mkdir(("Wheel" + wheelStr.str()).c_str());
214  edm::LogVerbatim("Calibration") << "Accessing " << ("Station" + stationStr.str());
215  TDirectory* stationDir = wheelDir->GetDirectory(("Station" + stationStr.str()).c_str());
216  if(!stationDir) stationDir = wheelDir->mkdir(("Station" + stationStr.str()).c_str());
217  edm::LogVerbatim("Calibration") << "Accessing " << ("Sector" + sectorStr.str());
218  TDirectory* sectorDir = stationDir->GetDirectory(("Sector" + sectorStr.str()).c_str());
219  if(!sectorDir) sectorDir = stationDir->mkdir(("Sector" + sectorStr.str()).c_str());
220 
221  /*std::string dirName = rootBaseDir_ + "/Wheel" + wheelStr.str() +
222  "/Station" + stationStr.str() +
223  "/Sector" + sectorStr.str();
224 
225  TDirectory* dir = rootFile_->GetDirectory(dirName.c_str());
226  if(!dir) dir = rootFile_->mkdir(dirName.c_str());
227  dir->cd();*/
228  sectorDir->cd();
229  // Create the monitor elements
230  std::vector<TH1F*> histosTH1F;
231  histosTH1F.push_back(new TH1F(("hResDist"+slHistoName).c_str(),
232  "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
233  200, -0.4, 0.4));
234  std::vector<TH2F*> histosTH2F;
235  histosTH2F.push_back(new TH2F(("hResDistVsDist"+slHistoName).c_str(),
236  "Residuals on the dist. (cm) from wire (rec_hit - segm_extr) vs dist. (cm)",
237  100, 0, 2.5, 200, -0.4, 0.4));
238  histoMapTH1F_[slId] = histosTH1F;
239  histoMapTH2F_[slId] = histosTH2F;
240 }
241 
242 // Fill a set of histograms for a given SL
244  float distance,
245  float residualOnDistance) {
246  std::vector<TH1F*> const& histosTH1F = histoMapTH1F_[slId];
247  std::vector<TH2F*> const& histosTH2F = histoMapTH2F_[slId];
248  histosTH1F[0]->Fill(residualOnDistance);
249  histosTH2F[0]->Fill(distance, residualOnDistance);
250 }
T getUntrackedParameter(std::string const &, T const &) const
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:73
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:53
void fillHistos(DTSuperLayerId slId, float distance, float residualOnDistance)
void bookHistos(DTSuperLayerId slId)
list step
Definition: launcher.py:15
DTResidualCalibration(const edm::ParameterSet &pset)
Constructor.
DTChamberId chamberId() const
Return the corresponding ChamberId.
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
virtual ~DTResidualCalibration()
Destructor.
T y() const
Definition: PV3DBase.h:57
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:64
identifier iterator
Definition: RangeMap.h:139
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:61
Geom::Theta< T > theta() const
Definition: PV3DBase.h:69
void analyze(const edm::Event &event, const edm::EventSetup &setup)
const DTLayer * layer(DTLayerId id) const
Return a layer given its id.
Definition: DTGeometry.cc:115
virtual LocalVector localDirection() const
Local direction in Chamber frame.
void beginRun(const edm::Run &, const edm::EventSetup &)
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:46
const DTTopology & specificTopology() const
Definition: DTLayer.cc:44
std::map< DTSuperLayerId, std::vector< TH1F * > > histoMapTH1F_
T z() const
Definition: PV3DBase.h:58
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::map< DTSuperLayerId, std::vector< TH2F * > > histoMapTH2F_
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
virtual LocalPoint localPosition() const
Local position in Chamber frame.
#define LogTrace(id)
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
const std::vector< DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:90
virtual LocalPoint localPosition() const
Return the 3-dimensional local position.
Definition: DTRecHit1D.h:62
int wire() const
Return the wire number.
Definition: DTWireId.h:58
int superlayer() const
Return the superlayer number (deprecated method name)
float compute(const DTGeometry *, const DTRecHit1D &, const DTRecSegment4D &)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
const DTChamber * chamber(DTChamberId id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:105
const DTGeometry * dtGeom_
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:64
int sector() const
Definition: DTChamberId.h:63
float segmentToWireDistance(const DTRecHit1D &recHit1D, const DTRecSegment4D &segment)
int station() const
Return the station number.
Definition: DTChamberId.h:53
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
T x() const
Definition: PV3DBase.h:56
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
Definition: Run.h:32
DTWireId wireId() const
Return the wireId.
Definition: DTRecHit1D.h:109