CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ErrorsAnalyzer.cc
Go to the documentation of this file.
1 #ifndef ERRORSANALYZER_CC
2 #define ERRORSANALYZER_CC
3 
4 #include "ErrorsAnalyzer.h"
5 
7  : treeFileName_(iConfig.getParameter<std::string>("InputFileName")),
8  resolFitType_(iConfig.getParameter<int>("ResolFitType")),
9  maxEvents_(iConfig.getParameter<int>("MaxEvents")),
10  outputFileName_(iConfig.getParameter<std::string>("OutputFileName")),
11  ptBins_(iConfig.getParameter<int>("PtBins")),
12  ptMin_(iConfig.getParameter<double>("PtMin")),
13  ptMax_(iConfig.getParameter<double>("PtMax")),
14  etaBins_(iConfig.getParameter<int>("EtaBins")),
15  etaMin_(iConfig.getParameter<double>("EtaMin")),
16  etaMax_(iConfig.getParameter<double>("EtaMax")),
17  debug_(iConfig.getParameter<bool>("Debug")) {
18  parameters_ = iConfig.getParameter<std::vector<double> >("Parameters");
19  errors_ = iConfig.getParameter<std::vector<double> >("Errors");
20  errorFactors_ = iConfig.getParameter<std::vector<int> >("ErrorFactors");
21 
22  if ((parameters_.size() != errors_.size()) || (parameters_.size() != errorFactors_.size())) {
23  std::cout << "Error: parameters and errors have different number of values" << std::endl;
24  exit(1);
25  }
26 
28 
29  sigmaPtVsPt_ = new TProfile("sigmaPtVsPtProfile", "sigmaPtVsPt", ptBins_, ptMin_, ptMax_);
30  sigmaPtVsPtPlusErr_ = new TProfile("sigmaPtVsPtPlusErrProfile", "sigmaPtVsPtPlusErr", ptBins_, ptMin_, ptMax_);
31  sigmaPtVsPtMinusErr_ = new TProfile("sigmaPtVsPtMinusErrProfile", "sigmaPtVsPtMinusErr", ptBins_, ptMin_, ptMax_);
32 
33  sigmaPtVsEta_ = new TProfile("sigmaPtVsEtaProfile", "sigmaPtVsEta", etaBins_, etaMin_, etaMax_);
34  sigmaPtVsEtaPlusErr_ = new TProfile("sigmaPtVsEtaPlusErrProfile", "sigmaPtVsEtaPlusErr", etaBins_, etaMin_, etaMax_);
36  new TProfile("sigmaPtVsEtaMinusErrProfile", "sigmaPtVsEtaMinusErr", etaBins_, etaMin_, etaMax_);
37 
38  sigmaMassVsPt_ = new TProfile("sigmaMassVsPtProfile", "sigmaMassVsPt", ptBins_, ptMin_, ptMax_);
39  sigmaMassVsPtPlusErr_ = new TProfile("sigmaMassVsPtPlusErrProfile", "sigmaMassVsPtPlusErr", ptBins_, ptMin_, ptMax_);
41  new TProfile("sigmaMassVsPtMinusErrProfile", "sigmaMassVsPtMinusErr", ptBins_, ptMin_, ptMax_);
42 
43  sigmaMassVsEta_ = new TProfile("sigmaMassVsEtaProfile", "sigmaMassVsEta", etaBins_, etaMin_, etaMax_);
45  new TProfile("sigmaMassVsEtaPlusErrProfile", "sigmaMassVsEtaPlusErr", etaBins_, etaMin_, etaMax_);
47  new TProfile("sigmaMassVsEtaMinusErrProfile", "sigmaMassVsEtaMinusErr", etaBins_, etaMin_, etaMax_);
48 }
49 
51  valuePlusError_.resize(parameters_.size());
52  valueMinusError_.resize(parameters_.size());
53 
54  std::vector<double>::const_iterator parIt = parameters_.begin();
55  std::vector<double>::const_iterator errIt = errors_.begin();
56  std::vector<int>::const_iterator errFactorIt = errorFactors_.begin();
57  int i = 0;
58  for (; parIt != parameters_.end(); ++parIt, ++errIt, ++errFactorIt, ++i) {
59  valuePlusError_[i] = *parIt + (*errIt) * (*errFactorIt);
60  valueMinusError_[i] = *parIt - (*errIt) * (*errFactorIt);
61  }
62 }
63 
65  gROOT->SetStyle("Plain");
66 
68 
69  TFile* outputFile = new TFile(outputFileName_, "RECREATE");
70 
73 
76 
77  outputFile->Write();
78  outputFile->Close();
79 }
80 
82  const TProfile* histoPlusErr,
83  const TProfile* histoMinusErr,
84  const TString& type) {
85  TH1D* sigmaPtVsEtaTH1D =
86  new TH1D(type, type, histo->GetNbinsX(), histo->GetXaxis()->GetXmin(), histo->GetXaxis()->GetXmax());
87  TH1D* sigmaPtVsEtaPlusErrTH1D = new TH1D(type + "PlusErr",
88  type + "PlusErr",
89  histo->GetNbinsX(),
90  histo->GetXaxis()->GetXmin(),
91  histo->GetXaxis()->GetXmax());
92  TH1D* sigmaPtVsEtaMinusErrTH1D = new TH1D(type + "MinusErr",
93  type + "MinusErr",
94  histo->GetNbinsX(),
95  histo->GetXaxis()->GetXmin(),
96  histo->GetXaxis()->GetXmax());
97 
98  TH1D* sigmaMassVsEtaTH1D =
99  new TH1D(type, type, histo->GetNbinsX(), histo->GetXaxis()->GetXmin(), histo->GetXaxis()->GetXmax());
100  TH1D* sigmaMassVsEtaPlusErrTH1D = new TH1D(type + "PlusErr",
101  type + "PlusErr",
102  histo->GetNbinsX(),
103  histo->GetXaxis()->GetXmin(),
104  histo->GetXaxis()->GetXmax());
105  TH1D* sigmaMassVsEtaMinusErrTH1D = new TH1D(type + "MinusErr",
106  type + "MinusErr",
107  histo->GetNbinsX(),
108  histo->GetXaxis()->GetXmin(),
109  histo->GetXaxis()->GetXmax());
110 
111  TCanvas* canvas = new TCanvas("canvas" + type, "canvas" + type, 1000, 800);
112  for (int iBin = 1; iBin <= histo->GetNbinsX(); ++iBin) {
113  sigmaPtVsEtaTH1D->SetBinContent(iBin, histo->GetBinContent(iBin));
114  sigmaPtVsEtaPlusErrTH1D->SetBinContent(iBin, histoPlusErr->GetBinContent(iBin));
115  sigmaPtVsEtaMinusErrTH1D->SetBinContent(iBin, histoMinusErr->GetBinContent(iBin));
116  }
117  int numBins = sigmaPtVsEtaTH1D->GetNbinsX();
118  // Draw TGraph with asymmetric errors
119  Double_t* values = sigmaPtVsEtaTH1D->GetArray();
120  Double_t* valuesPlus = sigmaPtVsEtaPlusErrTH1D->GetArray();
121  Double_t* valuesMinus = sigmaPtVsEtaMinusErrTH1D->GetArray();
122  double* posErrors = new double[numBins];
123  double* negErrors = new double[numBins];
124 
125  TGraphAsymmErrors* graphAsymmErrors = new TGraphAsymmErrors(sigmaPtVsEtaTH1D);
126  TGraph* graph = new TGraph(sigmaPtVsEtaTH1D);
127 
128  for (int i = 1; i <= numBins; ++i) {
129  // std::cout << "filling " << i << std::endl;
130  posErrors[i - 1] = valuesPlus[i] - values[i];
131  if (valuesMinus[i] < 0)
132  negErrors[i - 1] = values[i];
133  else
134  negErrors[i - 1] = values[i] - valuesMinus[i];
135 
136  graphAsymmErrors->SetPointEYlow(i - 1, negErrors[i - 1]);
137  graphAsymmErrors->SetPointEYhigh(i - 1, posErrors[i - 1]);
138  }
139 
140  canvas->Draw();
141 
142  delete[] posErrors;
143  delete[] negErrors;
144 
145  graphAsymmErrors->SetTitle("");
146  graphAsymmErrors->SetFillColor(kGray);
147  graphAsymmErrors->Draw("A2");
148  TString title(type);
149  if (type == "Eta")
150  title = "#eta";
151  graphAsymmErrors->GetXaxis()->SetTitle(title);
152  graphAsymmErrors->GetYaxis()->SetTitle("#sigmaPt/Pt");
153  graph->Draw();
154 
155  // graph->SetLineColor(kGray);
156  // graph->SetMarkerColor(kBlack);
157  // graph->Draw("AP");
158 
159  // if( debug_ ) {
160  // sigmaPtVsEtaPlusErrTH1D->SetLineColor(kRed);
161  // sigmaPtVsEtaMinusErrTH1D->SetLineColor(kBlue);
162  // }
163  // else {
164  // sigmaPtVsEtaPlusErrTH1D->SetFillColor(kGray);
165  // sigmaPtVsEtaPlusErrTH1D->SetLineColor(kWhite);
166  // sigmaPtVsEtaMinusErrTH1D->SetFillColor(kWhite);
167  // sigmaPtVsEtaMinusErrTH1D->SetLineColor(kWhite);
168  // }
169  // sigmaPtVsEtaPlusErrTH1D->Draw();
170  // sigmaPtVsEtaTH1D->Draw("SAME");
171  // sigmaPtVsEtaMinusErrTH1D->Draw("SAME");
172 
173  sigmaPtVsEtaPlusErrTH1D->Write();
174  sigmaPtVsEtaTH1D->Write();
175  sigmaPtVsEtaMinusErrTH1D->Write();
176 
177  sigmaPtVsEtaPlusErr_->Write();
178  sigmaPtVsEta_->Write();
179  sigmaPtVsEtaMinusErr_->Write();
180 
181  // Mass
182  sigmaMassVsEtaPlusErrTH1D->Write();
183  sigmaMassVsEtaTH1D->Write();
184  sigmaMassVsEtaMinusErrTH1D->Write();
185 
186  sigmaMassVsEtaPlusErr_->Write();
187  sigmaMassVsEta_->Write();
188  sigmaMassVsEtaMinusErr_->Write();
189 
190  canvas->Write();
191 }
192 
193 // ------------ method called to for each event ------------
195 
197  std::cout << "Reading muon pairs from Root Tree in " << treeFileName_ << std::endl;
198  RootTreeHandler rootTreeHandler;
199 
200  typedef std::vector<std::pair<lorentzVector, lorentzVector> > MuonPairVector;
201  MuonPairVector savedPair;
202  std::vector<std::pair<unsigned int, unsigned long long> > evtRun;
203  rootTreeHandler.readTree(maxEvents_, treeFileName_, &savedPair, 0, &evtRun);
204  // rootTreeHandler.readTree(maxEvents, inputRootTreeFileName_, &savedPair, &(MuScleFitUtils::genPair));
205 
210  MuScleFitUtils::resfind = std::vector<int>(6, 0);
211 
212  // Loop on all the pairs
213  unsigned int i = 0;
214  MuonPairVector::iterator it = savedPair.begin();
215  std::cout << "Starting loop on " << savedPair.size() << " muons" << std::endl;
216  for (; it != savedPair.end(); ++it, ++i) {
217  double pt1 = it->first.pt();
218  double eta1 = it->first.eta();
219  double pt2 = it->second.pt();
220  double eta2 = it->second.eta();
221 
222  if (debug_) {
223  std::cout << "pt1 = " << pt1 << ", eta1 = " << eta1 << ", pt2 = " << pt2 << ", eta2 = " << eta2 << std::endl;
224  }
225 
226  if (pt1 == 0 && pt2 == 0 && eta1 == 0 && eta2 == 0)
227  continue;
228 
229  double sigmaPt1 = resolutionFunctionForVec->sigmaPt(pt1, eta1, parameters_);
230  double sigmaPt2 = resolutionFunctionForVec->sigmaPt(pt2, eta2, parameters_);
231  double sigmaPtPlusErr1 = resolutionFunctionForVec->sigmaPt(pt1, eta1, valuePlusError_);
232  double sigmaPtPlusErr2 = resolutionFunctionForVec->sigmaPt(pt2, eta2, valuePlusError_);
233  double sigmaPtMinusErr1 = resolutionFunctionForVec->sigmaPt(pt1, eta1, valueMinusError_);
234  double sigmaPtMinusErr2 = resolutionFunctionForVec->sigmaPt(pt2, eta2, valueMinusError_);
235 
236  double sigmaMass = MuScleFitUtils::massResolution(it->first, it->second, parameters_);
237  double sigmaMassPlusErr = MuScleFitUtils::massResolution(it->first, it->second, valuePlusError_);
238  double sigmaMassMinusErr = MuScleFitUtils::massResolution(it->first, it->second, valueMinusError_);
239 
240  if (debug_) {
241  std::cout << "sigmaPt1 = " << sigmaPt1 << " + " << sigmaPtPlusErr1 << " - " << sigmaPtMinusErr1 << std::endl;
242  std::cout << "sigmaPt2 = " << sigmaPt2 << " + " << sigmaPtPlusErr2 << " - " << sigmaPtMinusErr2 << std::endl;
243  std::cout << "sigmaMass = " << sigmaMass << " + " << sigmaMassPlusErr << " - " << sigmaMassMinusErr << std::endl;
244  }
245 
246  // Protections from nans
247  if (pt1 != pt1)
248  continue;
249  if (pt2 != pt2)
250  continue;
251  if (eta1 != eta1)
252  continue;
253  if (eta2 != eta2)
254  continue;
255  if (sigmaPt1 != sigmaPt1)
256  continue;
257  if (sigmaPt2 != sigmaPt2)
258  continue;
259  if (sigmaPtPlusErr1 != sigmaPtPlusErr1)
260  continue;
261  if (sigmaPtPlusErr2 != sigmaPtPlusErr2)
262  continue;
263  if (sigmaPtMinusErr1 != sigmaPtMinusErr1)
264  continue;
265  if (sigmaPtMinusErr2 != sigmaPtMinusErr2)
266  continue;
267  if (sigmaMass != sigmaMass)
268  continue;
269  if (sigmaMassPlusErr != sigmaMassPlusErr)
270  continue;
271  if (sigmaMassMinusErr != sigmaMassMinusErr)
272  continue;
273 
274  std::cout << "Muon pair number " << i << std::endl;
275 
276  // std::cout << "sigmaPt1 = " << sigmaPt1 << ", sigmaPt2 = " << sigmaPt2 << std::endl;
277  // std::cout << "sigmaPtPlusErr1 = " << sigmaPtPlusErr1 << ", sigmaPtMinusErr2 = " << sigmaPtMinusErr2 << std::endl;
278  // std::cout << "sigmaPtMinusErr1 = " << sigmaPtPlusErr1 << ", sigmaPtMinusErr2 = " << sigmaPtMinusErr2 << std::endl;
279 
280  sigmaPtVsPt_->Fill(pt1, sigmaPt1);
281  sigmaPtVsPt_->Fill(pt2, sigmaPt2);
282  sigmaPtVsPtPlusErr_->Fill(pt1, sigmaPtPlusErr1);
283  sigmaPtVsPtPlusErr_->Fill(pt2, sigmaPtPlusErr2);
284  sigmaPtVsPtMinusErr_->Fill(pt1, sigmaPtMinusErr1);
285  sigmaPtVsPtMinusErr_->Fill(pt2, sigmaPtMinusErr2);
286 
287  sigmaPtVsEta_->Fill(eta1, sigmaPt1);
288  sigmaPtVsEta_->Fill(eta2, sigmaPt2);
289  sigmaPtVsEtaPlusErr_->Fill(eta1, sigmaPtPlusErr1);
290  sigmaPtVsEtaPlusErr_->Fill(eta2, sigmaPtPlusErr2);
291  sigmaPtVsEtaMinusErr_->Fill(eta1, sigmaPtMinusErr1);
292  sigmaPtVsEtaMinusErr_->Fill(eta2, sigmaPtMinusErr2);
293 
294  sigmaMassVsPt_->Fill(pt1, sigmaMass);
295  sigmaMassVsPt_->Fill(pt2, sigmaMass);
296  sigmaMassVsPtPlusErr_->Fill(pt1, sigmaMassPlusErr);
297  sigmaMassVsPtPlusErr_->Fill(pt2, sigmaMassPlusErr);
298  sigmaMassVsPtMinusErr_->Fill(pt1, sigmaMassMinusErr);
299  sigmaMassVsPtMinusErr_->Fill(pt2, sigmaMassMinusErr);
300 
301  sigmaMassVsEta_->Fill(eta1, sigmaMass);
302  sigmaMassVsEta_->Fill(eta2, sigmaMass);
303  sigmaMassVsEtaPlusErr_->Fill(eta1, sigmaMassPlusErr);
304  sigmaMassVsEtaPlusErr_->Fill(eta2, sigmaMassPlusErr);
305  sigmaMassVsEtaMinusErr_->Fill(eta1, sigmaMassMinusErr);
306  sigmaMassVsEtaMinusErr_->Fill(eta2, sigmaMassMinusErr);
307  }
308 }
309 
310 //define this as a plug-in
312 
313 #endif
TProfile * sigmaMassVsPtPlusErr_
TProfile * sigmaMassVsEta_
~ErrorsAnalyzer() override
std::vector< double > valueMinusError_
std::vector< double > errors_
void drawHistograms(const TProfile *histo, const TProfile *histoPlusErr, const TProfile *histoMinusErr, const TString &type)
TProfile * sigmaPtVsPtMinusErr_
TProfile * sigmaMassVsPtMinusErr_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static bool debugMassResol_
TProfile * sigmaMassVsEtaPlusErr_
void analyze(const edm::Event &, const edm::EventSetup &) override
static int debug
TString outputFileName_
void readTree(const int maxEvents, const TString &fileName, MuonPairVector *savedPair, const int muonType, std::vector< std::pair< unsigned int, unsigned long long > > *evtRun, MuonPairVector *genPair=nullptr)
etaMax_(conf.getParameter< double >("etaMax"))
TProfile * sigmaPtVsEtaMinusErr_
def canvas
Definition: svgfig.py:482
resolutionFunctionBase< double * > * resolutionFunctionService(const int identifier)
Service to build the resolution functor corresponding to the passed identifier.
Definition: Functions.cc:70
virtual double sigmaPt(const double &pt, const double &eta, const T &parval)=0
ptMin_(conf.getParameter< double >("ptMin"))
int iEvent
Definition: GenABIO.cc:224
resolutionFunctionBase< std::vector< double > > * resolutionFunctionVecService(const int identifier)
Service to build the resolution functor corresponding to the passed identifier when receiving a std::...
Definition: Functions.cc:92
static double massResolution(const lorentzVector &mu1, const lorentzVector &mu2)
std::vector< int > errorFactors_
uint32_t maxEvents_
std::vector< double > valuePlusError_
TProfile * sigmaPtVsPtPlusErr_
TProfile * sigmaPtVsPt_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< std::pair< lorentzVector, lorentzVector > > MuonPairVector
TProfile * sigmaMassVsEtaMinusErr_
TProfile * sigmaMassVsPt_
std::vector< double > parameters_
tuple cout
Definition: gather_cfg.py:144
ErrorsAnalyzer(const edm::ParameterSet &)
static std::vector< int > resfind
TString treeFileName_
TProfile * sigmaPtVsEta_
etaMin_(conf.getParameter< double >("etaMin"))
TProfile * sigmaPtVsEtaPlusErr_
static resolutionFunctionBase< double * > * resolutionFunction