Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "DTTTrigOffsetCalibration.h"
00011
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/Utilities/interface/InputTag.h"
00018
00019 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00020 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00021
00022 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00023 #include "CondFormats/DTObjects/interface/DTTtrig.h"
00024 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
00025 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00026 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00027
00028 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
00029 #include "CalibMuon/DTCalibration/interface/DTSegmentSelector.h"
00030
00031 #include <string>
00032 #include <sstream>
00033 #include "TFile.h"
00034 #include "TH1F.h"
00035
00036 using namespace std;
00037 using namespace edm;
00038
00039 DTTTrigOffsetCalibration::DTTTrigOffsetCalibration(const ParameterSet& pset):
00040 select_(pset),
00041 theRecHits4DLabel_(pset.getParameter<InputTag>("recHits4DLabel")),
00042 doTTrigCorrection_(pset.getUntrackedParameter<bool>("doT0SegCorrection", false)),
00043 theCalibChamber_(pset.getUntrackedParameter<string>("calibChamber", "All")),
00044 dbLabel_(pset.getUntrackedParameter<string>("dbLabel", "")) {
00045
00046 LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Constructor called!";
00047
00048
00049 string rootFileName = pset.getUntrackedParameter<string>("rootFileName","DTT0SegHistos.root");
00050 rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
00051 rootFile_->cd();
00052 }
00053
00054 void DTTTrigOffsetCalibration::beginRun(const edm::Run& run, const edm::EventSetup& setup) {
00055 if(doTTrigCorrection_){
00056 ESHandle<DTTtrig> tTrig;
00057 setup.get<DTTtrigRcd>().get(dbLabel_,tTrig);
00058 tTrigMap_ = &*tTrig;
00059 LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration]: TTrig version: " << tTrig->version() << endl;
00060 }
00061 }
00062
00063 DTTTrigOffsetCalibration::~DTTTrigOffsetCalibration(){
00064 rootFile_->Close();
00065 LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Destructor called!";
00066 }
00067
00068 void DTTTrigOffsetCalibration::analyze(const Event & event, const EventSetup& eventSetup) {
00069 rootFile_->cd();
00070 DTChamberId chosenChamberId;
00071
00072 if(theCalibChamber_ != "All") {
00073 stringstream linestr;
00074 int selWheel, selStation, selSector;
00075 linestr << theCalibChamber_;
00076 linestr >> selWheel >> selStation >> selSector;
00077 chosenChamberId = DTChamberId(selWheel, selStation, selSector);
00078 LogVerbatim("Calibration") << " chosen chamber " << chosenChamberId << endl;
00079 }
00080
00081
00082 ESHandle<DTGeometry> dtGeom;
00083 eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00084
00085
00086 Handle<DTRecSegment4DCollection> all4DSegments;
00087 event.getByLabel(theRecHits4DLabel_, all4DSegments);
00088
00089
00090 DTRecSegment4DCollection::id_iterator chamberIdIt;
00091 for(chamberIdIt = all4DSegments->id_begin(); chamberIdIt != all4DSegments->id_end(); ++chamberIdIt){
00092
00093
00094 const DTChamber* chamber = dtGeom->chamber(*chamberIdIt);
00095 LogTrace("Calibration") << "Chamber Id: " << *chamberIdIt;
00096
00097
00098 if(theT0SegHistoMap_.find(*chamberIdIt) == theT0SegHistoMap_.end()){
00099 bookHistos(*chamberIdIt);
00100 }
00101
00102
00103 if((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId)) continue;
00104
00105
00106 DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt));
00107
00108
00109 for(DTRecSegment4DCollection::const_iterator segment = range.first;
00110 segment != range.second; ++segment){
00111
00112 LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
00113 << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
00114
00115 if( !select_(*segment, event, eventSetup) ) continue;
00116
00117
00118 if( (*segment).hasPhi() ) {
00119
00120 if( (segment->phiSegment()->t0()) != 0.00 ){
00121 (theT0SegHistoMap_[*chamberIdIt])[0]->Fill(segment->phiSegment()->t0());
00122 }
00123 }
00124 if( (*segment).hasZed() ){
00125
00126 if( (segment->zSegment()->t0()) != 0.00 ){
00127 (theT0SegHistoMap_[*chamberIdIt])[1]->Fill(segment->zSegment()->t0());
00128 }
00129 }
00130 }
00131 }
00132 }
00133
00134 void DTTTrigOffsetCalibration::endJob() {
00135 rootFile_->cd();
00136
00137 LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Writing histos to file!" << endl;
00138
00139 for(ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end(); ++itChHistos){
00140 for(vector<TH1F*>::const_iterator itHist = (*itChHistos).second.begin();
00141 itHist != (*itChHistos).second.end(); ++itHist) (*itHist)->Write();
00142 }
00143
00144 if(doTTrigCorrection_){
00145
00146 DTTtrig* tTrig = new DTTtrig();
00147
00148 for(ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end(); ++itChHistos){
00149 DTChamberId chId = itChHistos->first;
00150
00151 vector<DTSuperLayerId> slIds;
00152 slIds.push_back(DTSuperLayerId(chId,1));
00153 slIds.push_back(DTSuperLayerId(chId,3));
00154 if(chId.station() != 4) slIds.push_back(DTSuperLayerId(chId,2));
00155
00156 for(vector<DTSuperLayerId>::const_iterator itSl = slIds.begin(); itSl != slIds.end(); ++itSl){
00157
00158 float ttrigMean = 0;
00159 float ttrigSigma = 0;
00160 float kFactor = 0;
00161 tTrigMap_->get(*itSl,ttrigMean,ttrigSigma,kFactor,DTTimeUnits::ns);
00162
00163
00164 float ttrigMeanNew = ttrigMean;
00165 float ttrigSigmaNew = ttrigSigma;
00166 float t0SegMean = (itSl->superLayer() != 2)?itChHistos->second[0]->GetMean():itChHistos->second[1]->GetMean();
00167
00168 float kFactorNew = (kFactor*ttrigSigma+t0SegMean)/ttrigSigma;
00169
00170 tTrig->set(*itSl,ttrigMeanNew,ttrigSigmaNew,kFactorNew,DTTimeUnits::ns);
00171 }
00172 }
00173 LogVerbatim("Calibration")<< "[DTTTrigOffsetCalibration] Writing ttrig object to DB!" << endl;
00174
00175 string tTrigRecord = "DTTtrigRcd";
00176 DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
00177 }
00178 }
00179
00180
00181 void DTTTrigOffsetCalibration::bookHistos(DTChamberId chId) {
00182
00183 LogTrace("Calibration") << " Booking histos for Chamber: " << chId;
00184
00185
00186 stringstream wheel; wheel << chId.wheel();
00187 stringstream station; station << chId.station();
00188 stringstream sector; sector << chId.sector();
00189
00190 string chHistoName =
00191 "_W" + wheel.str() +
00192 "_St" + station.str() +
00193 "_Sec" + sector.str();
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213 vector<TH1F*> histos;
00214
00215 histos.push_back(new TH1F(("hRPhiSegT0"+chHistoName).c_str(), "t0 from Phi segments", 500, -60., 60.));
00216 if(chId.station() != 4) histos.push_back(new TH1F(("hRZSegT0"+chHistoName).c_str(), "t0 from Z segments", 500, -60., 60.));
00217
00218 theT0SegHistoMap_[chId] = histos;
00219 }