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(), LogDebug, nevent, rootFile_, DTAnalyzerDetailed_cfi::rootFileName, segmbad, segment4DLabel_, segmok, select_, and AlCaHLTBitMon_QueryRunRegistry::string.

34  :
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)) {
39 
41  select_ = new DTSegmentSelector(pset,collector);
42 
43  LogDebug("Calibration") << "[DTResidualCalibration] Constructor called.";
44  consumes< DTRecSegment4DCollection >(edm::InputTag(segment4DLabel_));
45  std::string rootFileName = pset.getUntrackedParameter<std::string>("rootFileName","residuals.root");
46  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
47  rootFile_->cd();
48 
49  segmok=0;
50  segmbad=0;
51  nevent=0;
52 }
#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 54 of file DTResidualCalibration.cc.

References nevent, segmbad, segmok, and select_.

54  {
55  delete select_;
56  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Destructor called.";
57  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Analyzed events: " << nevent;
58  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Good segments: " << segmok;
59  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Bad segments: " << segmbad;
60 }
DTSegmentSelector * select_

Member Function Documentation

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

Definition at line 91 of file DTResidualCalibration.cc.

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

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

Reimplemented from edm::EDAnalyzer.

Definition at line 62 of file DTResidualCalibration.cc.

62  {
63  TH1::SetDefaultSumw2(true);
64 }
void DTResidualCalibration::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)
override

Definition at line 66 of file DTResidualCalibration.cc.

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

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

Definition at line 200 of file DTResidualCalibration.cc.

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

Referenced by beginRun().

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

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

247  {
248  TH1AddDirectorySentry addDir;
249  rootFile_->cd();
250 
251  LogDebug("Calibration") << "[DTResidualCalibration] Booking histos for layer: " << layerId;
252 
253  // Compose the chamber name
254  std::string wheelStr = std::to_string(layerId.wheel());
255  std::string stationStr = std::to_string(layerId.station());
256  std::string sectorStr = std::to_string(layerId.sector());
257  std::string superLayerStr = std::to_string(layerId.superlayer());
258  std::string layerStr = std::to_string(layerId.layer());
259  // Define the step
260  int step = 3;
261 
262  std::string layerHistoName =
263  "_STEP" + std::to_string(step) +
264  "_W" + wheelStr +
265  "_St" + stationStr +
266  "_Sec" + sectorStr +
267  "_SL" + superLayerStr +
268  "_Layer" + layerStr;
269 
270  LogDebug("Calibration") << "Accessing " << rootBaseDir_;
271  TDirectory* baseDir = rootFile_->GetDirectory(rootBaseDir_.c_str());
272  if(!baseDir) baseDir = rootFile_->mkdir(rootBaseDir_.c_str());
273  LogDebug("Calibration") << "Accessing " << ("Wheel" + wheelStr);
274  TDirectory* wheelDir = baseDir->GetDirectory(("Wheel" + wheelStr).c_str());
275  if(!wheelDir) wheelDir = baseDir->mkdir(("Wheel" + wheelStr).c_str());
276  LogDebug("Calibration") << "Accessing " << ("Station" + stationStr);
277  TDirectory* stationDir = wheelDir->GetDirectory(("Station" + stationStr).c_str());
278  if(!stationDir) stationDir = wheelDir->mkdir(("Station" + stationStr).c_str());
279  LogDebug("Calibration") << "Accessing " << ("Sector" + sectorStr);
280  TDirectory* sectorDir = stationDir->GetDirectory(("Sector" + sectorStr).c_str());
281  if(!sectorDir) sectorDir = stationDir->mkdir(("Sector" + sectorStr).c_str());
282  LogDebug("Calibration") << "Accessing " << ("SL" + superLayerStr);
283  TDirectory* superLayerDir = sectorDir->GetDirectory(("SL" + superLayerStr).c_str());
284  if(!superLayerDir) superLayerDir = sectorDir->mkdir(("SL" + superLayerStr).c_str());
285 
286  superLayerDir->cd();
287  // Create histograms
288  TH1F* histosTH1F = new TH1F(("hResDist"+layerHistoName).c_str(),
289  "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
290  200, -histRange_, histRange_);
291  TH2F* histosTH2F = 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 }
#define LogDebug(id)
std::map< DTLayerId, TH1F * > histoMapPerLayerTH1F_
std::map< DTLayerId, TH2F * > histoMapPerLayerTH2F_
step
Definition: StallMonitor.cc:94
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 
181  LogDebug("Calibration") << "[DTResidualCalibration] Writing histos to file.";
182  rootFile_->cd();
183  rootFile_->Write();
184  rootFile_->Close();
185 
186  /*std::map<DTSuperLayerId, TH1F* >::const_iterator itSlHistos = histoMapTH1F_.begin();
187  std::map<DTSuperLayerId, TH1F* >::const_iterator itSlHistos_end = histoMapTH1F_.end();
188  for(; itSlHistos != itSlHistos_end; ++itSlHistos){
189  std::vector<TH1F*>::const_iterator itHistTH1F = (*itSlHistos).second.begin();
190  std::vector<TH1F*>::const_iterator itHistTH1F_end = (*itSlHistos).second.end();
191  for(; itHistTH1F != itHistTH1F_end; ++itHistTH1F) (*itHistTH1F)->Write();
192 
193  std::vector<TH2F*>::const_iterator itHistTH2F = histoMapTH2F_[(*itSlHistos).first].begin();
194  std::vector<TH2F*>::const_iterator itHistTH2F_end = histoMapTH2F_[(*itSlHistos).first].end();
195  for(; itHistTH2F != itHistTH2F_end; ++itHistTH2F) (*itHistTH2F)->Write();
196  }*/
197 
198 }
#define LogDebug(id)
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  histoMapTH1F_[slId]->Fill(residualOnDistance);
303  histoMapTH2F_[slId]->Fill(distance, residualOnDistance);
304 }
std::map< DTSuperLayerId, TH1F * > histoMapTH1F_
std::map< DTSuperLayerId, TH2F * > histoMapTH2F_
void DTResidualCalibration::fillHistos ( DTLayerId  slId,
float  distance,
float  residualOnDistance 
)
private

Definition at line 307 of file DTResidualCalibration.cc.

References histoMapPerLayerTH1F_, and histoMapPerLayerTH2F_.

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

Definition at line 152 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().

152  {
153 
154  // Get the layer and the wire position
155  const DTWireId wireId = recHit1D.wireId();
156  const DTLayer* layer = dtGeom_->layer(wireId);
157  float wireX = layer->specificTopology().wirePosition(wireId.wire());
158 
159  // Extrapolate the segment to the z of the wire
160  // Get wire position in chamber RF
161  // (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)
162  LocalPoint wirePosInLay(wireX,recHit1D.localPosition().y(),recHit1D.localPosition().z());
163  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
164  const DTChamber* chamber = dtGeom_->chamber(wireId.layerId().chamberId());
165  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
166 
167  // Segment position at Wire z in chamber local frame
168  LocalPoint segPosAtZWire = segment.localPosition() + segment.localDirection()*wirePosInChamber.z()/cos(segment.localDirection().theta());
169 
170  // Compute the distance of the segment from the wire
171  int sl = wireId.superlayer();
172  float segmDistance = -1;
173  if(sl == 1 || sl == 3) segmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
174  else if(sl == 2) segmDistance = fabs(segPosAtZWire.y() - wirePosInChamber.y());
175 
176  return segmDistance;
177 }
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:117
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:127
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, TH1F* > DTResidualCalibration::histoMapPerLayerTH1F_
private

Definition at line 68 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

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

Definition at line 69 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

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

Definition at line 65 of file DTResidualCalibration.h.

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

std::map<DTSuperLayerId, 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(), DTResidualCalibration(), and ~DTResidualCalibration().