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_
 
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.3
Author
A. Vilela Pereira

Definition at line 26 of file DTTTrigResidualCorrection.h.

Constructor & Destructor Documentation

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

Definition at line 40 of file DTTTrigResidualCorrection.cc.

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

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

Definition at line 82 of file DTTTrigResidualCorrection.cc.

82  {
83  delete rootFile_;
84 }

Member Function Documentation

DTTTrigData DTTTrigResidualCorrection::correction ( const DTSuperLayerId slId)
virtual

Implements DTTTrigBaseCorrection.

Definition at line 100 of file DTTTrigResidualCorrection.cc.

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

100  {
101 
102  float tTrigMean,tTrigSigma,kFactor;
103  int status = tTrigMap_->get(slId,tTrigMean,tTrigSigma,kFactor,DTTimeUnits::ns);
104  if(status != 0) throw cms::Exception("[DTTTrigResidualCorrection]") << "Could not find tTrig entry in DB for"
105  << slId << endl;
106 
107  float vDrift,hitResolution;
108  status = mTimeMap_->get(slId,vDrift,hitResolution,DTVelocityUnits::cm_per_ns);
109  if(status != 0) throw cms::Exception("[DTTTrigResidualCorrection]") << "Could not find vDrift entry in DB for"
110  << slId << endl;
111  TH1F residualHisto = *(getHisto(slId));
112  double fitMean = -1.;
113 
114  if(useFit_){
115 
116  RooRealVar x("x","residual",-1.,1.);
117  RooRealVar mean("mean","mean",residualHisto.GetMean(),-0.3,0.3);
118  RooRealVar sigma1("sigma1","sigma1",0.,0.5);
119  RooRealVar sigma2("sigma2","sigma2",0.,0.5);
120 
121  RooRealVar frac("frac","frac",0.,1.);
122 
123  RooGaussian myg1("myg1","Gaussian distribution",x,mean,sigma1);
124  RooGaussian myg2("myg2","Gaussian distribution",x,mean,sigma2);
125 
126  RooAddPdf myg("myg","myg",RooArgList(myg1,myg2),RooArgList(frac));
127 
128  RooDataHist hdata("hdata","Binned data",RooArgList(x),&residualHisto);
129  myg.fitTo(hdata,RooFit::Minos(0),RooFit::Range(-0.2,0.2));
130 
131  fitMean = mean.getVal();
132  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: \n"
133  << " Mean, Fit Mean = " << residualHisto.GetMean()
134  << ", " << fitMean << "\n"
135  << " RMS, Fit RMS = " << residualHisto.GetRMS()
136  << ", " << sigma1.getVal();
137 
138  //static int count = 0;
139  /*
140  if(count == 0){
141  RooPlot *xframe = x.frame();
142  hdata.plotOn(xframe);
143  myg.plotOn(xframe);
144  TCanvas c1;
145  c1.cd();
146  xframe->Draw();
147  c1.SaveAs("prova.eps");
148  count++;
149  }
150  */
151 
152  }
153 
154  int wheel = slId.wheel();
155  int sector = slId.sector();
156  int station = slId.station();
157  int superLayer = slId.superLayer();
158 
159  double resTime = 0.;
160  if(useSlopesCalib_){
161  double vdrift_eff = vDriftEff_[wheel+2][sector-1][station-1][superLayer-1];
162  if(vdrift_eff == 0) vdrift_eff = vDrift;
163 
164  if(vdrift_eff) resTime = ( (useFit_) ? fitMean : residualHisto.GetMean() )/vdrift_eff;
165 
166  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Effective vDrift, correction to tTrig = "
167  << vdrift_eff << ", " << resTime;
168  } else{
169  if(vDrift) resTime = ( (useFit_) ? fitMean : residualHisto.GetMean() )/vDrift;
170 
171  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: vDrift from DB, correction to tTrig = "
172  << vDrift << ", " << resTime;
173  }
174 
175  double corrMean = tTrigMean;
176  double corrSigma = (tTrigSigma != 0.) ? tTrigSigma : 1.;
177  double corrKFact = (tTrigSigma != 0.) ? (kFactor + resTime/tTrigSigma) : resTime;
178 
179  return DTTTrigData(corrMean,corrSigma,corrKFact);
180 }
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
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
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
const TH1F * DTTTrigResidualCorrection::getHisto ( const DTSuperLayerId slId)
private

Definition at line 182 of file DTTTrigResidualCorrection.cc.

References edm::hlt::Exception, trackerHits::histo, and LogTrace.

182  {
183  string histoName = getHistoName(slId);
184  LogTrace("Calibration") << "[DTTTrigResidualCorrection]: Accessing histogram " << histoName.c_str();
185  TH1F* histo = static_cast<TH1F*>(rootFile_->Get(histoName.c_str()));
186  if(!histo) throw cms::Exception("[DTTTrigResidualCorrection]") << "residual histogram not found:"
187  << histoName << endl;
188  return histo;
189 }
tuple histo
Definition: trackerHits.py:12
std::string getHistoName(const DTSuperLayerId &slID)
#define LogTrace(id)
string DTTTrigResidualCorrection::getHistoName ( const DTSuperLayerId slID)
private

Definition at line 191 of file DTTTrigResidualCorrection.cc.

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

191  {
192 
193  int step = 3;
194  stringstream wheel; wheel << slId.wheel();
195  stringstream station; station << slId.station();
196  stringstream sector; sector << slId.sector();
197  stringstream superLayer; superLayer << slId.superlayer();
198  stringstream Step; Step << step;
199 
200  string histoName =
201  rootBaseDir_ + "/Wheel" + wheel.str() +
202  "/Station" + station.str() +
203  "/Sector" + sector.str() +
204  "/hResDist_STEP" + Step.str() +
205  "_W" + wheel.str() +
206  "_St" + station.str() +
207  "_Sec" + sector.str() +
208  "_SL" + superLayer.str();
209 
210  return histoName;
211 }
void DTTTrigResidualCorrection::setES ( const edm::EventSetup setup)
virtual

Implements DTTTrigBaseCorrection.

Definition at line 86 of file DTTTrigResidualCorrection.cc.

References edm::EventSetup::get().

86  {
87  // Get tTrig record from DB
88  ESHandle<DTTtrig> tTrig;
89  //setup.get<DTTtrigRcd>().get(tTrig);
90  setup.get<DTTtrigRcd>().get(dbLabel_,tTrig);
91  tTrigMap_ = &*tTrig;
92 
93  // Get vDrift record
94  ESHandle<DTMtime> mTimeHandle;
95  setup.get<DTMtimeRcd>().get(mTimeHandle);
96  mTimeMap_ = &*mTimeHandle;
97 
98 }
const T & get() const
Definition: EventSetup.h:55

Member Data Documentation

std::string DTTTrigResidualCorrection::dbLabel_
private

Definition at line 45 of file DTTTrigResidualCorrection.h.

const DTMtime* DTTTrigResidualCorrection::mTimeMap_
private

Definition at line 51 of file DTTTrigResidualCorrection.h.

std::string DTTTrigResidualCorrection::rootBaseDir_
private

Definition at line 43 of file DTTTrigResidualCorrection.h.

TFile* DTTTrigResidualCorrection::rootFile_
private

Definition at line 41 of file DTTTrigResidualCorrection.h.

const DTTtrig* DTTTrigResidualCorrection::tTrigMap_
private

Definition at line 50 of file DTTTrigResidualCorrection.h.

bool DTTTrigResidualCorrection::useFit_
private

Definition at line 44 of file DTTTrigResidualCorrection.h.

bool DTTTrigResidualCorrection::useSlopesCalib_
private

Definition at line 46 of file DTTTrigResidualCorrection.h.

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

Definition at line 48 of file DTTTrigResidualCorrection.h.