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
dtCalibration::DTVDriftSegment Class Reference

#include <DTVDriftSegment.h>

Inheritance diagram for dtCalibration::DTVDriftSegment:
dtCalibration::DTVDriftBaseAlgo

Public Member Functions

virtual DTVDriftData compute (const DTSuperLayerId &)
 
 DTVDriftSegment (edm::ParameterSet const &)
 
virtual void setES (const edm::EventSetup &setup)
 
virtual ~DTVDriftSegment ()
 
- Public Member Functions inherited from dtCalibration::DTVDriftBaseAlgo
 DTVDriftBaseAlgo ()
 
virtual ~DTVDriftBaseAlgo ()
 

Private Member Functions

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

Private Attributes

DTResidualFitterfitter_
 
const DTMtimemTimeMap_
 
unsigned int nSigmas_
 
TFile * rootFile_
 

Detailed Description

Definition at line 24 of file DTVDriftSegment.h.

Constructor & Destructor Documentation

DTVDriftSegment::DTVDriftSegment ( edm::ParameterSet const &  pset)

Definition at line 36 of file DTVDriftSegment.cc.

References debug, fitter_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), rootFile_, and dtT0Analyzer_cfg::rootFileName.

36  :
37  nSigmas_( pset.getUntrackedParameter<unsigned int>("nSigmasFitRange", 1) ) {
38 
39  string rootFileName = pset.getParameter<string>("rootFileName");
40  rootFile_ = new TFile(rootFileName.c_str(), "READ");
41 
42  bool debug = pset.getUntrackedParameter<bool>("debug",false);
43  fitter_ = new DTResidualFitter(debug);
44  //bool debug = pset.getUntrackedParameter<bool>("debug", false);
45  //if(debug) fitter_->setVerbosity(1);
46 }
#define debug
Definition: MEtoEDMFormat.h:34
DTVDriftSegment::~DTVDriftSegment ( )
virtual

Definition at line 48 of file DTVDriftSegment.cc.

References fitter_, and rootFile_.

48  {
49  rootFile_->Close();
50  delete fitter_;
51 }

Member Function Documentation

DTVDriftData DTVDriftSegment::compute ( const DTSuperLayerId slId)
virtual

Implements dtCalibration::DTVDriftBaseAlgo.

Definition at line 60 of file DTVDriftSegment.cc.

References DTVelocityUnits::cm_per_ns, edm::hlt::Exception, DTResidualFitResult::fitMean, DTResidualFitResult::fitMeanError, DTResidualFitter::fitResiduals(), DTResidualFitResult::fitSigma, DTResidualFitResult::fitSigmaError, fitter_, DTMtime::get(), getHisto(), LogTrace, mTimeMap_, nSigmas_, dtDQMClient_cfg::resolution, ntuplemaker::status, and DTSuperLayerId::superLayer().

60  {
61 
62  // Get original value from DB; vdrift is cm/ns , resolution is cm
63  float vDrift = 0., resolution = 0.;
65 
66  if(status != 0) throw cms::Exception("DTCalibration") << "Could not find vDrift entry in DB for"
67  << slId << endl;
68 
69  // For RZ superlayers use original value
70  if(slId.superLayer() == 2){
71  LogTrace("Calibration") << "[DTVDriftSegment]: RZ superlayer\n"
72  << " Will use original vDrift and resolution.";
73  return DTVDriftData(vDrift,resolution);
74  } else{
75  TH1F* vDriftCorrHisto = getHisto(slId);
76  // If empty histogram
77  if(vDriftCorrHisto->GetEntries() == 0){
78  LogError("Calibration") << "[DTVDriftSegment]: Histogram " << vDriftCorrHisto->GetName() << " is empty.\n"
79  << " Will use original vDrift and resolution.";
80  return DTVDriftData(vDrift,resolution);
81  }
82 
83  LogTrace("Calibration") << "[DTVDriftSegment]: Fitting histogram " << vDriftCorrHisto->GetName();
84  DTResidualFitResult fitResult = fitter_->fitResiduals(*vDriftCorrHisto,nSigmas_);
85  LogTrace("Calibration") << "[DTVDriftSegment]: \n"
86  << " Fit Mean = " << fitResult.fitMean << " +/- "
87  << fitResult.fitMeanError << "\n"
88  << " Fit Sigma = " << fitResult.fitSigma << " +/- "
89  << fitResult.fitSigmaError;
90 
91  float vDriftCorr = fitResult.fitMean;
92  float vDriftNew = vDrift*(1. - vDriftCorr);
93  float resolutionNew = resolution;
94  return DTVDriftData(vDriftNew,resolutionNew);
95  }
96 }
int superLayer() const
Return the superlayer number.
#define LogTrace(id)
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
tuple status
Definition: ntuplemaker.py:245
TH1F * getHisto(const DTSuperLayerId &)
TH1F * DTVDriftSegment::getHisto ( const DTSuperLayerId slId)
private

Definition at line 98 of file DTVDriftSegment.cc.

References edm::hlt::Exception, getHistoName(), timingPdfMaker::histo, and rootFile_.

Referenced by compute().

98  {
99  string histoName = getHistoName(slId);
100  TH1F* histo = static_cast<TH1F*>(rootFile_->Get(histoName.c_str()));
101  if(!histo) throw cms::Exception("DTCalibration") << "v-drift correction histogram not found:"
102  << histoName << endl;
103  return histo;
104 }
std::string getHistoName(const DTSuperLayerId &)
string DTVDriftSegment::getHistoName ( const DTSuperLayerId slId)
private

Definition at line 106 of file DTVDriftSegment.cc.

References DTSuperLayerId::chamberId(), DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, DTSuperLayerId::superLayer(), and DTChamberId::wheel().

Referenced by getHisto().

106  {
107  DTChamberId chId = slId.chamberId();
108 
109  // Compose the chamber name
110  stringstream wheel; wheel << chId.wheel();
111  stringstream station; station << chId.station();
112  stringstream sector; sector << chId.sector();
113 
114  string chHistoName =
115  "_W" + wheel.str() +
116  "_St" + station.str() +
117  "_Sec" + sector.str();
118 
119  return (slId.superLayer() != 2)?("hRPhiVDriftCorr" + chHistoName):("hRZVDriftCorr" + chHistoName);
120 }
DTChamberId chamberId() const
Return the corresponding ChamberId.
int superLayer() const
Return the superlayer number.
int sector() const
Definition: DTChamberId.h:63
int station() const
Return the station number.
Definition: DTChamberId.h:53
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
void DTVDriftSegment::setES ( const edm::EventSetup setup)
virtual

Implements dtCalibration::DTVDriftBaseAlgo.

Definition at line 53 of file DTVDriftSegment.cc.

References edm::EventSetup::get(), and mTimeMap_.

53  {
54  // Get the map of vdrift from the setup
55  ESHandle<DTMtime> mTime;
56  setup.get<DTMtimeRcd>().get(mTime);
57  mTimeMap_ = &*mTime;
58 }
const T & get() const
Definition: EventSetup.h:55

Member Data Documentation

DTResidualFitter* dtCalibration::DTVDriftSegment::fitter_
private

Definition at line 39 of file DTVDriftSegment.h.

Referenced by compute(), DTVDriftSegment(), and ~DTVDriftSegment().

const DTMtime* dtCalibration::DTVDriftSegment::mTimeMap_
private

Definition at line 37 of file DTVDriftSegment.h.

Referenced by compute(), and setES().

unsigned int dtCalibration::DTVDriftSegment::nSigmas_
private

Definition at line 35 of file DTVDriftSegment.h.

Referenced by compute().

TFile* dtCalibration::DTVDriftSegment::rootFile_
private

Definition at line 38 of file DTVDriftSegment.h.

Referenced by DTVDriftSegment(), getHisto(), and ~DTVDriftSegment().