CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTResidualCalibration Class Reference

#include <DTResidualCalibration.h>

Inheritance diagram for DTResidualCalibration:
edm::one::EDAnalyzer< edm::one::WatchRuns > edm::one::EDAnalyzerBase 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
 
void endRun (const edm::Run &, const edm::EventSetup &) override
 
 ~DTResidualCalibration () override
 Destructor. More...
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::WatchRuns >
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () 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
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex >
const & 
esGetTokenRecordIndicesVector (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::array< std::vector< ModuleDescription const * > *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, 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
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
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_
 
const edm::ESGetToken
< DTGeometry,
MuonGeometryRecord
dtGeomToken_
 
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::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Extracts DT segment residuals

Definition at line 29 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_FULL_cff::InputTag, LogDebug, nevent, rootFile_, dtT0Analyzer_cfg::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)),
39  dtGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
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 }
T getUntrackedParameter(std::string const &, T const &) const
DTSegmentSelector * select_
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
#define LogDebug(id)
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 }
Log< level::Info, true > LogVerbatim
DTSegmentSelector * select_

Member Function Documentation

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

Implements edm::one::EDAnalyzerBase.

Definition at line 88 of file DTResidualCalibration.cc.

References DTGeometry::chamber(), DTRecHitSegmentResidual::compute(), filterCSVwithJSON::copy, detailedAnalysis_, dtGeom_, edmPickEvents::event, fillHistos(), DTWireId::layerId(), LogTrace, nevent, sistrip::SpyUtilities::range(), rootFile_, segmbad, segment4DLabel_, segmentToWireDistance(), segmok, select_, GeneralSetup::setup(), 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
98  for (chamberIdIt = segments4D->id_begin(); chamberIdIt != segments4D->id_end(); ++chamberIdIt) {
99  const DTChamber* chamber = dtGeom_->chamber(*chamberIdIt);
100 
101  // Get the range for the corresponding ChamberId
102  DTRecSegment4DCollection::range range = segments4D->get((*chamberIdIt));
103  // Loop over the rechits of this DetUnit
104  for (DTRecSegment4DCollection::const_iterator segment = range.first; segment != range.second; ++segment) {
105  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
106  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
107 
108  if (!(*select_)(*segment, event, setup)) {
109  segmbad++;
110  continue;
111  }
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(); recHit1D != recHits1D_S3.end();
131  ++recHit1D) {
132  const DTWireId wireId = recHit1D->wireId();
133 
134  float segmDistance = segmentToWireDistance(*recHit1D, *segment);
135  if (segmDistance > 2.1)
136  LogTrace("Calibration") << "WARNING: segment-wire distance: " << segmDistance;
137  else
138  LogTrace("Calibration") << "segment-wire distance: " << segmDistance;
139 
140  float residualOnDistance = DTRecHitSegmentResidual().compute(dtGeom_, *recHit1D, *segment);
141  LogTrace("Calibration") << "Wire Id " << wireId << " residual on distance: " << residualOnDistance;
142 
143  fillHistos(wireId.superlayerId(), segmDistance, residualOnDistance);
144  if (detailedAnalysis_)
145  fillHistos(wireId.layerId(), segmDistance, residualOnDistance);
146  }
147  }
148  }
149 }
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
#define LogTrace(id)
const uint16_t range(const Frame &aFrame)
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
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 ( )
overridevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 62 of file DTResidualCalibration.cc.

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

Definition at line 64 of file DTResidualCalibration.cc.

References bookHistos(), DTGeometry::chambers(), detailedAnalysis_, dtGeom_, dtGeomH, dtGeomToken_, edm::EventSetup::getHandle(), histoMapTH1F_, LayerTriplets::layers(), and edm::ESHandle< class >::product().

64  {
65  // get the geometry
67  dtGeomH = setup.getHandle(dtGeomToken_);
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 }
std::map< DTSuperLayerId, TH1F * > histoMapTH1F_
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)
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
T const * product() const
Definition: ESHandle.h:86
const DTGeometry * dtGeom_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
void DTResidualCalibration::bookHistos ( DTSuperLayerId  slId)
private

Definition at line 199 of file DTResidualCalibration.cc.

References histoMapTH1F_, histoMapTH2F_, histRange_, LogDebug, rootBaseDir_, rootFile_, DTChamberId::sector(), mergeVDriftHistosByStation::sectorStr, DTChamberId::station(), mergeVDriftHistosByStation::stationStr, AlCaHLTBitMon_QueryRunRegistry::string, DTSuperLayerId::superlayer(), cond::impl::to_string(), DTChamberId::wheel(), and mergeVDriftHistosByStation::wheelStr.

Referenced by beginRun().

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

Definition at line 253 of file DTResidualCalibration.cc.

References histoMapPerLayerTH1F_, histoMapPerLayerTH2F_, histRange_, DTLayerId::layer(), LogDebug, rootBaseDir_, rootFile_, DTChamberId::sector(), mergeVDriftHistosByStation::sectorStr, DTChamberId::station(), mergeVDriftHistosByStation::stationStr, AlCaHLTBitMon_QueryRunRegistry::string, DTSuperLayerId::superlayer(), cond::impl::to_string(), DTChamberId::wheel(), and mergeVDriftHistosByStation::wheelStr.

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

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 180 of file DTResidualCalibration.cc.

References LogDebug, and rootFile_.

Referenced by o2olib.O2ORunMgr::executeJob().

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 }
#define LogDebug(id)
void DTResidualCalibration::endRun ( const edm::Run ,
const edm::EventSetup  
)
inlineoverride

Definition at line 39 of file DTResidualCalibration.h.

39 {};
void DTResidualCalibration::fillHistos ( DTSuperLayerId  slId,
float  distance,
float  residualOnDistance 
)
private

Definition at line 312 of file DTResidualCalibration.cc.

References histoMapTH1F_, and histoMapTH2F_.

Referenced by analyze().

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

Definition at line 318 of file DTResidualCalibration.cc.

References histoMapPerLayerTH1F_, and histoMapPerLayerTH2F_.

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

Definition at line 151 of file DTResidualCalibration.cc.

References DTGeometry::chamber(), DTSuperLayerId::chamberId(), funct::cos(), dtGeom_, DTGeometry::layer(), phase1PixelTopology::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().

151  {
152  // Get the layer and the wire position
153  const DTWireId wireId = recHit1D.wireId();
154  const DTLayer* layer = dtGeom_->layer(wireId);
155  float wireX = layer->specificTopology().wirePosition(wireId.wire());
156 
157  // Extrapolate the segment to the z of the wire
158  // Get wire position in chamber RF
159  // (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)
160  LocalPoint wirePosInLay(wireX, recHit1D.localPosition().y(), recHit1D.localPosition().z());
161  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
162  const DTChamber* chamber = dtGeom_->chamber(wireId.layerId().chamberId());
163  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
164 
165  // Segment position at Wire z in chamber local frame
166  LocalPoint segPosAtZWire =
167  segment.localPosition() + segment.localDirection() * wirePosInChamber.z() / cos(segment.localDirection().theta());
168 
169  // Compute the distance of the segment from the wire
170  int sl = wireId.superlayer();
171  float segmDistance = -1;
172  if (sl == 1 || sl == 3)
173  segmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
174  else if (sl == 2)
175  segmDistance = fabs(segPosAtZWire.y() - wirePosInChamber.y());
176 
177  return segmDistance;
178 }
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:59
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90
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
LocalVector localDirection() const override
Local direction in Chamber frame.
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
LocalPoint localPosition() const override
Local position in Chamber frame.
constexpr std::array< uint8_t, layerIndexSize > layer
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
LocalPoint localPosition() const override
Return the 3-dimensional local position.
Definition: DTRecHit1D.h:47
DTWireId wireId() const
Return the wireId.
Definition: DTRecHit1D.h:76

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().

const edm::ESGetToken<DTGeometry, MuonGeometryRecord> DTResidualCalibration::dtGeomToken_
private

Definition at line 64 of file DTResidualCalibration.h.

Referenced by beginRun().

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

Definition at line 70 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

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

Definition at line 71 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

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

Definition at line 67 of file DTResidualCalibration.h.

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

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

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