#include <DTVDriftSegment.h>
Public Member Functions | |
virtual DTVDriftData | compute (const DTSuperLayerId &) |
DTVDriftSegment (edm::ParameterSet const &) | |
virtual void | setES (const edm::EventSetup &setup) |
virtual | ~DTVDriftSegment () |
Private Member Functions | |
TH1F * | getHisto (const DTSuperLayerId &) |
std::string | getHistoName (const DTSuperLayerId &) |
Private Attributes | |
DTResidualFitter * | fitter_ |
const DTMtime * | mTimeMap_ |
unsigned int | nSigmas_ |
TFile * | rootFile_ |
Concrete implementation of a DTVDriftBaseAlgo. Computes vDrift using fit result segment by segment.
Definition at line 22 of file DTVDriftSegment.h.
DTVDriftSegment::DTVDriftSegment | ( | edm::ParameterSet const & | pset | ) |
Definition at line 34 of file DTVDriftSegment.cc.
References debug, fitter_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), rootFile_, and dtTPAnalyzer_cfg::rootFileName.
: nSigmas_( pset.getUntrackedParameter<unsigned int>("nSigmasFitRange", 1) ) { string rootFileName = pset.getParameter<string>("rootFileName"); rootFile_ = new TFile(rootFileName.c_str(), "READ"); bool debug = pset.getUntrackedParameter<bool>("debug",false); fitter_ = new DTResidualFitter(debug); //bool debug = pset.getUntrackedParameter<bool>("debug", false); //if(debug) fitter_->setVerbosity(1); }
DTVDriftSegment::~DTVDriftSegment | ( | ) | [virtual] |
DTVDriftData DTVDriftSegment::compute | ( | const DTSuperLayerId & | slId | ) | [virtual] |
Implements 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_, dtDQMClient_cfg::resolution, ntuplemaker::status, and DTSuperLayerId::superLayer().
{ // Get original value from DB; vdrift is cm/ns , resolution is cm float vDrift = 0., resolution = 0.; int status = mTimeMap_->get(slId,vDrift,resolution,DTVelocityUnits::cm_per_ns); if(status != 0) throw cms::Exception("DTCalibration") << "Could not find vDrift entry in DB for" << slId << endl; // For RZ superlayers use original value if(slId.superLayer() == 2){ LogTrace("Calibration") << "[DTVDriftSegment]: RZ superlayer\n" << " Will use original vDrift and resolution."; return DTVDriftData(vDrift,resolution); } else{ TH1F* vDriftCorrHisto = getHisto(slId); // If empty histogram if(vDriftCorrHisto->GetEntries() == 0){ LogError("Calibration") << "[DTVDriftSegment]: Histogram " << vDriftCorrHisto->GetName() << " is empty.\n" << " Will use original vDrift and resolution."; return DTVDriftData(vDrift,resolution); } LogTrace("Calibration") << "[DTVDriftSegment]: Fitting histogram " << vDriftCorrHisto->GetName(); DTResidualFitResult fitResult = fitter_->fitResiduals(*vDriftCorrHisto,nSigmas_); LogTrace("Calibration") << "[DTVDriftSegment]: \n" << " Fit Mean = " << fitResult.fitMean << " +/- " << fitResult.fitMeanError << "\n" << " Fit Sigma = " << fitResult.fitSigma << " +/- " << fitResult.fitSigmaError; float vDriftCorr = fitResult.fitMean; float vDriftNew = vDrift*(1. - vDriftCorr); float resolutionNew = resolution; return DTVDriftData(vDriftNew,resolutionNew); } }
TH1F * DTVDriftSegment::getHisto | ( | const DTSuperLayerId & | slId | ) | [private] |
Definition at line 96 of file DTVDriftSegment.cc.
References Exception, getHistoName(), timingPdfMaker::histo, and rootFile_.
Referenced by compute().
{ string histoName = getHistoName(slId); TH1F* histo = static_cast<TH1F*>(rootFile_->Get(histoName.c_str())); if(!histo) throw cms::Exception("DTCalibration") << "v-drift correction histogram not found:" << histoName << endl; return histo; }
string DTVDriftSegment::getHistoName | ( | const DTSuperLayerId & | slId | ) | [private] |
Definition at line 104 of file DTVDriftSegment.cc.
References DTSuperLayerId::chamberId(), DTChamberId::sector(), relativeConstraints::station, DTChamberId::station(), DTSuperLayerId::superLayer(), and DTChamberId::wheel().
Referenced by getHisto().
{ DTChamberId chId = slId.chamberId(); // Compose the chamber name stringstream wheel; wheel << chId.wheel(); stringstream station; station << chId.station(); stringstream sector; sector << chId.sector(); string chHistoName = "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str(); return (slId.superLayer() != 2)?("hRPhiVDriftCorr" + chHistoName):("hRZVDriftCorr" + chHistoName); }
void DTVDriftSegment::setES | ( | const edm::EventSetup & | setup | ) | [virtual] |
Implements DTVDriftBaseAlgo.
Definition at line 51 of file DTVDriftSegment.cc.
References edm::EventSetup::get(), and mTimeMap_.
{ // Get the map of vdrift from the setup ESHandle<DTMtime> mTime; setup.get<DTMtimeRcd>().get(mTime); mTimeMap_ = &*mTime; }
DTResidualFitter* DTVDriftSegment::fitter_ [private] |
Definition at line 37 of file DTVDriftSegment.h.
Referenced by compute(), DTVDriftSegment(), and ~DTVDriftSegment().
const DTMtime* DTVDriftSegment::mTimeMap_ [private] |
Definition at line 35 of file DTVDriftSegment.h.
unsigned int DTVDriftSegment::nSigmas_ [private] |
Definition at line 33 of file DTVDriftSegment.h.
Referenced by compute().
TFile* DTVDriftSegment::rootFile_ [private] |
Definition at line 36 of file DTVDriftSegment.h.
Referenced by DTVDriftSegment(), getHisto(), and ~DTVDriftSegment().