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
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
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)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
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, TH1F * > histoMapPerLayerTH1F_
 
std::map< DTLayerId, TH2F * > histoMapPerLayerTH2F_
 
std::map< DTSuperLayerId, TH1F * > histoMapTH1F_
 
std::map< DTSuperLayerId, TH2F * > histoMapTH2F_
 
double histRange_
 
unsigned int nevent
 
std::string rootBaseDir_
 
TFile * rootFile_
 
unsigned int segmbad
 
edm::InputTag segment4DLabel_
 
unsigned int segmok
 
DTSegmentSelectorselect_
 

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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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::EDConsumerBase::consumesCollector(), edm::ParameterSet::getUntrackedParameter(), HLT_2018_cff::InputTag, LogDebug, nevent, rootFile_, CSCSkim_cfi::rootFileName, segmbad, segment4DLabel_, segmok, select_, and AlCaHLTBitMon_QueryRunRegistry::string.

35  : histRange_(pset.getParameter<double>("histogramRange")),
36  segment4DLabel_(pset.getParameter<edm::InputTag>("segment4DLabel")),
37  rootBaseDir_(pset.getUntrackedParameter<std::string>("rootBaseDir", "DT/Residuals")),
38  detailedAnalysis_(pset.getUntrackedParameter<bool>("detailedAnalysis", false)) {
40  select_ = new DTSegmentSelector(pset, collector);
41 
42  LogDebug("Calibration") << "[DTResidualCalibration] Constructor called.";
43  consumes<DTRecSegment4DCollection>(edm::InputTag(segment4DLabel_));
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 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTSegmentSelector * select_
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
DTResidualCalibration::~DTResidualCalibration ( )
override

Destructor.

Definition at line 53 of file DTResidualCalibration.cc.

References nevent, segmbad, segmok, and select_.

53  {
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 }
DTSegmentSelector * select_

Member Function Documentation

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

Definition at line 87 of file DTResidualCalibration.cc.

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

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

Reimplemented from edm::EDAnalyzer.

Definition at line 61 of file DTResidualCalibration.cc.

61 { TH1::SetDefaultSumw2(true); }
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_, hgcalTopologyTester_cfi::layers, and edm::ESHandle< T >::product().

63  {
64  // get the geometry
66  setup.get<MuonGeometryRecord>().get(dtGeomH);
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 }
std::map< DTSuperLayerId, TH1F * > histoMapTH1F_
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
void bookHistos(DTSuperLayerId slId)
const DTGeometry * dtGeom_
T get() const
Definition: EventSetup.h:73
T const * product() const
Definition: ESHandle.h:86
void DTResidualCalibration::bookHistos ( DTSuperLayerId  slId)
private

Definition at line 198 of file DTResidualCalibration.cc.

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

Referenced by beginRun().

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

Definition at line 252 of file DTResidualCalibration.cc.

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

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

Reimplemented from edm::EDAnalyzer.

Definition at line 179 of file DTResidualCalibration.cc.

References LogDebug, and rootFile_.

Referenced by o2olib.O2ORunMgr::executeJob().

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

Definition at line 311 of file DTResidualCalibration.cc.

References histoMapTH1F_, and histoMapTH2F_.

Referenced by analyze().

311  {
312  histoMapTH1F_[slId]->Fill(residualOnDistance);
313  histoMapTH2F_[slId]->Fill(distance, residualOnDistance);
314 }
std::map< DTSuperLayerId, TH1F * > histoMapTH1F_
std::map< DTSuperLayerId, TH2F * > histoMapTH2F_
void DTResidualCalibration::fillHistos ( DTLayerId  slId,
float  distance,
float  residualOnDistance 
)
private

Definition at line 317 of file DTResidualCalibration.cc.

References histoMapPerLayerTH1F_, and histoMapPerLayerTH2F_.

317  {
318  histoMapPerLayerTH1F_[layerId]->Fill(residualOnDistance);
319  histoMapPerLayerTH2F_[layerId]->Fill(distance, residualOnDistance);
320 }
std::map< DTLayerId, TH2F * > histoMapPerLayerTH2F_
std::map< DTLayerId, TH1F * > histoMapPerLayerTH1F_
float DTResidualCalibration::segmentToWireDistance ( const DTRecHit1D recHit1D,
const DTRecSegment4D segment 
)
private

Definition at line 150 of file DTResidualCalibration.cc.

References relativeConstraints::chamber, DTGeometry::chamber(), DTSuperLayerId::chamberId(), funct::cos(), dtGeom_, DTGeometry::layer(), DTWireId::layerId(), DTRecSegment4D::localDirection(), DTRecHit1D::localPosition(), DTRecSegment4D::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().

150  {
151  // Get the layer and the wire position
152  const DTWireId wireId = recHit1D.wireId();
153  const DTLayer* layer = dtGeom_->layer(wireId);
154  float wireX = layer->specificTopology().wirePosition(wireId.wire());
155 
156  // Extrapolate the segment to the z of the wire
157  // Get wire position in chamber RF
158  // (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)
159  LocalPoint wirePosInLay(wireX, recHit1D.localPosition().y(), recHit1D.localPosition().z());
160  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
161  const DTChamber* chamber = dtGeom_->chamber(wireId.layerId().chamberId());
162  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
163 
164  // Segment position at Wire z in chamber local frame
165  LocalPoint segPosAtZWire =
166  segment.localPosition() + segment.localDirection() * wirePosInChamber.z() / cos(segment.localDirection().theta());
167 
168  // Compute the distance of the segment from the wire
169  int sl = wireId.superlayer();
170  float segmDistance = -1;
171  if (sl == 1 || sl == 3)
172  segmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
173  else if (sl == 2)
174  segmDistance = fabs(segPosAtZWire.y() - wirePosInChamber.y());
175 
176  return segmDistance;
177 }
LocalPoint localPosition() const override
Return the 3-dimensional local position.
Definition: DTRecHit1D.h:47
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:59
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:90
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:49
T y() const
Definition: PV3DBase.h:60
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
const DTTopology & specificTopology() const
Definition: DTLayer.cc:37
T z() const
Definition: PV3DBase.h:61
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int wire() const
Return the wire number.
Definition: DTWireId.h:42
int superlayer() const
Return the superlayer number (deprecated method name)
const DTGeometry * dtGeom_
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:45
const DTLayer * layer(const DTLayerId &id) const
Return a layer given its id.
Definition: DTGeometry.cc:96
T x() const
Definition: PV3DBase.h:59
DTWireId wireId() const
Return the wireId.
Definition: DTRecHit1D.h:76

Member Data Documentation

bool DTResidualCalibration::detailedAnalysis_
private

Definition at line 58 of file DTResidualCalibration.h.

Referenced by analyze(), and beginRun().

const DTGeometry* DTResidualCalibration::dtGeom_
private

Definition at line 61 of file DTResidualCalibration.h.

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

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

Definition at line 66 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

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

Definition at line 67 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

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

Definition at line 63 of file DTResidualCalibration.h.

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

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

Definition at line 64 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

double DTResidualCalibration::histRange_
private

Definition at line 54 of file DTResidualCalibration.h.

Referenced by bookHistos().

unsigned int DTResidualCalibration::nevent
private

Definition at line 42 of file DTResidualCalibration.h.

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

std::string DTResidualCalibration::rootBaseDir_
private

Definition at line 56 of file DTResidualCalibration.h.

Referenced by bookHistos().

TFile* DTResidualCalibration::rootFile_
private

Definition at line 59 of file DTResidualCalibration.h.

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

unsigned int DTResidualCalibration::segmbad
private

Definition at line 43 of file DTResidualCalibration.h.

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

edm::InputTag DTResidualCalibration::segment4DLabel_
private

Definition at line 55 of file DTResidualCalibration.h.

Referenced by analyze(), and DTResidualCalibration().

unsigned int DTResidualCalibration::segmok
private

Definition at line 43 of file DTResidualCalibration.h.

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

DTSegmentSelector* DTResidualCalibration::select_
private

Definition at line 53 of file DTResidualCalibration.h.

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