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