CMS 3D CMS Logo

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