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  string rootFileName = pset.getParameter<string>("rootFileName");
37  rootFile_ = new TFile(rootFileName.c_str(), "READ");
38 
39  bool debug = pset.getUntrackedParameter<bool>("debug", false);
40  fitter_ = new DTResidualFitter(debug);
41  //bool debug = pset.getUntrackedParameter<bool>("debug", false);
42  //if(debug) fitter_->setVerbosity(1);
43  }
44 
46  rootFile_->Close();
47  delete fitter_;
48  }
49 
51  // Get the map of vdrift from the setup
52  ESHandle<DTMtime> mTime;
53  setup.get<DTMtimeRcd>().get(mTime);
54  mTimeMap_ = &*mTime;
55  }
56 
58  // Get original value from DB; vdrift is cm/ns , resolution is cm
59  float vDrift = 0., resolution = 0.;
61 
62  if (status != 0)
63  throw cms::Exception("DTCalibration") << "Could not find vDrift entry in DB for" << slId << endl;
64 
65  // For RZ superlayers use original value
66  if (slId.superLayer() == 2) {
67  LogTrace("Calibration") << "[DTVDriftSegment]: RZ superlayer\n"
68  << " Will use original vDrift and resolution.";
69  return DTVDriftData(vDrift, resolution);
70  } else {
71  TH1F* vDriftCorrHisto = getHisto(slId);
72  // If empty histogram
73  if (vDriftCorrHisto->GetEntries() == 0) {
74  LogError("Calibration") << "[DTVDriftSegment]: Histogram " << vDriftCorrHisto->GetName() << " is empty.\n"
75  << " Will use original vDrift and resolution.";
76  return DTVDriftData(vDrift, resolution);
77  }
78 
79  LogTrace("Calibration") << "[DTVDriftSegment]: Fitting histogram " << vDriftCorrHisto->GetName();
80  DTResidualFitResult fitResult = fitter_->fitResiduals(*vDriftCorrHisto, nSigmas_);
81  LogTrace("Calibration") << "[DTVDriftSegment]: \n"
82  << " Fit Mean = " << fitResult.fitMean << " +/- " << fitResult.fitMeanError << "\n"
83  << " Fit Sigma = " << fitResult.fitSigma << " +/- " << fitResult.fitSigmaError;
84 
85  float vDriftCorr = fitResult.fitMean;
86  float vDriftNew = vDrift * (1. - vDriftCorr);
87  float resolutionNew = resolution;
88  return DTVDriftData(vDriftNew, resolutionNew);
89  }
90  }
91 
93  string histoName = getHistoName(slId);
94  TH1F* histo = static_cast<TH1F*>(rootFile_->Get(histoName.c_str()));
95  if (!histo)
96  throw cms::Exception("DTCalibration") << "v-drift correction histogram not found:" << histoName << endl;
97  return histo;
98  }
99 
101  DTChamberId chId = slId.chamberId();
102 
103  // Compose the chamber name
104  std::string wheel = std::to_string(chId.wheel());
105  std::string station = std::to_string(chId.station());
106  std::string sector = std::to_string(chId.sector());
107 
108  string chHistoName = "_W" + wheel + "_St" + station + "_Sec" + sector;
109 
110  return (slId.superLayer() != 2) ? ("hRPhiVDriftCorr" + chHistoName) : ("hRZVDriftCorr" + chHistoName);
111  }
112 
113 } // namespace dtCalibration
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTChamberId chamberId() const
Return the corresponding ChamberId.
std::string getHistoName(const DTSuperLayerId &)
void setES(const edm::EventSetup &setup) override
int superLayer() const
Return the superlayer number.
#define LogTrace(id)
DTResidualFitResult fitResiduals(TH1F &histo, int nSigmas=1)
#define debug
Definition: HDRShower.cc:19
int get(int wheelId, int stationId, int sectorId, int slId, float &mTime, float &mTrms, DTTimeUnits::type unit) const
Definition: DTMtime.cc:56
HLT enums.
int sector() const
Definition: DTChamberId.h:49
T get() const
Definition: EventSetup.h:73
int station() const
Return the station number.
Definition: DTChamberId.h:42
TH1F * getHisto(const DTSuperLayerId &)
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
DTVDriftData compute(const DTSuperLayerId &) override