CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTTTrigOffsetCalibration.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * $Date: 2010/11/18 20:33:11 $
6  * $Revision: 1.9 $
7  * \author A. Vilela Pereira
8  */
9 
11 
18 
21 
27 
30 
31 #include <string>
32 #include <sstream>
33 #include "TFile.h"
34 #include "TH1F.h"
35 
36 using namespace std;
37 using namespace edm;
38 
40  select_(pset),
41  theRecHits4DLabel_(pset.getParameter<InputTag>("recHits4DLabel")),
42  doTTrigCorrection_(pset.getUntrackedParameter<bool>("doT0SegCorrection", false)),
43  theCalibChamber_(pset.getUntrackedParameter<string>("calibChamber", "All")),
44  dbLabel_(pset.getUntrackedParameter<string>("dbLabel", "")) {
45 
46  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Constructor called!";
47 
48  // the root file which will contain the histos
49  string rootFileName = pset.getUntrackedParameter<string>("rootFileName","DTT0SegHistos.root");
50  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
51  rootFile_->cd();
52 }
53 
56  ESHandle<DTTtrig> tTrig;
57  setup.get<DTTtrigRcd>().get(dbLabel_,tTrig);
58  tTrigMap_ = &*tTrig;
59  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration]: TTrig version: " << tTrig->version() << endl;
60  }
61 }
62 
64  rootFile_->Close();
65  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Destructor called!";
66 }
67 
68 void DTTTrigOffsetCalibration::analyze(const Event & event, const EventSetup& eventSetup) {
69  rootFile_->cd();
70  DTChamberId chosenChamberId;
71 
72  if(theCalibChamber_ != "All") {
73  stringstream linestr;
74  int selWheel, selStation, selSector;
75  linestr << theCalibChamber_;
76  linestr >> selWheel >> selStation >> selSector;
77  chosenChamberId = DTChamberId(selWheel, selStation, selSector);
78  LogVerbatim("Calibration") << " chosen chamber " << chosenChamberId << endl;
79  }
80 
81  // Get the DT Geometry
82  ESHandle<DTGeometry> dtGeom;
83  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
84 
85  // Get the rechit collection from the event
87  event.getByLabel(theRecHits4DLabel_, all4DSegments);
88 
89  // Loop over segments by chamber
91  for(chamberIdIt = all4DSegments->id_begin(); chamberIdIt != all4DSegments->id_end(); ++chamberIdIt){
92 
93  // Get the chamber from the setup
94  const DTChamber* chamber = dtGeom->chamber(*chamberIdIt);
95  LogTrace("Calibration") << "Chamber Id: " << *chamberIdIt;
96 
97  // Book histos
98  if(theT0SegHistoMap_.find(*chamberIdIt) == theT0SegHistoMap_.end()){
99  bookHistos(*chamberIdIt);
100  }
101 
102  // Calibrate just the chosen chamber/s
103  if((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId)) continue;
104 
105  // Get the range for the corresponding ChamberId
106  DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt));
107 
108  // Loop over the rechits of this DetUnit
109  for(DTRecSegment4DCollection::const_iterator segment = range.first;
110  segment != range.second; ++segment){
111 
112  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
113  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
114 
115  if( !select_(*segment, event, eventSetup) ) continue;
116 
117  // Fill t0-seg values
118  if( (*segment).hasPhi() ) {
119  //if( segment->phiSegment()->ist0Valid() ){
120  if( (segment->phiSegment()->t0()) != 0.00 ){
121  (theT0SegHistoMap_[*chamberIdIt])[0]->Fill(segment->phiSegment()->t0());
122  }
123  }
124  if( (*segment).hasZed() ){
125  //if( segment->zSegment()->ist0Valid() ){
126  if( (segment->zSegment()->t0()) != 0.00 ){
127  (theT0SegHistoMap_[*chamberIdIt])[1]->Fill(segment->zSegment()->t0());
128  }
129  }
130  } // DTRecSegment4DCollection::const_iterator segment
131  } // DTRecSegment4DCollection::id_iterator chamberIdIt
132 } // DTTTrigOffsetCalibration::analyze
133 
135  rootFile_->cd();
136 
137  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Writing histos to file!" << endl;
138 
139  for(ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end(); ++itChHistos){
140  for(vector<TH1F*>::const_iterator itHist = (*itChHistos).second.begin();
141  itHist != (*itChHistos).second.end(); ++itHist) (*itHist)->Write();
142  }
143 
144  if(doTTrigCorrection_){
145  // Create the object to be written to DB
146  DTTtrig* tTrig = new DTTtrig();
147 
148  for(ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end(); ++itChHistos){
149  DTChamberId chId = itChHistos->first;
150  // Get SuperLayerId's for each ChamberId
151  vector<DTSuperLayerId> slIds;
152  slIds.push_back(DTSuperLayerId(chId,1));
153  slIds.push_back(DTSuperLayerId(chId,3));
154  if(chId.station() != 4) slIds.push_back(DTSuperLayerId(chId,2));
155 
156  for(vector<DTSuperLayerId>::const_iterator itSl = slIds.begin(); itSl != slIds.end(); ++itSl){
157  // Get old values from DB
158  float ttrigMean = 0;
159  float ttrigSigma = 0;
160  float kFactor = 0;
161  tTrigMap_->get(*itSl,ttrigMean,ttrigSigma,kFactor,DTTimeUnits::ns);
162  //FIXME: verify if values make sense
163  // Set new values
164  float ttrigMeanNew = ttrigMean;
165  float ttrigSigmaNew = ttrigSigma;
166  float t0SegMean = (itSl->superLayer() != 2)?itChHistos->second[0]->GetMean():itChHistos->second[1]->GetMean();
167 
168  float kFactorNew = (kFactor*ttrigSigma+t0SegMean)/ttrigSigma;
169 
170  tTrig->set(*itSl,ttrigMeanNew,ttrigSigmaNew,kFactorNew,DTTimeUnits::ns);
171  }
172  }
173  LogVerbatim("Calibration")<< "[DTTTrigOffsetCalibration] Writing ttrig object to DB!" << endl;
174  // Write the object to DB
175  string tTrigRecord = "DTTtrigRcd";
176  DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
177  }
178 }
179 
180 // Book a set of histograms for a given Chamber
182 
183  LogTrace("Calibration") << " Booking histos for Chamber: " << chId;
184 
185  // Compose the chamber name
186  stringstream wheel; wheel << chId.wheel();
187  stringstream station; station << chId.station();
188  stringstream sector; sector << chId.sector();
189 
190  string chHistoName =
191  "_W" + wheel.str() +
192  "_St" + station.str() +
193  "_Sec" + sector.str();
194 
195  /*// Define the step
196  stringstream Step; Step << step;
197 
198  string chHistoName =
199  "_STEP" + Step.str() +
200  "_W" + wheel.str() +
201  "_St" + station.str() +
202  "_Sec" + sector.str();
203 
204  theDbe->setCurrentFolder("DT/DTCalibValidation/Wheel" + wheel.str() +
205  "/Station" + station.str() +
206  "/Sector" + sector.str());
207  // Create the monitor elements
208  vector<MonitorElement *> histos;
209  // Note hte order matters
210  histos.push_back(theDbe->book1D("hRPhiSegT0"+chHistoName, "t0 from Phi segments", 200, -25., 25.));
211  histos.push_back(theDbe->book1D("hRZSegT0"+chHistoName, "t0 from Z segments", 200, -25., 25.));*/
212 
213  vector<TH1F*> histos;
214  // Note the order matters
215  histos.push_back(new TH1F(("hRPhiSegT0"+chHistoName).c_str(), "t0 from Phi segments", 500, -60., 60.));
216  if(chId.station() != 4) histos.push_back(new TH1F(("hRZSegT0"+chHistoName).c_str(), "t0 from Z segments", 500, -60., 60.));
217 
218  theT0SegHistoMap_[chId] = histos;
219 }
T getUntrackedParameter(std::string const &, T const &) const
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:52
int set(int wheelId, int stationId, int sectorId, int slId, float tTrig, float tTrms, float kFact, DTTimeUnits::type unit)
Definition: DTTtrig.cc:262
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
void beginRun(const edm::Run &run, const edm::EventSetup &setup)
DTSuperLayerId
identifier iterator
Definition: RangeMap.h:138
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
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
#define LogTrace(id)
int get(int wheelId, int stationId, int sectorId, int slId, float &tTrig, float &tTrms, float &kFact, DTTimeUnits::type unit) const
get content
Definition: DTTtrig.cc:87
const T & get() const
Definition: EventSetup.h:55
int sector() const
Definition: DTChamberId.h:63
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup)
int station() const
Return the station number.
Definition: DTChamberId.h:53
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
static void writeToDB(std::string record, T *payload)
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
DTTTrigOffsetCalibration(const edm::ParameterSet &pset)
Definition: Run.h:33