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