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  dtGeomToken_(esConsumes()) {
38  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Constructor called!";
39 
41  select_ = new DTSegmentSelector(pset, collector);
42  consumes<DTRecSegment4DCollection>(edm::InputTag(theRecHits4DLabel_));
43  // the root file which will contain the histos
44  string rootFileName = pset.getUntrackedParameter<string>("rootFileName", "DTVDriftHistos.root");
45  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
46  rootFile_->cd();
47 }
48 
49 void DTVDriftSegmentCalibration::beginJob() { TH1::SetDefaultSumw2(true); }
50 
52 
54  rootFile_->Close();
55  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Destructor called!";
56 }
57 
59  rootFile_->cd();
60 
61  // Get the DT Geometry
62  ESHandle<DTGeometry> dtGeom;
63  dtGeom = eventSetup.getHandle(dtGeomToken_);
64 
65  // Get the rechit collection from the event
67  event.getByLabel(theRecHits4DLabel_, all4DSegments);
68 
69  DTChamberId chosenChamberId;
70  if (theCalibChamber_ != "All") {
71  stringstream linestr;
72  int selWheel, selStation, selSector;
73  linestr << theCalibChamber_;
74  linestr >> selWheel >> selStation >> selSector;
75  chosenChamberId = DTChamberId(selWheel, selStation, selSector);
76  LogVerbatim("Calibration") << " Chosen chamber: " << chosenChamberId << endl;
77  }
78  // Loop over segments by chamber
80  for (chamberIdIt = all4DSegments->id_begin(); chamberIdIt != all4DSegments->id_end(); ++chamberIdIt) {
81  // Calibrate just the chosen chamber/s
82  if ((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId))
83  continue;
84 
85  // Book histos
86  if (theVDriftHistoMapTH1F_.find(*chamberIdIt) == theVDriftHistoMapTH1F_.end()) {
87  LogTrace("Calibration") << " Booking histos for Chamber: " << *chamberIdIt;
88  bookHistos(*chamberIdIt);
89  }
90 
91  // Get the chamber from the setup
92  const DTChamber* chamber = dtGeom->chamber(*chamberIdIt);
93  // Get the range for the corresponding ChamberId
94  DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt));
95  // Loop over the rechits of this DetUnit
96  for (DTRecSegment4DCollection::const_iterator segment = range.first; segment != range.second; ++segment) {
97  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
98  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
99 
100  if (!(*select_)(*segment, event, eventSetup))
101  continue;
102 
103  // Fill v-drift values
104  if ((*segment).hasPhi()) {
105  //if( segment->phiSegment()->ist0Valid() ){
106  double segmentVDrift = segment->phiSegment()->vDrift();
107  if (segmentVDrift != 0.00) {
108  (theVDriftHistoMapTH1F_[*chamberIdIt])[0]->Fill(segmentVDrift);
109  (theVDriftHistoMapTH2F_[*chamberIdIt])[0]->Fill(segment->localPosition().x(), segmentVDrift);
110  (theVDriftHistoMapTH2F_[*chamberIdIt])[1]->Fill(segment->localPosition().y(), segmentVDrift);
111  }
112  }
113  // Probably not meaningful
114  if ((*segment).hasZed()) {
115  //if( segment->zSegment()->ist0Valid() ){
116  double segmentVDrift = segment->zSegment()->vDrift();
117  if (segmentVDrift != 0.00) {
118  (theVDriftHistoMapTH1F_[*chamberIdIt])[1]->Fill(segmentVDrift);
119  }
120  }
121  } // DTRecSegment4DCollection::const_iterator segment
122  } // DTRecSegment4DCollection::id_iterator chamberIdIt
123 } // DTVDriftSegmentCalibration::analyze
124 
126  rootFile_->cd();
127 
128  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Writing histos to file!" << endl;
129 
130  for (ChamberHistosMapTH1F::const_iterator itChHistos = theVDriftHistoMapTH1F_.begin();
131  itChHistos != theVDriftHistoMapTH1F_.end();
132  ++itChHistos) {
133  vector<TH1F*>::const_iterator itHistTH1F = (*itChHistos).second.begin();
134  vector<TH1F*>::const_iterator itHistTH1F_end = (*itChHistos).second.end();
135  for (; itHistTH1F != itHistTH1F_end; ++itHistTH1F)
136  (*itHistTH1F)->Write();
137 
138  vector<TH2F*>::const_iterator itHistTH2F = theVDriftHistoMapTH2F_[(*itChHistos).first].begin();
139  vector<TH2F*>::const_iterator itHistTH2F_end = theVDriftHistoMapTH2F_[(*itChHistos).first].end();
140  for (; itHistTH2F != itHistTH2F_end; ++itHistTH2F)
141  (*itHistTH2F)->Write();
142  }
143 
144  /*if(writeVDriftDB_){
145  // ...
146  }*/
147 }
148 
149 // Book a set of histograms for a given Chamber
151  // Compose the chamber name
154  std::string sector = std::to_string(chId.sector());
155 
156  string chHistoName = "_W" + wheel + "_St" + station + "_Sec" + sector;
157 
158  vector<TH1F*> histosTH1F;
159  histosTH1F.push_back(
160  new TH1F(("hRPhiVDriftCorr" + chHistoName).c_str(), "v-drift corr. from Phi segments", 200, -0.4, 0.4));
161  if (chId.station() != 4)
162  histosTH1F.push_back(
163  new TH1F(("hRZVDriftCorr" + chHistoName).c_str(), "v-drift corr. from Z segments", 200, -0.4, 0.4));
164 
165  vector<TH2F*> histosTH2F;
166  histosTH2F.push_back(new TH2F(("hRPhiVDriftCorrVsSegmPosX" + chHistoName).c_str(),
167  "v-drift corr. vs. segment x position",
168  250,
169  -125.,
170  125.,
171  200,
172  -0.4,
173  0.4));
174  histosTH2F.push_back(new TH2F(("hRPhiVDriftCorrVsSegmPosY" + chHistoName).c_str(),
175  "v-drift corr. vs. segment y position",
176  250,
177  -125.,
178  125.,
179  200,
180  -0.4,
181  0.4));
182  //if(chId.station() != 4) ...
183 
184  theVDriftHistoMapTH1F_[chId] = histosTH1F;
185  theVDriftHistoMapTH2F_[chId] = histosTH2F;
186 }
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:71
#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)
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)