CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
21 
24 
25 #include <string>
26 #include <sstream>
27 #include "TFile.h"
28 #include "TH1F.h"
29 #include "TH2F.h"
30 
31 using namespace std;
32 using namespace edm;
33 
35  select_(pset),
36  theRecHits4DLabel_(pset.getParameter<InputTag>("recHits4DLabel")),
37  //writeVDriftDB_(pset.getUntrackedParameter<bool>("writeVDriftDB", false)),
38  theCalibChamber_(pset.getUntrackedParameter<string>("calibChamber", "All")) {
39 
40  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Constructor called!";
41 
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 
49  TH1::SetDefaultSumw2(true);
50 }
51 
53 
55  rootFile_->Close();
56  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Destructor called!";
57 }
58 
59 void DTVDriftSegmentCalibration::analyze(const Event & event, const EventSetup& eventSetup) {
60  rootFile_->cd();
61 
62  // Get the DT Geometry
63  ESHandle<DTGeometry> dtGeom;
64  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
65 
66  // Get the rechit collection from the event
68  event.getByLabel(theRecHits4DLabel_, all4DSegments);
69 
70  DTChamberId chosenChamberId;
71  if(theCalibChamber_ != "All") {
72  stringstream linestr;
73  int selWheel, selStation, selSector;
74  linestr << theCalibChamber_;
75  linestr >> selWheel >> selStation >> selSector;
76  chosenChamberId = DTChamberId(selWheel, selStation, selSector);
77  LogVerbatim("Calibration") << " Chosen chamber: " << chosenChamberId << endl;
78  }
79  // Loop over segments by chamber
80  DTRecSegment4DCollection::id_iterator chamberIdIt;
81  for(chamberIdIt = all4DSegments->id_begin(); chamberIdIt != all4DSegments->id_end(); ++chamberIdIt){
82 
83  // Calibrate just the chosen chamber/s
84  if((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId)) continue;
85 
86  // Book histos
87  if(theVDriftHistoMapTH1F_.find(*chamberIdIt) == theVDriftHistoMapTH1F_.end()){
88  LogTrace("Calibration") << " Booking histos for Chamber: " << *chamberIdIt;
89  bookHistos(*chamberIdIt);
90  }
91 
92  // Get the chamber from the setup
93  const DTChamber* chamber = dtGeom->chamber(*chamberIdIt);
94  // Get the range for the corresponding ChamberId
95  DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt));
96  // Loop over the rechits of this DetUnit
97  for(DTRecSegment4DCollection::const_iterator segment = range.first;
98  segment != range.second; ++segment){
99 
100 
101  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
102  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
103 
104  if( !select_(*segment, event, eventSetup) ) continue;
105 
106  // Fill v-drift values
107  if( (*segment).hasPhi() ) {
108  //if( segment->phiSegment()->ist0Valid() ){
109  double segmentVDrift = segment->phiSegment()->vDrift();
110  if( segmentVDrift != 0.00 ){
111  (theVDriftHistoMapTH1F_[*chamberIdIt])[0]->Fill(segmentVDrift);
112  (theVDriftHistoMapTH2F_[*chamberIdIt])[0]->Fill(segment->localPosition().x(),segmentVDrift);
113  (theVDriftHistoMapTH2F_[*chamberIdIt])[1]->Fill(segment->localPosition().y(),segmentVDrift);
114  }
115  }
116  // Probably not meaningful
117  if( (*segment).hasZed() ){
118  //if( segment->zSegment()->ist0Valid() ){
119  double segmentVDrift = segment->zSegment()->vDrift();
120  if( segmentVDrift != 0.00 ){
121  (theVDriftHistoMapTH1F_[*chamberIdIt])[1]->Fill(segmentVDrift);
122  }
123  }
124  } // DTRecSegment4DCollection::const_iterator segment
125  } // DTRecSegment4DCollection::id_iterator chamberIdIt
126 } // DTVDriftSegmentCalibration::analyze
127 
129  rootFile_->cd();
130 
131  LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Writing histos to file!" << endl;
132 
133  for(ChamberHistosMapTH1F::const_iterator itChHistos = theVDriftHistoMapTH1F_.begin(); itChHistos != theVDriftHistoMapTH1F_.end(); ++itChHistos){
134  vector<TH1F*>::const_iterator itHistTH1F = (*itChHistos).second.begin();
135  vector<TH1F*>::const_iterator itHistTH1F_end = (*itChHistos).second.end();
136  for(; itHistTH1F != itHistTH1F_end; ++itHistTH1F) (*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) (*itHistTH2F)->Write();
141  }
142 
143  /*if(writeVDriftDB_){
144  // ...
145  }*/
146 }
147 
148 // Book a set of histograms for a given Chamber
150 
151  // Compose the chamber name
152  stringstream wheel; wheel << chId.wheel();
153  stringstream station; station << chId.station();
154  stringstream sector; sector << chId.sector();
155 
156  string chHistoName =
157  "_W" + wheel.str() +
158  "_St" + station.str() +
159  "_Sec" + sector.str();
160 
161  vector<TH1F*> histosTH1F;
162  histosTH1F.push_back(new TH1F(("hRPhiVDriftCorr" + chHistoName).c_str(), "v-drift corr. from Phi segments", 200, -0.4, 0.4));
163  if(chId.station() != 4) histosTH1F.push_back(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(), "v-drift corr. vs. segment x position", 250, -125., 125., 200, -0.4, 0.4));
167  histosTH2F.push_back(new TH2F(("hRPhiVDriftCorrVsSegmPosY" + chHistoName).c_str(), "v-drift corr. vs. segment y position", 250, -125., 125., 200, -0.4, 0.4));
168  //if(chId.station() != 4) ...
169 
170  theVDriftHistoMapTH1F_[chId] = histosTH1F;
171  theVDriftHistoMapTH2F_[chId] = histosTH2F;
172 }
T getUntrackedParameter(std::string const &, T const &) const
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:52
void bookHistos()
Definition: Histogram.h:33
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
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
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup)
#define LogTrace(id)
const T & get() const
Definition: EventSetup.h:55
int sector() const
Definition: DTChamberId.h:61
void beginRun(const edm::Run &run, const edm::EventSetup &setup)
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
Definition: Run.h:41
DTVDriftSegmentCalibration(const edm::ParameterSet &pset)