CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTResidualCalibration Class Reference

#include <DTResidualCalibration.h>

Inheritance diagram for DTResidualCalibration:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &setup) override
 
void beginJob () override
 
void beginRun (const edm::Run &, const edm::EventSetup &) override
 
 DTResidualCalibration (const edm::ParameterSet &pset)
 Constructor. More...
 
void endJob () override
 
 ~DTResidualCalibration () override
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void bookHistos (DTSuperLayerId slId)
 
void bookHistos (DTLayerId slId)
 
void fillHistos (DTSuperLayerId slId, float distance, float residualOnDistance)
 
void fillHistos (DTLayerId slId, float distance, float residualOnDistance)
 
float segmentToWireDistance (const DTRecHit1D &recHit1D, const DTRecSegment4D &segment)
 

Private Attributes

bool detailedAnalysis_
 
const DTGeometrydtGeom_
 
std::map< DTLayerId, std::vector< TH1F * > > histoMapPerLayerTH1F_
 
std::map< DTLayerId, std::vector< TH2F * > > histoMapPerLayerTH2F_
 
std::map< DTSuperLayerId, std::vector< TH1F * > > histoMapTH1F_
 
std::map< DTSuperLayerId, std::vector< TH2F * > > histoMapTH2F_
 
double histRange_
 
unsigned int nevent
 
std::string rootBaseDir_
 
TFile * rootFile_
 
unsigned int segmbad
 
edm::InputTag segment4DLabel_
 
unsigned int segmok
 
DTSegmentSelector select_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Extracts DT segment residuals

Definition at line 28 of file DTResidualCalibration.h.

Constructor & Destructor Documentation

DTResidualCalibration::DTResidualCalibration ( const edm::ParameterSet pset)

Constructor.

Definition at line 34 of file DTResidualCalibration.cc.

References edm::ParameterSet::getUntrackedParameter(), nevent, rootFile_, DTAnalyzerDetailed_cfi::rootFileName, segmbad, segment4DLabel_, segmok, and AlCaHLTBitMon_QueryRunRegistry::string.

34  :
35  select_(pset),
36  histRange_(pset.getParameter<double>("histogramRange")),
37  segment4DLabel_(pset.getParameter<edm::InputTag>("segment4DLabel")),
38  rootBaseDir_(pset.getUntrackedParameter<std::string>("rootBaseDir","DT/Residuals")),
39  detailedAnalysis_(pset.getUntrackedParameter<bool>("detailedAnalysis",false)) {
40 
41  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Constructor called.";
42  consumes< DTRecSegment4DCollection >(edm::InputTag(segment4DLabel_));
43  std::string rootFileName = pset.getUntrackedParameter<std::string>("rootFileName","residuals.root");
44  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
45  rootFile_->cd();
46 
47  segmok=0;
48  segmbad=0;
49  nevent=0;
50 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTResidualCalibration::~DTResidualCalibration ( )
override

Destructor.

Definition at line 52 of file DTResidualCalibration.cc.

References nevent, segmbad, and segmok.

52  {
53  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Destructor called.";
54  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Analyzed events: " << nevent;
55  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Good segments: " << segmok;
56  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Bad segments: " << segmbad;
57 }

Member Function Documentation

void DTResidualCalibration::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
override

Definition at line 88 of file DTResidualCalibration.cc.

References relativeConstraints::chamber, DTGeometry::chamber(), DTRecHitSegmentResidual::compute(), popcon2dropbox::copy(), detailedAnalysis_, dtGeom_, fillHistos(), DTWireId::layerId(), LogTrace, nevent, rootFile_, segmbad, segment4DLabel_, segmentToWireDistance(), segmok, select_, DTRecSegment2D::specificRecHits(), DTLayerId::superlayerId(), and GeomDet::toGlobal().

88  {
89  rootFile_->cd();
90  ++nevent;
91 
92  // Get the 4D rechits from the event
94  event.getByLabel(segment4DLabel_, segments4D);
95 
96  // Loop over segments by chamber
97  DTRecSegment4DCollection::id_iterator chamberIdIt;
98  for(chamberIdIt = segments4D->id_begin(); chamberIdIt != segments4D->id_end(); ++chamberIdIt){
99 
100  const DTChamber* chamber = dtGeom_->chamber(*chamberIdIt);
101 
102  // Get the range for the corresponding ChamberId
103  DTRecSegment4DCollection::range range = segments4D->get((*chamberIdIt));
104  // Loop over the rechits of this DetUnit
105  for(DTRecSegment4DCollection::const_iterator segment = range.first;
106  segment != range.second; ++segment){
107 
108  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
109  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
110 
111  if( !select_(*segment, event, setup) ) { segmbad++; continue; }
112  segmok++;
113 
114  // Get all 1D RecHits at step 3 within the 4D segment
115  std::vector<DTRecHit1D> recHits1D_S3;
116 
117  if( (*segment).hasPhi() ){
118  const DTChamberRecSegment2D* phiSeg = (*segment).phiSegment();
119  const std::vector<DTRecHit1D>& phiRecHits = phiSeg->specificRecHits();
120  std::copy(phiRecHits.begin(), phiRecHits.end(), back_inserter(recHits1D_S3));
121  }
122 
123  if( (*segment).hasZed() ){
124  const DTSLRecSegment2D* zSeg = (*segment).zSegment();
125  const std::vector<DTRecHit1D>& zRecHits = zSeg->specificRecHits();
126  std::copy(zRecHits.begin(), zRecHits.end(), back_inserter(recHits1D_S3));
127  }
128 
129  // Loop over 1D RecHit inside 4D segment
130  for(std::vector<DTRecHit1D>::const_iterator recHit1D = recHits1D_S3.begin();
131  recHit1D != recHits1D_S3.end(); ++recHit1D) {
132  const DTWireId wireId = recHit1D->wireId();
133 
134  float segmDistance = segmentToWireDistance(*recHit1D,*segment);
135  if(segmDistance > 2.1) LogTrace("Calibration") << "WARNING: segment-wire distance: " << segmDistance;
136  else 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_) fillHistos(wireId.layerId(), segmDistance, residualOnDistance);
143  }
144  }
145  }
146 }
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
def copy(args, dbName)
void fillHistos(DTSuperLayerId slId, float distance, float residualOnDistance)
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:99
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:59
#define LogTrace(id)
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
float compute(const DTGeometry *, const DTRecHit1D &, const DTRecSegment4D &)
const DTGeometry * dtGeom_
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:62
float segmentToWireDistance(const DTRecHit1D &recHit1D, const DTRecSegment4D &segment)
void DTResidualCalibration::beginJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 59 of file DTResidualCalibration.cc.

59  {
60  TH1::SetDefaultSumw2(true);
61 }
void DTResidualCalibration::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)
override

Definition at line 63 of file DTResidualCalibration.cc.

References bookHistos(), DTGeometry::chambers(), detailedAnalysis_, dtGeom_, edm::EventSetup::get(), histoMapTH1F_, LayerTriplets::layers(), and edm::ESHandle< T >::product().

63  {
64 
65  // get the geometry
67  setup.get<MuonGeometryRecord>().get(dtGeomH);
68  dtGeom_ = dtGeomH.product();
69 
70  // Loop over all the chambers
71  if(histoMapTH1F_.empty()) {
72  for (auto ch_it : dtGeom_->chambers()) {
73  // Loop over the SLs
74  for (auto sl_it : ch_it->superLayers()) {
75  DTSuperLayerId slId = (sl_it)->id();
76  bookHistos(slId);
77  if(detailedAnalysis_) {
78  for (auto layer_it : (sl_it)->layers()) {
79  DTLayerId layerId = (layer_it)->id();
80  bookHistos(layerId);
81  }
82  }
83  }
84  }
85  }
86 }
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
void bookHistos(DTSuperLayerId slId)
std::map< DTSuperLayerId, std::vector< TH1F * > > histoMapTH1F_
const T & get() const
Definition: EventSetup.h:59
const DTGeometry * dtGeom_
T const * product() const
Definition: ESHandle.h:86
void DTResidualCalibration::bookHistos ( DTSuperLayerId  slId)
private

Definition at line 196 of file DTResidualCalibration.cc.

References histoMapTH1F_, histoMapTH2F_, histRange_, rootBaseDir_, rootFile_, DTChamberId::sector(), DTChamberId::station(), AlCaHLTBitMon_QueryRunRegistry::string, DTSuperLayerId::superlayer(), and DTChamberId::wheel().

Referenced by beginRun().

196  {
197  TH1AddDirectorySentry addDir;
198  rootFile_->cd();
199 
200  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Booking histos for SL: " << slId;
201 
202  // Compose the chamber name
203  // Define the step
204  int step = 3;
205 
206  std::string wheelStr = std::to_string(slId.wheel());
207  std::string stationStr = std::to_string(slId.station());
208  std::string sectorStr = std::to_string(slId.sector());
209 
210  std::string slHistoName =
211  "_STEP" + std::to_string(step) +
212  "_W" + wheelStr +
213  "_St" + stationStr +
214  "_Sec" + sectorStr +
215  "_SL" + std::to_string(slId.superlayer());
216 
217  edm::LogVerbatim("Calibration") << "Accessing " << rootBaseDir_;
218  TDirectory* baseDir = rootFile_->GetDirectory(rootBaseDir_.c_str());
219  if(!baseDir) baseDir = rootFile_->mkdir(rootBaseDir_.c_str());
220  edm::LogVerbatim("Calibration") << "Accessing " << ("Wheel" + wheelStr);
221  TDirectory* wheelDir = baseDir->GetDirectory(("Wheel" + wheelStr).c_str());
222  if(!wheelDir) wheelDir = baseDir->mkdir(("Wheel" + wheelStr).c_str());
223  edm::LogVerbatim("Calibration") << "Accessing " << ("Station" + stationStr);
224  TDirectory* stationDir = wheelDir->GetDirectory(("Station" + stationStr).c_str());
225  if(!stationDir) stationDir = wheelDir->mkdir(("Station" + stationStr).c_str());
226  edm::LogVerbatim("Calibration") << "Accessing " << ("Sector" + sectorStr);
227  TDirectory* sectorDir = stationDir->GetDirectory(("Sector" + sectorStr).c_str());
228  if(!sectorDir) sectorDir = stationDir->mkdir(("Sector" + sectorStr).c_str());
229 
230  sectorDir->cd();
231 
232  // Create the monitor elements
233  std::vector<TH1F*> histosTH1F;
234  histosTH1F.push_back(new TH1F(("hResDist"+slHistoName).c_str(),
235  "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
236  200, -histRange_, histRange_));
237  std::vector<TH2F*> histosTH2F;
238  histosTH2F.push_back(new TH2F(("hResDistVsDist"+slHistoName).c_str(),
239  "Residuals on the dist. (cm) from wire (rec_hit - segm_extr) vs dist. (cm)",
240  100, 0, 2.5, 200, -histRange_, histRange_));
241  histoMapTH1F_[slId] = histosTH1F;
242  histoMapTH2F_[slId] = histosTH2F;
243 }
std::map< DTSuperLayerId, std::vector< TH1F * > > histoMapTH1F_
std::map< DTSuperLayerId, std::vector< TH2F * > > histoMapTH2F_
int superlayer() const
Return the superlayer number (deprecated method name)
int sector() const
Definition: DTChamberId.h:61
step
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void DTResidualCalibration::bookHistos ( DTLayerId  slId)
private

Definition at line 245 of file DTResidualCalibration.cc.

References histoMapPerLayerTH1F_, histoMapPerLayerTH2F_, histRange_, DTLayerId::layer(), rootBaseDir_, rootFile_, DTChamberId::sector(), DTChamberId::station(), AlCaHLTBitMon_QueryRunRegistry::string, DTSuperLayerId::superlayer(), and DTChamberId::wheel().

245  {
246  TH1AddDirectorySentry addDir;
247  rootFile_->cd();
248 
249  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Booking histos for layer: " << layerId;
250 
251  // Compose the chamber name
252  std::string wheelStr = std::to_string(layerId.wheel());
253  std::string stationStr = std::to_string(layerId.station());
254  std::string sectorStr = std::to_string(layerId.sector());
255  std::string superLayerStr = std::to_string(layerId.superlayer());
256  std::string layerStr = std::to_string(layerId.layer());
257  // Define the step
258  int step = 3;
259 
260  std::string layerHistoName =
261  "_STEP" + std::to_string(step) +
262  "_W" + wheelStr +
263  "_St" + stationStr +
264  "_Sec" + sectorStr +
265  "_SL" + superLayerStr +
266  "_Layer" + layerStr;
267 
268  edm::LogVerbatim("Calibration") << "Accessing " << rootBaseDir_;
269  TDirectory* baseDir = rootFile_->GetDirectory(rootBaseDir_.c_str());
270  if(!baseDir) baseDir = rootFile_->mkdir(rootBaseDir_.c_str());
271  edm::LogVerbatim("Calibration") << "Accessing " << ("Wheel" + wheelStr);
272  TDirectory* wheelDir = baseDir->GetDirectory(("Wheel" + wheelStr).c_str());
273  if(!wheelDir) wheelDir = baseDir->mkdir(("Wheel" + wheelStr).c_str());
274  edm::LogVerbatim("Calibration") << "Accessing " << ("Station" + stationStr);
275  TDirectory* stationDir = wheelDir->GetDirectory(("Station" + stationStr).c_str());
276  if(!stationDir) stationDir = wheelDir->mkdir(("Station" + stationStr).c_str());
277  edm::LogVerbatim("Calibration") << "Accessing " << ("Sector" + sectorStr);
278  TDirectory* sectorDir = stationDir->GetDirectory(("Sector" + sectorStr).c_str());
279  if(!sectorDir) sectorDir = stationDir->mkdir(("Sector" + sectorStr).c_str());
280  edm::LogVerbatim("Calibration") << "Accessing " << ("SL" + superLayerStr);
281  TDirectory* superLayerDir = sectorDir->GetDirectory(("SL" + superLayerStr).c_str());
282  if(!superLayerDir) superLayerDir = sectorDir->mkdir(("SL" + superLayerStr).c_str());
283 
284  superLayerDir->cd();
285  // Create histograms
286  std::vector<TH1F*> histosTH1F;
287  histosTH1F.push_back(new TH1F(("hResDist"+layerHistoName).c_str(),
288  "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
289  200, -histRange_, histRange_));
290  std::vector<TH2F*> histosTH2F;
291  histosTH2F.push_back(new TH2F(("hResDistVsDist"+layerHistoName).c_str(),
292  "Residuals on the dist. (cm) from wire (rec_hit - segm_extr) vs dist. (cm)",
293  100, 0, 2.5, 200, -histRange_, histRange_));
294  histoMapPerLayerTH1F_[layerId] = histosTH1F;
295  histoMapPerLayerTH2F_[layerId] = histosTH2F;
296 }
std::map< DTLayerId, std::vector< TH1F * > > histoMapPerLayerTH1F_
std::map< DTLayerId, std::vector< TH2F * > > histoMapPerLayerTH2F_
step
void DTResidualCalibration::endJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 175 of file DTResidualCalibration.cc.

References rootFile_.

Referenced by o2olib.O2ORunMgr::executeJob().

175  {
176 
177  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Writing histos to file.";
178  rootFile_->cd();
179  rootFile_->Write();
180  rootFile_->Close();
181 
182  /*std::map<DTSuperLayerId, std::vector<TH1F*> >::const_iterator itSlHistos = histoMapTH1F_.begin();
183  std::map<DTSuperLayerId, std::vector<TH1F*> >::const_iterator itSlHistos_end = histoMapTH1F_.end();
184  for(; itSlHistos != itSlHistos_end; ++itSlHistos){
185  std::vector<TH1F*>::const_iterator itHistTH1F = (*itSlHistos).second.begin();
186  std::vector<TH1F*>::const_iterator itHistTH1F_end = (*itSlHistos).second.end();
187  for(; itHistTH1F != itHistTH1F_end; ++itHistTH1F) (*itHistTH1F)->Write();
188 
189  std::vector<TH2F*>::const_iterator itHistTH2F = histoMapTH2F_[(*itSlHistos).first].begin();
190  std::vector<TH2F*>::const_iterator itHistTH2F_end = histoMapTH2F_[(*itSlHistos).first].end();
191  for(; itHistTH2F != itHistTH2F_end; ++itHistTH2F) (*itHistTH2F)->Write();
192  }*/
193 
194 }
void DTResidualCalibration::fillHistos ( DTSuperLayerId  slId,
float  distance,
float  residualOnDistance 
)
private

Definition at line 299 of file DTResidualCalibration.cc.

References histoMapTH1F_, and histoMapTH2F_.

Referenced by analyze().

301  {
302  std::vector<TH1F*> const& histosTH1F = histoMapTH1F_[slId];
303  std::vector<TH2F*> const& histosTH2F = histoMapTH2F_[slId];
304  histosTH1F[0]->Fill(residualOnDistance);
305  histosTH2F[0]->Fill(distance, residualOnDistance);
306 }
std::map< DTSuperLayerId, std::vector< TH1F * > > histoMapTH1F_
std::map< DTSuperLayerId, std::vector< TH2F * > > histoMapTH2F_
void DTResidualCalibration::fillHistos ( DTLayerId  slId,
float  distance,
float  residualOnDistance 
)
private

Definition at line 309 of file DTResidualCalibration.cc.

References histoMapPerLayerTH1F_, and histoMapPerLayerTH2F_.

311  {
312  std::vector<TH1F*> const& histosTH1F = histoMapPerLayerTH1F_[layerId];
313  std::vector<TH2F*> const& histosTH2F = histoMapPerLayerTH2F_[layerId];
314  histosTH1F[0]->Fill(residualOnDistance);
315  histosTH2F[0]->Fill(distance, residualOnDistance);
316 }
std::map< DTLayerId, std::vector< TH1F * > > histoMapPerLayerTH1F_
std::map< DTLayerId, std::vector< TH2F * > > histoMapPerLayerTH2F_
float DTResidualCalibration::segmentToWireDistance ( const DTRecHit1D recHit1D,
const DTRecSegment4D segment 
)
private

Definition at line 148 of file DTResidualCalibration.cc.

References relativeConstraints::chamber, DTGeometry::chamber(), DTSuperLayerId::chamberId(), funct::cos(), dtGeom_, DTGeometry::layer(), DTWireId::layerId(), DTRecSegment4D::localDirection(), DTRecSegment4D::localPosition(), DTRecHit1D::localPosition(), DTLayer::specificTopology(), DTSuperLayerId::superlayer(), PV3DBase< T, PVType, FrameType >::theta(), GeomDet::toGlobal(), GeomDet::toLocal(), DTWireId::wire(), DTRecHit1D::wireId(), DTTopology::wirePosition(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by analyze().

148  {
149 
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 = segment.localPosition() + segment.localDirection()*wirePosInChamber.z()/cos(segment.localDirection().theta());
165 
166  // Compute the distance of the segment from the wire
167  int sl = wireId.superlayer();
168  float segmDistance = -1;
169  if(sl == 1 || sl == 3) segmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
170  else if(sl == 2) segmDistance = fabs(segPosAtZWire.y() - wirePosInChamber.y());
171 
172  return segmDistance;
173 }
LocalPoint localPosition() const override
Return the 3-dimensional local position.
Definition: DTRecHit1D.h:60
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:86
LocalPoint localPosition() const override
Local position in Chamber frame.
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:99
LocalVector localDirection() const override
Local direction in Chamber frame.
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:54
T y() const
Definition: PV3DBase.h:63
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:69
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
const DTTopology & specificTopology() const
Definition: DTLayer.cc:42
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int wire() const
Return the wire number.
Definition: DTWireId.h:56
int superlayer() const
Return the superlayer number (deprecated method name)
const DTGeometry * dtGeom_
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:62
const DTLayer * layer(const DTLayerId &id) const
Return a layer given its id.
Definition: DTGeometry.cc:109
T x() const
Definition: PV3DBase.h:62
DTWireId wireId() const
Return the wireId.
Definition: DTRecHit1D.h:107

Member Data Documentation

bool DTResidualCalibration::detailedAnalysis_
private

Definition at line 60 of file DTResidualCalibration.h.

Referenced by analyze(), and beginRun().

const DTGeometry* DTResidualCalibration::dtGeom_
private

Definition at line 63 of file DTResidualCalibration.h.

Referenced by analyze(), beginRun(), and segmentToWireDistance().

std::map<DTLayerId, std::vector<TH1F*> > DTResidualCalibration::histoMapPerLayerTH1F_
private

Definition at line 68 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

std::map<DTLayerId, std::vector<TH2F*> > DTResidualCalibration::histoMapPerLayerTH2F_
private

Definition at line 69 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

std::map<DTSuperLayerId, std::vector<TH1F*> > DTResidualCalibration::histoMapTH1F_
private

Definition at line 65 of file DTResidualCalibration.h.

Referenced by beginRun(), bookHistos(), and fillHistos().

std::map<DTSuperLayerId, std::vector<TH2F*> > DTResidualCalibration::histoMapTH2F_
private

Definition at line 66 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

double DTResidualCalibration::histRange_
private

Definition at line 56 of file DTResidualCalibration.h.

Referenced by bookHistos().

unsigned int DTResidualCalibration::nevent
private

Definition at line 44 of file DTResidualCalibration.h.

Referenced by analyze(), DTResidualCalibration(), and ~DTResidualCalibration().

std::string DTResidualCalibration::rootBaseDir_
private

Definition at line 58 of file DTResidualCalibration.h.

Referenced by bookHistos().

TFile* DTResidualCalibration::rootFile_
private

Definition at line 61 of file DTResidualCalibration.h.

Referenced by analyze(), bookHistos(), DTResidualCalibration(), and endJob().

unsigned int DTResidualCalibration::segmbad
private

Definition at line 45 of file DTResidualCalibration.h.

Referenced by analyze(), DTResidualCalibration(), and ~DTResidualCalibration().

edm::InputTag DTResidualCalibration::segment4DLabel_
private

Definition at line 57 of file DTResidualCalibration.h.

Referenced by analyze(), and DTResidualCalibration().

unsigned int DTResidualCalibration::segmok
private

Definition at line 45 of file DTResidualCalibration.h.

Referenced by analyze(), DTResidualCalibration(), and ~DTResidualCalibration().

DTSegmentSelector DTResidualCalibration::select_
private

Definition at line 55 of file DTResidualCalibration.h.

Referenced by analyze().