Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:34

0001 #include "multibsvspvplots.h"
0002 #include <vector>
0003 #include <algorithm>
0004 #include <string>
0005 #include <map>
0006 #include <iostream>
0007 #include "TPad.h"
0008 #include "TH1D.h"
0009 #include "TFile.h"
0010 #include "TH2F.h"
0011 #include "TH1F.h"
0012 #include "TF1.h"
0013 #include "TProfile.h"
0014 #include "TProfile2D.h"
0015 #include "DPGAnalysis/SiStripTools/interface/CommonAnalyzer.h"
0016 #include "TCanvas.h"
0017 #include "TSystem.h"
0018 #include "TStyle.h"
0019 #include "TText.h"
0020 #include "TLegend.h"
0021 
0022 void multibsvspvplots(const char* module, const char* label, const char* postfix, const char* /* family="" */) {
0023   TFile* file[5];
0024 
0025   file[0] = new TFile("rootfiles/Tracking_PFG_Run2012A_prompt_minbias_v1_190456-194076_muon_relumi_v55_fittedV0.root");
0026   file[1] = new TFile("rootfiles/Tracking_PFG_Run2012B_prompt_minbias_v1_190456-196531_muon_relumi_v55_fittedV0.root");
0027   file[2] = new TFile("rootfiles/Tracking_PFG_Run2012C_prompt_minbias_v1_190456-199011_muon_relumi_v55_fittedV0.root");
0028   file[3] = new TFile("rootfiles/Tracking_PFG_Run2012C_prompt_minbias_v2_190456-203002_muon_relumi_v57_fittedV0.root");
0029   file[4] = new TFile("rootfiles/Tracking_PFG_Run2012D_prompt_minbias_v1_190456-208686_muon_relumi_v57_fittedV0.root");
0030 
0031   char modfull[300];
0032   sprintf(modfull, "%s%s", module, postfix);
0033 
0034   char labfull[300];
0035   sprintf(labfull, "%s%s", label, postfix);
0036 
0037   // Summary histograms
0038   TH1D* deltaxsum = new TH1D("deltaxsum", "(PV-BS) Fitted X position vs run", 10, 0., 10.);
0039   deltaxsum->SetCanExtend(TH1::kAllAxes);
0040   TH1D* deltaysum = new TH1D("deltaysum", "(PV-BS) Fitted Y position vs run", 10, 0., 10.);
0041   deltaysum->SetCanExtend(TH1::kAllAxes);
0042 
0043   TH1D* deltaxmeansum = new TH1D("deltaxmeansum", "(PV-BS) Mean X position vs run", 10, 0., 10.);
0044   deltaxmeansum->SetCanExtend(TH1::kAllAxes);
0045   TH1D* deltaymeansum = new TH1D("deltaymeansum", "(PV-BS) Mean Y position vs run", 10, 0., 10.);
0046   deltaymeansum->SetCanExtend(TH1::kAllAxes);
0047 
0048   TH1D* deltazmeansum = new TH1D("deltazmeansum", "(PV-BS) Mean Z position vs run", 10, 0., 10.);
0049   deltazmeansum->SetCanExtend(TH1::kAllAxes);
0050 
0051   TF1* fdoubleg = new TF1("doubleg",
0052                           "[1]*exp(-0.5*((x-[0])/[2])**2)+[3]*exp(-0.5*((x-[0])/[4])**2)+[5]*exp(sqrt((x-[0])**2)*[6])",
0053                           -.2,
0054                           .2);
0055   fdoubleg->SetLineColor(kRed);
0056   fdoubleg->SetLineWidth(1);
0057 
0058   for (unsigned int ifile = 0; ifile < 5; ++ifile) {
0059     CommonAnalyzer castat(file[ifile], "", modfull);
0060 
0061     std::vector<unsigned int> runs = castat.getRunList();
0062     std::sort(runs.begin(), runs.end());
0063 
0064     {
0065       std::cout << "Found " << runs.size() << " runs" << std::endl;
0066 
0067       for (unsigned int i = 0; i < runs.size(); ++i) {
0068         char runlabel[100];
0069         sprintf(runlabel, "%d", runs[i]);
0070         char runpath[100];
0071         sprintf(runpath, "run_%d", runs[i]);
0072         castat.setPath(runpath);
0073         std::cout << runpath << std::endl;
0074 
0075         TH1F* deltax = (TH1F*)castat.getObject("deltaxrun");
0076         if (deltax && deltax->GetEntries() > 0) {
0077           fdoubleg->SetParameter(0, deltax->GetMean());
0078           fdoubleg->SetParameter(2, deltax->GetRMS());
0079           fdoubleg->SetParameter(4, deltax->GetRMS());
0080           fdoubleg->SetParameter(1, deltax->GetMaximum());
0081           fdoubleg->SetParameter(3, 0.1 * deltax->GetMaximum());
0082           fdoubleg->SetParameter(5, 0.1 * deltax->GetMaximum());
0083           /* const int result = */ deltax->Fit(fdoubleg, "q0b", "", -.05, .05);
0084 
0085           int bin = deltaxsum->Fill(runlabel, fdoubleg->GetParameter(0));
0086           deltaxsum->SetBinError(bin, fdoubleg->GetParError(0));
0087 
0088           bin = deltaxmeansum->Fill(runlabel, deltax->GetMean());
0089           deltaxmeansum->SetBinError(bin, deltax->GetMeanError());
0090         }
0091         TH1F* deltay = (TH1F*)castat.getObject("deltayrun");
0092         if (deltay && deltay->GetEntries() > 0) {
0093           fdoubleg->SetParameter(0, deltay->GetMean());
0094           fdoubleg->SetParameter(2, deltay->GetRMS());
0095           fdoubleg->SetParameter(4, deltay->GetRMS());
0096           fdoubleg->SetParameter(1, deltay->GetMaximum());
0097           fdoubleg->SetParameter(3, 0.1 * deltay->GetMaximum());
0098           fdoubleg->SetParameter(5, 0.1 * deltay->GetMaximum());
0099           /* const int result = */ deltay->Fit(fdoubleg, "q0b", "", -.05, .05);
0100           int bin = deltaysum->Fill(runlabel, fdoubleg->GetParameter(0));
0101           deltaysum->SetBinError(bin, fdoubleg->GetParError(0));
0102 
0103           bin = deltaymeansum->Fill(runlabel, deltay->GetMean());
0104           deltaymeansum->SetBinError(bin, deltay->GetMeanError());
0105         }
0106         TH1F* deltaz = (TH1F*)castat.getObject("deltazrun");
0107         if (deltaz && deltaz->GetEntries() > 0) {
0108           int bin = deltazmeansum->Fill(runlabel, deltaz->GetMean());
0109           deltazmeansum->SetBinError(bin, deltaz->GetMeanError());
0110         }
0111       }
0112     }
0113   }
0114 
0115   std::string plotfilename;
0116 
0117   plotfilename = "/afs/cern.ch/cms/tracking/output/";
0118   //    plotfilename += dirname;
0119   //    plotfilename += "/";
0120   //  plotfilename += labfull;
0121   plotfilename += "/deltaxsum_";
0122   plotfilename += labfull;
0123   //    plotfilename += "_";
0124   //    plotfilename += dirname;
0125   plotfilename += ".gif";
0126 
0127   /* TCanvas * cwidedeltax = */ new TCanvas(plotfilename.c_str(), plotfilename.c_str(), 1500, 500);
0128 
0129   deltaxsum->SetLineColor(kRed);
0130   deltaxsum->SetMarkerColor(kRed);
0131   deltaxsum->GetYaxis()->SetRangeUser(-.002, .002);
0132   deltaxsum->GetYaxis()->SetTitle("#Delta x (cm)");
0133   deltaxsum->Draw();
0134   deltaxmeansum->Draw("esame");
0135   TLegend deltaxleg(.7, .8, .85, .9, "#Delta(x)");
0136   deltaxleg.AddEntry(deltaxsum, "fitted mean", "l");
0137   deltaxleg.AddEntry(deltaxmeansum, "aritm. mean", "l");
0138   deltaxleg.Draw();
0139   gPad->Print(plotfilename.c_str());
0140   //  delete cwidedeltax;
0141 
0142   plotfilename = "/afs/cern.ch/cms/tracking/output/";
0143   //    plotfilename += dirname;
0144   //    plotfilename += "/";
0145   //  plotfilename += labfull;
0146   plotfilename += "/deltaysum_";
0147   plotfilename += labfull;
0148   //    plotfilename += "_";
0149   //    plotfilename += dirname;
0150   plotfilename += ".gif";
0151 
0152   /* TCanvas * cwidedeltay = */ new TCanvas(plotfilename.c_str(), plotfilename.c_str(), 1500, 500);
0153 
0154   deltaysum->SetLineColor(kRed);
0155   deltaysum->SetMarkerColor(kRed);
0156   deltaysum->GetYaxis()->SetRangeUser(-.002, .002);
0157   deltaysum->GetYaxis()->SetTitle("#Delta y (cm)");
0158   deltaysum->Draw();
0159   deltaymeansum->Draw("esame");
0160   TLegend deltayleg(.7, .8, .85, .9, "#Delta(y)");
0161   deltayleg.AddEntry(deltaysum, "fitted mean", "l");
0162   deltayleg.AddEntry(deltaymeansum, "aritm. mean", "l");
0163   deltayleg.Draw();
0164   gPad->Print(plotfilename.c_str());
0165   //  delete cwidedeltay;
0166 
0167   plotfilename = "/afs/cern.ch/cms/tracking/output/";
0168   //    plotfilename += dirname;
0169   //    plotfilename += "/";
0170   //  plotfilename += labfull;
0171   plotfilename += "/deltazsum_";
0172   plotfilename += labfull;
0173   //    plotfilename += "_";
0174   //    plotfilename += dirname;
0175   plotfilename += ".gif";
0176 
0177   /* TCanvas * cwidedeltaz = */ new TCanvas(plotfilename.c_str(), plotfilename.c_str(), 1500, 500);
0178 
0179   deltazmeansum->GetYaxis()->SetRangeUser(-2., 2.);
0180   deltazmeansum->GetYaxis()->SetTitle("#Delta z (cm)");
0181   deltazmeansum->Draw();
0182   gPad->Print(plotfilename.c_str());
0183   //  delete cwidedeltaz;
0184   /*
0185   delete deltaxsum;
0186   delete deltaysum;
0187   delete deltaxmeansum;
0188   delete deltaymeansum;
0189   delete deltazmeansum;
0190   */
0191   delete fdoubleg;
0192 }