CMS 3D CMS Logo

DTVDriftSegment.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * \author A. Vilela Pereira
6  */
7 
8 #include "DTVDriftSegment.h"
9 
14 
19 
22 
23 #include <string>
24 #include <vector>
25 
26 #include "TH1F.h"
27 #include "TFile.h"
28 
29 using namespace std;
30 using namespace edm;
31 
32 namespace dtCalibration {
33 
34 DTVDriftSegment::DTVDriftSegment(const ParameterSet& pset):
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 }
45 
47  rootFile_->Close();
48  delete fitter_;
49 }
50 
52  // Get the map of vdrift from the setup
53  ESHandle<DTMtime> mTime;
54  setup.get<DTMtimeRcd>().get(mTime);
55  mTimeMap_ = &*mTime;
56 }
57 
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 }
95 
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 }
103 
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 }
119 
120 } // namespace
virtual void setES(const edm::EventSetup &setup)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTChamberId chamberId() const
Return the corresponding ChamberId.
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
std::string getHistoName(const DTSuperLayerId &)
virtual DTVDriftData compute(const DTSuperLayerId &)
int superLayer() const
Return the superlayer number.
#define LogTrace(id)
DTResidualFitResult fitResiduals(TH1F &histo, int nSigmas=1)
#define debug
Definition: HDRShower.cc:19
const T & get() const
Definition: EventSetup.h:56
int get(int wheelId, int stationId, int sectorId, int slId, float &mTime, float &mTrms, DTTimeUnits::type unit) const
Definition: DTMtime.cc:82
HLT enums.
int sector() const
Definition: DTChamberId.h:61
int station() const
Return the station number.
Definition: DTChamberId.h:51
TH1F * getHisto(const DTSuperLayerId &)
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45