CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTTTrigResidualCorrection Class Reference

#include <DTTTrigResidualCorrection.h>

Inheritance diagram for DTTTrigResidualCorrection:
DTTTrigBaseCorrection

Public Member Functions

virtual DTTTrigData correction (const DTSuperLayerId &)
 
 DTTTrigResidualCorrection (const edm::ParameterSet &)
 
virtual void setES (const edm::EventSetup &setup)
 
virtual ~DTTTrigResidualCorrection ()
 
- Public Member Functions inherited from DTTTrigBaseCorrection
 DTTTrigBaseCorrection ()
 
virtual ~DTTTrigBaseCorrection ()
 

Private Member Functions

const TH1F * getHisto (const DTSuperLayerId &)
 
std::string getHistoName (const DTSuperLayerId &slID)
 

Private Attributes

std::string dbLabel_
 
DTResidualFitterfitter_
 
const DTMtimemTimeMap_
 
std::string rootBaseDir_
 
TFile * rootFile_
 
const DTTtrigtTrigMap_
 
bool useFit_
 
bool useSlopesCalib_
 
double vDriftEff_ [5][14][4][3]
 

Detailed Description

Concrete implementation of a DTTTrigBaseCorrection. Computes residual correction for tTrig

Revision:
1.5
Author
A. Vilela Pereira

Definition at line 27 of file DTTTrigResidualCorrection.h.

Constructor & Destructor Documentation

DTTTrigResidualCorrection::DTTTrigResidualCorrection ( const edm::ParameterSet pset)

Definition at line 42 of file DTTTrigResidualCorrection.cc.

References debug, edm::FileInPath::fullPath(), edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().

42  {
43  string residualsRootFile = pset.getParameter<string>("residualsRootFile");
44  rootFile_ = new TFile(residualsRootFile.c_str(),"READ");
45  rootBaseDir_ = pset.getUntrackedParameter<string>("rootBaseDir","/DQMData/DT/DTCalibValidation");
46  useFit_ = pset.getParameter<bool>("useFitToResiduals");
47  //useConstantvDrift_ = pset.getParameter<bool>("useConstantDriftVelocity");
48  dbLabel_ = pset.getUntrackedParameter<string>("dbLabel", "");
49  useSlopesCalib_ = pset.getUntrackedParameter<bool>("useSlopesCalib",false);
50 
51  // Load external slopes
52  if(useSlopesCalib_){
53  ifstream fileSlopes;
54  fileSlopes.open( pset.getParameter<FileInPath>("slopesFileName").fullPath().c_str() );
55 
56  int tmp_wheel = 0, tmp_sector = 0, tmp_station = 0, tmp_SL = 0;
57  double tmp_ttrig = 0., tmp_t0 = 0., tmp_kfact = 0.;
58  int tmp_a = 0, tmp_b = 0, tmp_c = 0, tmp_d = 0;
59  double tmp_v_eff = 0.;
60  while(!fileSlopes.eof()){
61  fileSlopes >> tmp_wheel >> tmp_sector >> tmp_station >> tmp_SL >> tmp_a >> tmp_b >>
62  tmp_ttrig >> tmp_t0 >> tmp_kfact >> tmp_c >> tmp_d >> tmp_v_eff;
63  vDriftEff_[tmp_wheel+2][tmp_sector-1][tmp_station-1][tmp_SL-1] = -tmp_v_eff;
64  }
65  fileSlopes.close();
66  }
67 
68  bool debug = pset.getUntrackedParameter<bool>("debug",false);
69  fitter_ = new DTResidualFitter(debug);
70 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::string fullPath() const
Definition: FileInPath.cc:171
#define debug
Definition: MEtoEDMFormat.h:34
DTTTrigResidualCorrection::~DTTTrigResidualCorrection ( )
virtual

Definition at line 72 of file DTTTrigResidualCorrection.cc.

72  {
73  delete rootFile_;
74  delete fitter_;
75 }

Member Function Documentation

DTTTrigData DTTTrigResidualCorrection::correction ( const DTSuperLayerId slId)
virtual

Implements DTTTrigBaseCorrection.

Definition at line 91 of file DTTTrigResidualCorrection.cc.

References DTVelocityUnits::cm_per_ns, edm::hlt::Exception, DTResidualFitResult::fitMean, DTResidualFitResult::fitSigma, cropTnPTrees::frac, getHisto(), LogTrace, timingPdfMaker::mean, DTTimeUnits::ns, DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, ntuplemaker::status, DTSuperLayerId::superLayer(), DTChamberId::wheel(), and vdt::x.

91  {
92 
93  float tTrigMean,tTrigSigma,kFactor;
94  int status = tTrigMap_->get(slId,tTrigMean,tTrigSigma,kFactor,DTTimeUnits::ns);
95  if(status != 0) throw cms::Exception("[DTTTrigResidualCorrection]") << "Could not find tTrig entry in DB for"
96  << slId << endl;
97 
98  float vDrift,hitResolution;
99  status = mTimeMap_->get(slId,vDrift,hitResolution,DTVelocityUnits::cm_per_ns);
100  if(status != 0) throw cms::Exception("[DTTTrigResidualCorrection]") << "Could not find vDrift entry in DB for"
101  << slId << endl;
102  TH1F residualHisto = *(getHisto(slId));
103  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: \n"
104  << " Mean, RMS = " << residualHisto.GetMean() << ", " << residualHisto.GetRMS();
105 
106  double fitMean = -1.;
107  double fitSigma = -1.;
108  if(useFit_){
109  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Fitting histogram " << residualHisto.GetName();
110  const bool originalFit = false; // FIXME: Include this option in fitter class
111  if(originalFit){
112  RooRealVar x("x","residual",-1.,1.);
113  RooRealVar mean("mean","mean",residualHisto.GetMean(),-0.3,0.3);
114  RooRealVar sigma1("sigma1","sigma1",0.,0.5);
115  RooRealVar sigma2("sigma2","sigma2",0.,0.5);
116 
117  RooRealVar frac("frac","frac",0.,1.);
118 
119  RooGaussian myg1("myg1","Gaussian distribution",x,mean,sigma1);
120  RooGaussian myg2("myg2","Gaussian distribution",x,mean,sigma2);
121 
122  RooAddPdf myg("myg","myg",RooArgList(myg1,myg2),RooArgList(frac));
123 
124  RooDataHist hdata("hdata","Binned data",RooArgList(x),&residualHisto);
125  myg.fitTo(hdata,RooFit::Minos(0),RooFit::Range(-0.2,0.2));
126 
127  fitMean = mean.getVal();
128  fitSigma = sigma1.getVal();
129  } else{
130  int nSigmas = 2;
131  DTResidualFitResult fitResult = fitter_->fitResiduals(residualHisto,nSigmas);
132  fitMean = fitResult.fitMean;
133  fitSigma = fitResult.fitSigma;
134  }
135  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: \n"
136  << " Fit Mean = " << fitMean << "\n"
137  << " Fit Sigma = " << fitSigma;
138  }
139  double resMean = (useFit_) ? fitMean : residualHisto.GetMean();
140 
141  int wheel = slId.wheel();
142  int sector = slId.sector();
143  int station = slId.station();
144  int superLayer = slId.superLayer();
145  double resTime = 0.;
146  if(useSlopesCalib_){
147  double vdrift_eff = vDriftEff_[wheel+2][sector-1][station-1][superLayer-1];
148  if(vdrift_eff == 0) vdrift_eff = vDrift;
149 
150  if(vdrift_eff) resTime = resMean/vdrift_eff;
151 
152  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Effective vDrift, correction to tTrig = "
153  << vdrift_eff << ", " << resTime;
154  } else{
155  if(vDrift) resTime = resMean/vDrift;
156 
157  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: vDrift from DB, correction to tTrig = "
158  << vDrift << ", " << resTime;
159  }
160 
161  double corrMean = tTrigMean;
162  double corrSigma = (tTrigSigma != 0.) ? tTrigSigma : 1.;
163  double corrKFact = (tTrigSigma != 0.) ? (kFactor + resTime/tTrigSigma) : resTime;
164 
165  return DTTTrigData(corrMean,corrSigma,corrKFact);
166 }
int superLayer() const
Return the superlayer number.
#define LogTrace(id)
PixelRecoRange< float > Range
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
DTResidualFitResult fitResiduals(TH1F &histo, int nSigmas=1)
int get(int wheelId, int stationId, int sectorId, int slId, float &mTime, float &mTrms, DTTimeUnits::type unit) const
Definition: DTMtime.cc:86
const TH1F * getHisto(const DTSuperLayerId &)
int sector() const
Definition: DTChamberId.h:63
tuple status
Definition: ntuplemaker.py:245
int station() const
Return the station number.
Definition: DTChamberId.h:53
x
Definition: VDTMath.h:216
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
const TH1F * DTTTrigResidualCorrection::getHisto ( const DTSuperLayerId slId)
private

Definition at line 168 of file DTTTrigResidualCorrection.cc.

References edm::hlt::Exception, mergeVDriftHistosByStation::getHistoName(), interpolateCardsSimple::histo, and LogTrace.

168  {
169  string histoName = getHistoName(slId);
170  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Accessing histogram " << histoName.c_str();
171  TH1F* histo = static_cast<TH1F*>(rootFile_->Get(histoName.c_str()));
172  if(!histo) throw cms::Exception("[DTTTrigResidualCorrection]") << "residual histogram not found:"
173  << histoName << endl;
174  return histo;
175 }
std::string getHistoName(const DTSuperLayerId &slID)
#define LogTrace(id)
string DTTTrigResidualCorrection::getHistoName ( const DTSuperLayerId slID)
private

Definition at line 177 of file DTTTrigResidualCorrection.cc.

References DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, launcher::step, cmsPerfCommons::Step, DTSuperLayerId::superlayer(), and DTChamberId::wheel().

177  {
178 
179  int step = 3;
180  stringstream wheel; wheel << slId.wheel();
181  stringstream station; station << slId.station();
182  stringstream sector; sector << slId.sector();
183  stringstream superLayer; superLayer << slId.superlayer();
184  stringstream Step; Step << step;
185 
186  string histoName =
187  rootBaseDir_ + "/Wheel" + wheel.str() +
188  "/Station" + station.str() +
189  "/Sector" + sector.str() +
190  "/hResDist_STEP" + Step.str() +
191  "_W" + wheel.str() +
192  "_St" + station.str() +
193  "_Sec" + sector.str() +
194  "_SL" + superLayer.str();
195 
196  return histoName;
197 }
list step
Definition: launcher.py:15
void DTTTrigResidualCorrection::setES ( const edm::EventSetup setup)
virtual

Implements DTTTrigBaseCorrection.

Definition at line 77 of file DTTTrigResidualCorrection.cc.

References edm::EventSetup::get().

77  {
78  // Get tTrig record from DB
79  ESHandle<DTTtrig> tTrig;
80  //setup.get<DTTtrigRcd>().get(tTrig);
81  setup.get<DTTtrigRcd>().get(dbLabel_,tTrig);
82  tTrigMap_ = &*tTrig;
83 
84  // Get vDrift record
85  ESHandle<DTMtime> mTimeHandle;
86  setup.get<DTMtimeRcd>().get(mTimeHandle);
87  mTimeMap_ = &*mTimeHandle;
88 
89 }
const T & get() const
Definition: EventSetup.h:55

Member Data Documentation

std::string DTTTrigResidualCorrection::dbLabel_
private

Definition at line 46 of file DTTTrigResidualCorrection.h.

DTResidualFitter* DTTTrigResidualCorrection::fitter_
private

Definition at line 53 of file DTTTrigResidualCorrection.h.

const DTMtime* DTTTrigResidualCorrection::mTimeMap_
private

Definition at line 52 of file DTTTrigResidualCorrection.h.

std::string DTTTrigResidualCorrection::rootBaseDir_
private

Definition at line 44 of file DTTTrigResidualCorrection.h.

TFile* DTTTrigResidualCorrection::rootFile_
private

Definition at line 42 of file DTTTrigResidualCorrection.h.

const DTTtrig* DTTTrigResidualCorrection::tTrigMap_
private

Definition at line 51 of file DTTTrigResidualCorrection.h.

bool DTTTrigResidualCorrection::useFit_
private

Definition at line 45 of file DTTTrigResidualCorrection.h.

bool DTTTrigResidualCorrection::useSlopesCalib_
private

Definition at line 47 of file DTTTrigResidualCorrection.h.

double DTTTrigResidualCorrection::vDriftEff_[5][14][4][3]
private

Definition at line 49 of file DTTTrigResidualCorrection.h.