00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "DTTTrigResidualCorrection.h"
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 #include "FWCore/Framework/interface/EventSetup.h"
00012 #include "FWCore/Framework/interface/ESHandle.h"
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014
00015 #include "DataFormats/MuonDetId/interface/DTSuperLayerId.h"
00016 #include "CondFormats/DTObjects/interface/DTTtrig.h"
00017 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
00018 #include "CondFormats/DTObjects/interface/DTMtime.h"
00019 #include "CondFormats/DataRecord/interface/DTMtimeRcd.h"
00020
00021 #include "CalibMuon/DTCalibration/interface/DTResidualFitter.h"
00022
00023 #include "TFile.h"
00024 #include "TH1F.h"
00025 #include "TF1.h"
00026 #include "TCanvas.h"
00027
00028 #include "RooPlot.h"
00029 #include "RooRealVar.h"
00030 #include "RooDataHist.h"
00031 #include "RooGaussian.h"
00032 #include "RooAddPdf.h"
00033 #include "RooFitResult.h"
00034 #include "RooGlobalFunc.h"
00035
00036 #include <string>
00037 #include <sstream>
00038
00039 using namespace std;
00040 using namespace edm;
00041
00042 DTTTrigResidualCorrection::DTTTrigResidualCorrection(const ParameterSet& pset) {
00043 string residualsRootFile = pset.getParameter<string>("residualsRootFile");
00044 rootFile_ = new TFile(residualsRootFile.c_str(),"READ");
00045 rootBaseDir_ = pset.getUntrackedParameter<string>("rootBaseDir","/DQMData/DT/DTCalibValidation");
00046 useFit_ = pset.getParameter<bool>("useFitToResiduals");
00047
00048 dbLabel_ = pset.getUntrackedParameter<string>("dbLabel", "");
00049 useSlopesCalib_ = pset.getUntrackedParameter<bool>("useSlopesCalib",false);
00050
00051
00052 if(useSlopesCalib_){
00053 ifstream fileSlopes;
00054 fileSlopes.open( pset.getParameter<FileInPath>("slopesFileName").fullPath().c_str() );
00055
00056 int tmp_wheel = 0, tmp_sector = 0, tmp_station = 0, tmp_SL = 0;
00057 double tmp_ttrig = 0., tmp_t0 = 0., tmp_kfact = 0.;
00058 int tmp_a = 0, tmp_b = 0, tmp_c = 0, tmp_d = 0;
00059 double tmp_v_eff = 0.;
00060 while(!fileSlopes.eof()){
00061 fileSlopes >> tmp_wheel >> tmp_sector >> tmp_station >> tmp_SL >> tmp_a >> tmp_b >>
00062 tmp_ttrig >> tmp_t0 >> tmp_kfact >> tmp_c >> tmp_d >> tmp_v_eff;
00063 vDriftEff_[tmp_wheel+2][tmp_sector-1][tmp_station-1][tmp_SL-1] = -tmp_v_eff;
00064 }
00065 fileSlopes.close();
00066 }
00067
00068 bool debug = pset.getUntrackedParameter<bool>("debug",false);
00069 fitter_ = new DTResidualFitter(debug);
00070 }
00071
00072 DTTTrigResidualCorrection::~DTTTrigResidualCorrection() {
00073 delete rootFile_;
00074 delete fitter_;
00075 }
00076
00077 void DTTTrigResidualCorrection::setES(const EventSetup& setup) {
00078
00079 ESHandle<DTTtrig> tTrig;
00080
00081 setup.get<DTTtrigRcd>().get(dbLabel_,tTrig);
00082 tTrigMap_ = &*tTrig;
00083
00084
00085 ESHandle<DTMtime> mTimeHandle;
00086 setup.get<DTMtimeRcd>().get(mTimeHandle);
00087 mTimeMap_ = &*mTimeHandle;
00088
00089 }
00090
00091 DTTTrigData DTTTrigResidualCorrection::correction(const DTSuperLayerId& slId) {
00092
00093 float tTrigMean,tTrigSigma,kFactor;
00094 int status = tTrigMap_->get(slId,tTrigMean,tTrigSigma,kFactor,DTTimeUnits::ns);
00095 if(status != 0) throw cms::Exception("[DTTTrigResidualCorrection]") << "Could not find tTrig entry in DB for"
00096 << slId << endl;
00097
00098 float vDrift,hitResolution;
00099 status = mTimeMap_->get(slId,vDrift,hitResolution,DTVelocityUnits::cm_per_ns);
00100 if(status != 0) throw cms::Exception("[DTTTrigResidualCorrection]") << "Could not find vDrift entry in DB for"
00101 << slId << endl;
00102 TH1F residualHisto = *(getHisto(slId));
00103 LogTrace("Calibration") << "[DTTTrigResidualCorrection]: \n"
00104 << " Mean, RMS = " << residualHisto.GetMean() << ", " << residualHisto.GetRMS();
00105
00106 double fitMean = -1.;
00107 double fitSigma = -1.;
00108 if(useFit_){
00109 LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Fitting histogram " << residualHisto.GetName();
00110 const bool originalFit = false;
00111 if(originalFit){
00112 RooRealVar x("x","residual",-1.,1.);
00113 RooRealVar mean("mean","mean",residualHisto.GetMean(),-0.3,0.3);
00114 RooRealVar sigma1("sigma1","sigma1",0.,0.5);
00115 RooRealVar sigma2("sigma2","sigma2",0.,0.5);
00116
00117 RooRealVar frac("frac","frac",0.,1.);
00118
00119 RooGaussian myg1("myg1","Gaussian distribution",x,mean,sigma1);
00120 RooGaussian myg2("myg2","Gaussian distribution",x,mean,sigma2);
00121
00122 RooAddPdf myg("myg","myg",RooArgList(myg1,myg2),RooArgList(frac));
00123
00124 RooDataHist hdata("hdata","Binned data",RooArgList(x),&residualHisto);
00125 myg.fitTo(hdata,RooFit::Minos(0),RooFit::Range(-0.2,0.2));
00126
00127 fitMean = mean.getVal();
00128 fitSigma = sigma1.getVal();
00129 } else{
00130 int nSigmas = 2;
00131 DTResidualFitResult fitResult = fitter_->fitResiduals(residualHisto,nSigmas);
00132 fitMean = fitResult.fitMean;
00133 fitSigma = fitResult.fitSigma;
00134 }
00135 LogTrace("Calibration") << "[DTTTrigResidualCorrection]: \n"
00136 << " Fit Mean = " << fitMean << "\n"
00137 << " Fit Sigma = " << fitSigma;
00138 }
00139 double resMean = (useFit_) ? fitMean : residualHisto.GetMean();
00140
00141 int wheel = slId.wheel();
00142 int sector = slId.sector();
00143 int station = slId.station();
00144 int superLayer = slId.superLayer();
00145 double resTime = 0.;
00146 if(useSlopesCalib_){
00147 double vdrift_eff = vDriftEff_[wheel+2][sector-1][station-1][superLayer-1];
00148 if(vdrift_eff == 0) vdrift_eff = vDrift;
00149
00150 if(vdrift_eff) resTime = resMean/vdrift_eff;
00151
00152 LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Effective vDrift, correction to tTrig = "
00153 << vdrift_eff << ", " << resTime;
00154 } else{
00155 if(vDrift) resTime = resMean/vDrift;
00156
00157 LogTrace("Calibration") << "[DTTTrigResidualCorrection]: vDrift from DB, correction to tTrig = "
00158 << vDrift << ", " << resTime;
00159 }
00160
00161 double corrMean = tTrigMean;
00162 double corrSigma = (tTrigSigma != 0.) ? tTrigSigma : 1.;
00163 double corrKFact = (tTrigSigma != 0.) ? (kFactor + resTime/tTrigSigma) : resTime;
00164
00165 return DTTTrigData(corrMean,corrSigma,corrKFact);
00166 }
00167
00168 const TH1F* DTTTrigResidualCorrection::getHisto(const DTSuperLayerId& slId) {
00169 string histoName = getHistoName(slId);
00170 LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Accessing histogram " << histoName.c_str();
00171 TH1F* histo = static_cast<TH1F*>(rootFile_->Get(histoName.c_str()));
00172 if(!histo) throw cms::Exception("[DTTTrigResidualCorrection]") << "residual histogram not found:"
00173 << histoName << endl;
00174 return histo;
00175 }
00176
00177 string DTTTrigResidualCorrection::getHistoName(const DTSuperLayerId& slId) {
00178
00179 int step = 3;
00180 stringstream wheel; wheel << slId.wheel();
00181 stringstream station; station << slId.station();
00182 stringstream sector; sector << slId.sector();
00183 stringstream superLayer; superLayer << slId.superlayer();
00184 stringstream Step; Step << step;
00185
00186 string histoName =
00187 rootBaseDir_ + "/Wheel" + wheel.str() +
00188 "/Station" + station.str() +
00189 "/Sector" + sector.str() +
00190 "/hResDist_STEP" + Step.str() +
00191 "_W" + wheel.str() +
00192 "_St" + station.str() +
00193 "_Sec" + sector.str() +
00194 "_SL" + superLayer.str();
00195
00196 return histoName;
00197 }