CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTTTrigResidualCorrection.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author A. Vilela Pereira
5  */
6 
12 
18 
20 
21 #include "TFile.h"
22 #include "TH1F.h"
23 #include "TF1.h"
24 #include "TCanvas.h"
25 
26 #include "RooPlot.h"
27 #include "RooRealVar.h"
28 #include "RooDataHist.h"
29 #include "RooGaussian.h"
30 #include "RooAddPdf.h"
31 #include "RooFitResult.h"
32 #include "RooGlobalFunc.h"
33 
34 #include <string>
35 #include <sstream>
36 #include <fstream>
37 
38 using namespace std;
39 using namespace edm;
40 
41 namespace dtCalibration {
42 
43 DTTTrigResidualCorrection::DTTTrigResidualCorrection(const ParameterSet& pset) {
44  string residualsRootFile = pset.getParameter<string>("residualsRootFile");
45  rootFile_ = new TFile(residualsRootFile.c_str(),"READ");
46  rootBaseDir_ = pset.getUntrackedParameter<string>("rootBaseDir","/DQMData/DT/DTCalibValidation");
47  useFit_ = pset.getParameter<bool>("useFitToResiduals");
48  //useConstantvDrift_ = pset.getParameter<bool>("useConstantDriftVelocity");
49  dbLabel_ = pset.getUntrackedParameter<string>("dbLabel", "");
50  useSlopesCalib_ = pset.getUntrackedParameter<bool>("useSlopesCalib",false);
51 
52  // Load external slopes
53  if(useSlopesCalib_){
54  ifstream fileSlopes;
55  fileSlopes.open( pset.getParameter<FileInPath>("slopesFileName").fullPath().c_str() );
56 
57  int tmp_wheel = 0, tmp_sector = 0, tmp_station = 0, tmp_SL = 0;
58  double tmp_ttrig = 0., tmp_t0 = 0., tmp_kfact = 0.;
59  int tmp_a = 0, tmp_b = 0, tmp_c = 0, tmp_d = 0;
60  double tmp_v_eff = 0.;
61  while(!fileSlopes.eof()){
62  fileSlopes >> tmp_wheel >> tmp_sector >> tmp_station >> tmp_SL >> tmp_a >> tmp_b >>
63  tmp_ttrig >> tmp_t0 >> tmp_kfact >> tmp_c >> tmp_d >> tmp_v_eff;
64  vDriftEff_[tmp_wheel+2][tmp_sector-1][tmp_station-1][tmp_SL-1] = -tmp_v_eff;
65  }
66  fileSlopes.close();
67  }
68 
69  bool debug = pset.getUntrackedParameter<bool>("debug",false);
70  fitter_ = new DTResidualFitter(debug);
71 }
72 
73 DTTTrigResidualCorrection::~DTTTrigResidualCorrection() {
74  delete rootFile_;
75  delete fitter_;
76 }
77 
78 void DTTTrigResidualCorrection::setES(const EventSetup& setup) {
79  // Get tTrig record from DB
80  ESHandle<DTTtrig> tTrig;
81  //setup.get<DTTtrigRcd>().get(tTrig);
82  setup.get<DTTtrigRcd>().get(dbLabel_,tTrig);
83  tTrigMap_ = &*tTrig;
84 
85  // Get vDrift record
86  ESHandle<DTMtime> mTimeHandle;
87  setup.get<DTMtimeRcd>().get(mTimeHandle);
88  mTimeMap_ = &*mTimeHandle;
89 
90 }
91 
92 DTTTrigData DTTTrigResidualCorrection::correction(const DTSuperLayerId& slId) {
93 
94  float tTrigMean,tTrigSigma,kFactor;
95  int status = tTrigMap_->get(slId,tTrigMean,tTrigSigma,kFactor,DTTimeUnits::ns);
96  if(status != 0) throw cms::Exception("[DTTTrigResidualCorrection]") << "Could not find tTrig entry in DB for"
97  << slId << endl;
98 
99  float vDrift,hitResolution;
100  status = mTimeMap_->get(slId,vDrift,hitResolution,DTVelocityUnits::cm_per_ns);
101  if(status != 0) throw cms::Exception("[DTTTrigResidualCorrection]") << "Could not find vDrift entry in DB for"
102  << slId << endl;
103  TH1F residualHisto = *(getHisto(slId));
104  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: \n"
105  << " Mean, RMS = " << residualHisto.GetMean() << ", " << residualHisto.GetRMS();
106 
107  double fitMean = -1.;
108  double fitSigma = -1.;
109  if(useFit_){
110  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Fitting histogram " << residualHisto.GetName();
111  const bool originalFit = false; // FIXME: Include this option in fitter class
112  if(originalFit){
113  RooRealVar x("x","residual",-1.,1.);
114  RooRealVar mean("mean","mean",residualHisto.GetMean(),-0.3,0.3);
115  RooRealVar sigma1("sigma1","sigma1",0.,0.5);
116  RooRealVar sigma2("sigma2","sigma2",0.,0.5);
117 
118  RooRealVar frac("frac","frac",0.,1.);
119 
120  RooGaussian myg1("myg1","Gaussian distribution",x,mean,sigma1);
121  RooGaussian myg2("myg2","Gaussian distribution",x,mean,sigma2);
122 
123  RooAddPdf myg("myg","myg",RooArgList(myg1,myg2),RooArgList(frac));
124 
125  RooDataHist hdata("hdata","Binned data",RooArgList(x),&residualHisto);
126  myg.fitTo(hdata,RooFit::Minos(0),RooFit::Range(-0.2,0.2));
127 
128  fitMean = mean.getVal();
129  fitSigma = sigma1.getVal();
130  } else{
131  int nSigmas = 2;
132  DTResidualFitResult fitResult = fitter_->fitResiduals(residualHisto,nSigmas);
133  fitMean = fitResult.fitMean;
134  fitSigma = fitResult.fitSigma;
135  }
136  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: \n"
137  << " Fit Mean = " << fitMean << "\n"
138  << " Fit Sigma = " << fitSigma;
139  }
140  double resMean = (useFit_) ? fitMean : residualHisto.GetMean();
141 
142  int wheel = slId.wheel();
143  int sector = slId.sector();
144  int station = slId.station();
145  int superLayer = slId.superLayer();
146  double resTime = 0.;
147  if(useSlopesCalib_){
148  double vdrift_eff = vDriftEff_[wheel+2][sector-1][station-1][superLayer-1];
149  if(vdrift_eff == 0) vdrift_eff = vDrift;
150 
151  if(vdrift_eff) resTime = resMean/vdrift_eff;
152 
153  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Effective vDrift, correction to tTrig = "
154  << vdrift_eff << ", " << resTime;
155  } else{
156  if(vDrift) resTime = resMean/vDrift;
157 
158  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: vDrift from DB, correction to tTrig = "
159  << vDrift << ", " << resTime;
160  }
161 
162  double corrMean = tTrigMean;
163  double corrSigma = (tTrigSigma != 0.) ? tTrigSigma : 1.;
164  double corrKFact = (tTrigSigma != 0.) ? (kFactor + resTime/tTrigSigma) : resTime;
165 
166  return DTTTrigData(corrMean,corrSigma,corrKFact);
167 }
168 
170  string histoName = getHistoName(slId);
171  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Accessing histogram " << histoName.c_str();
172  TH1F* histo = static_cast<TH1F*>(rootFile_->Get(histoName.c_str()));
173  if(!histo) throw cms::Exception("[DTTTrigResidualCorrection]") << "residual histogram not found:"
174  << histoName << endl;
175  return histo;
176 }
177 
179 
180  int step = 3;
181  stringstream wheel; wheel << slId.wheel();
182  stringstream station; station << slId.station();
183  stringstream sector; sector << slId.sector();
184  stringstream superLayer; superLayer << slId.superlayer();
185  stringstream Step; Step << step;
186 
187  string histoName =
188  rootBaseDir_ + "/Wheel" + wheel.str() +
189  "/Station" + station.str() +
190  "/Sector" + sector.str() +
191  "/hResDist_STEP" + Step.str() +
192  "_W" + wheel.str() +
193  "_St" + station.str() +
194  "_Sec" + sector.str() +
195  "_SL" + superLayer.str();
196 
197  return histoName;
198 }
199 
200 } // namespace
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
T x() const
Cartesian x coordinate.
int superLayer() const
Return the superlayer number.
#define LogTrace(id)
int superlayer() const
Return the superlayer number (deprecated method name)
PixelRecoRange< float > Range
#define debug
Definition: HDRShower.cc:19
const T & get() const
Definition: EventSetup.h:56
int sector() const
Definition: DTChamberId.h:61
std::string fullPath() const
Definition: FileInPath.cc:165
tuple status
Definition: ntuplemaker.py:245
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
TH1F * getHisto(std::string name, std::string process, DQMStore *dbe_, bool verb=false, bool clone=false)
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")