CMS 3D CMS Logo

DTVDriftSegmentCalibration.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * \author A. Vilela Pereira
6  */
7 
9 
16 
20 
23 
24 #include <string>
25 #include <sstream>
26 #include "TFile.h"
27 #include "TH1F.h"
28 #include "TH2F.h"
29 
30 using namespace std;
31 using namespace edm;
32 
34  : theRecHits4DLabel_(pset.getParameter<InputTag>("recHits4DLabel")),
35  //writeVDriftDB_(pset.getUntrackedParameter<bool>("writeVDriftDB", false)),
36  theCalibChamber_(pset.getUntrackedParameter<string>("calibChamber", "All")) {
37  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Constructor called!";
38 
40  select_ = new DTSegmentSelector(pset, collector);
41  consumes<DTRecSegment4DCollection>(edm::InputTag(theRecHits4DLabel_));
42  // the root file which will contain the histos
43  string rootFileName = pset.getUntrackedParameter<string>("rootFileName", "DTVDriftHistos.root");
44  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
45  rootFile_->cd();
46 }
47 
48 void DTVDriftSegmentCalibration::beginJob() { TH1::SetDefaultSumw2(true); }
49 
51 
53  rootFile_->Close();
54  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Destructor called!";
55 }
56 
57 void DTVDriftSegmentCalibration::analyze(const Event& event, const EventSetup& eventSetup) {
58  rootFile_->cd();
59 
60  // Get the DT Geometry
61  ESHandle<DTGeometry> dtGeom;
62  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
63 
64  // Get the rechit collection from the event
66  event.getByLabel(theRecHits4DLabel_, all4DSegments);
67 
68  DTChamberId chosenChamberId;
69  if (theCalibChamber_ != "All") {
70  stringstream linestr;
71  int selWheel, selStation, selSector;
72  linestr << theCalibChamber_;
73  linestr >> selWheel >> selStation >> selSector;
74  chosenChamberId = DTChamberId(selWheel, selStation, selSector);
75  LogVerbatim("Calibration") << " Chosen chamber: " << chosenChamberId << endl;
76  }
77  // Loop over segments by chamber
79  for (chamberIdIt = all4DSegments->id_begin(); chamberIdIt != all4DSegments->id_end(); ++chamberIdIt) {
80  // Calibrate just the chosen chamber/s
81  if ((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId))
82  continue;
83 
84  // Book histos
85  if (theVDriftHistoMapTH1F_.find(*chamberIdIt) == theVDriftHistoMapTH1F_.end()) {
86  LogTrace("Calibration") << " Booking histos for Chamber: " << *chamberIdIt;
87  bookHistos(*chamberIdIt);
88  }
89 
90  // Get the chamber from the setup
91  const DTChamber* chamber = dtGeom->chamber(*chamberIdIt);
92  // Get the range for the corresponding ChamberId
93  DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt));
94  // Loop over the rechits of this DetUnit
95  for (DTRecSegment4DCollection::const_iterator segment = range.first; segment != range.second; ++segment) {
96  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
97  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
98 
99  if (!(*select_)(*segment, event, eventSetup))
100  continue;
101 
102  // Fill v-drift values
103  if ((*segment).hasPhi()) {
104  //if( segment->phiSegment()->ist0Valid() ){
105  double segmentVDrift = segment->phiSegment()->vDrift();
106  if (segmentVDrift != 0.00) {
107  (theVDriftHistoMapTH1F_[*chamberIdIt])[0]->Fill(segmentVDrift);
108  (theVDriftHistoMapTH2F_[*chamberIdIt])[0]->Fill(segment->localPosition().x(), segmentVDrift);
109  (theVDriftHistoMapTH2F_[*chamberIdIt])[1]->Fill(segment->localPosition().y(), segmentVDrift);
110  }
111  }
112  // Probably not meaningful
113  if ((*segment).hasZed()) {
114  //if( segment->zSegment()->ist0Valid() ){
115  double segmentVDrift = segment->zSegment()->vDrift();
116  if (segmentVDrift != 0.00) {
117  (theVDriftHistoMapTH1F_[*chamberIdIt])[1]->Fill(segmentVDrift);
118  }
119  }
120  } // DTRecSegment4DCollection::const_iterator segment
121  } // DTRecSegment4DCollection::id_iterator chamberIdIt
122 } // DTVDriftSegmentCalibration::analyze
123 
125  rootFile_->cd();
126 
127  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Writing histos to file!" << endl;
128 
129  for (ChamberHistosMapTH1F::const_iterator itChHistos = theVDriftHistoMapTH1F_.begin();
130  itChHistos != theVDriftHistoMapTH1F_.end();
131  ++itChHistos) {
132  vector<TH1F*>::const_iterator itHistTH1F = (*itChHistos).second.begin();
133  vector<TH1F*>::const_iterator itHistTH1F_end = (*itChHistos).second.end();
134  for (; itHistTH1F != itHistTH1F_end; ++itHistTH1F)
135  (*itHistTH1F)->Write();
136 
137  vector<TH2F*>::const_iterator itHistTH2F = theVDriftHistoMapTH2F_[(*itChHistos).first].begin();
138  vector<TH2F*>::const_iterator itHistTH2F_end = theVDriftHistoMapTH2F_[(*itChHistos).first].end();
139  for (; itHistTH2F != itHistTH2F_end; ++itHistTH2F)
140  (*itHistTH2F)->Write();
141  }
142 
143  /*if(writeVDriftDB_){
144  // ...
145  }*/
146 }
147 
148 // Book a set of histograms for a given Chamber
150  // Compose the chamber name
151  std::string wheel = std::to_string(chId.wheel());
152  std::string station = std::to_string(chId.station());
153  std::string sector = std::to_string(chId.sector());
154 
155  string chHistoName = "_W" + wheel + "_St" + station + "_Sec" + sector;
156 
157  vector<TH1F*> histosTH1F;
158  histosTH1F.push_back(
159  new TH1F(("hRPhiVDriftCorr" + chHistoName).c_str(), "v-drift corr. from Phi segments", 200, -0.4, 0.4));
160  if (chId.station() != 4)
161  histosTH1F.push_back(
162  new TH1F(("hRZVDriftCorr" + chHistoName).c_str(), "v-drift corr. from Z segments", 200, -0.4, 0.4));
163 
164  vector<TH2F*> histosTH2F;
165  histosTH2F.push_back(new TH2F(("hRPhiVDriftCorrVsSegmPosX" + chHistoName).c_str(),
166  "v-drift corr. vs. segment x position",
167  250,
168  -125.,
169  125.,
170  200,
171  -0.4,
172  0.4));
173  histosTH2F.push_back(new TH2F(("hRPhiVDriftCorrVsSegmPosY" + chHistoName).c_str(),
174  "v-drift corr. vs. segment y position",
175  250,
176  -125.,
177  125.,
178  200,
179  -0.4,
180  0.4));
181  //if(chId.station() != 4) ...
182 
183  theVDriftHistoMapTH1F_[chId] = histosTH1F;
184  theVDriftHistoMapTH2F_[chId] = histosTH2F;
185 }
DTVDriftSegmentCalibration.h
FastTimerService_cff.range
range
Definition: FastTimerService_cff.py:34
DTVDriftSegmentCalibration::rootFile_
TFile * rootFile_
Definition: DTVDriftSegmentCalibration.h:45
DTVDriftSegmentCalibration::bookHistos
void bookHistos(DTChamberId)
Definition: DTVDriftSegmentCalibration.cc:149
MessageLogger.h
ESHandle.h
DTVDriftSegmentCalibration::theVDriftHistoMapTH2F_
ChamberHistosMapTH2F theVDriftHistoMapTH2F_
Definition: DTVDriftSegmentCalibration.h:47
DTVDriftSegmentCalibration::beginJob
void beginJob() override
Definition: DTVDriftSegmentCalibration.cc:48
edm::Run
Definition: Run.h:45
relativeConstraints.station
station
Definition: relativeConstraints.py:67
edm
HLT enums.
Definition: AlignableModifier.h:19
DTChamber
Definition: DTChamber.h:24
DTVDriftSegmentCalibration::theCalibChamber_
std::string theCalibChamber_
Definition: DTVDriftSegmentCalibration.h:43
DTVDriftSegmentCalibration::analyze
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
Definition: DTVDriftSegmentCalibration.cc:57
edm::EDConsumerBase::consumesCollector
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
Definition: EDConsumerBase.cc:47
edm::Handle< DTRecSegment4DCollection >
edm::RangeMap::id_iterator
identifier iterator
Definition: RangeMap.h:130
singleTopDQM_cfi.setup
setup
Definition: singleTopDQM_cfi.py:37
DTGeometry::chamber
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90
edm::EventSetup::get
T get() const
Definition: EventSetup.h:87
DTVDriftSegmentCalibration::theVDriftHistoMapTH1F_
ChamberHistosMapTH1F theVDriftHistoMapTH1F_
Definition: DTVDriftSegmentCalibration.h:46
edm::ESHandle< DTGeometry >
DTVDriftSegmentCalibration::~DTVDriftSegmentCalibration
~DTVDriftSegmentCalibration() override
Definition: DTVDriftSegmentCalibration.cc:52
DTSegmentSelector
Definition: DTSegmentSelector.h:24
DTGeometry.h
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
edm::RangeMap::const_iterator
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
DTVDriftSegmentCalibration::endJob
void endJob() override
Definition: DTVDriftSegmentCalibration.cc:124
edmPickEvents.event
event
Definition: edmPickEvents.py:273
makeMuonMisalignmentScenario.wheel
wheel
Definition: makeMuonMisalignmentScenario.py:319
DTCalibDBUtils.h
edm::EventSetup
Definition: EventSetup.h:58
get
#define get
DTVDriftSegmentCalibration::theRecHits4DLabel_
edm::InputTag theRecHits4DLabel_
Definition: DTVDriftSegmentCalibration.h:41
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HcalObjRepresent::Fill
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
Definition: HcalObjRepresent.h:1053
InputTag.h
DTChamberId::sector
int sector() const
Definition: DTChamberId.h:49
edm::RangeMap::range
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
std
Definition: JetResolutionObject.h:76
writedatasetfile.run
run
Definition: writedatasetfile.py:27
DTSegmentSelector.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
relativeConstraints.chamber
chamber
Definition: relativeConstraints.py:53
EventSetup.h
CSCSkim_cfi.rootFileName
rootFileName
Definition: CSCSkim_cfi.py:9
DTVDriftSegmentCalibration::beginRun
void beginRun(const edm::Run &run, const edm::EventSetup &setup) override
Definition: DTVDriftSegmentCalibration.cc:50
DTVDriftSegmentCalibration::select_
DTSegmentSelector * select_
Definition: DTVDriftSegmentCalibration.h:39
DTChamberId
Definition: DTChamberId.h:14
LogTrace
#define LogTrace(id)
Definition: MessageLogger.h:234
ParameterSet.h
MuonGeometryRecord.h
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
DTVDriftSegmentCalibration::DTVDriftSegmentCalibration
DTVDriftSegmentCalibration(const edm::ParameterSet &pset)
Definition: DTVDriftSegmentCalibration.cc:33
MuonGeometryRecord
Definition: MuonGeometryRecord.h:34
DTChamberId::wheel
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
edm::InputTag
Definition: InputTag.h:15
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
DTChamberId::station
int station() const
Return the station number.
Definition: DTChamberId.h:42
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
DTRecSegment4DCollection.h