CMS 3D CMS Logo

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 23 of file DTVDriftSegment.h.

Constructor & Destructor Documentation

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

Definition at line 34 of file DTVDriftSegment.cc.

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

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

Definition at line 46 of file DTVDriftSegment.cc.

References fitter_, and rootFile_.

46  {
47  rootFile_->Close();
48  delete fitter_;
49 }

Member Function Documentation

DTVDriftData DTVDriftSegment::compute ( const DTSuperLayerId )
virtual

Implements dtCalibration::DTVDriftBaseAlgo.

Definition at line 58 of file DTVDriftSegment.cc.

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

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

Definition at line 96 of file DTVDriftSegment.cc.

References Exception, getHistoName(), trackerHits::histo, and rootFile_.

Referenced by compute().

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

Definition at line 104 of file DTVDriftSegment.cc.

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

Referenced by getHisto().

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

Implements dtCalibration::DTVDriftBaseAlgo.

Definition at line 51 of file DTVDriftSegment.cc.

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

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

Member Data Documentation

DTResidualFitter* dtCalibration::DTVDriftSegment::fitter_
private

Definition at line 38 of file DTVDriftSegment.h.

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

const DTMtime* dtCalibration::DTVDriftSegment::mTimeMap_
private

Definition at line 36 of file DTVDriftSegment.h.

Referenced by compute(), and setES().

unsigned int dtCalibration::DTVDriftSegment::nSigmas_
private

Definition at line 34 of file DTVDriftSegment.h.

Referenced by compute().

TFile* dtCalibration::DTVDriftSegment::rootFile_
private

Definition at line 37 of file DTVDriftSegment.h.

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