Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:02

0001 void CSCSegmentVisualise(int nHisto){
0002 
0003   /* Macro to plot histograms produced by CSCRecHitVisualise.cc
0004    * You may need to update the TFile name, and will need to
0005    * input the segtype as shown below.
0006    *
0007    * Author:  Dominique Fortin - UCR
0008    */
0009 
0010 
0011 // Files for histogram output --> set suffixps to desired file type:  e.g. .eps, .jpg, ...
0012 
0013 TString suffixps = ".gif";
0014 
0015  TString endfile = ".root";
0016  TString tfile = "cscsegments_plot"+endfile;
0017 
0018  TFile *file = TFile::Open(tfile);
0019 
0020 for (int i = nHisto; i < nHisto+1; ++i ) {
0021 
0022   // ********************************************************************
0023   // Pointers to histograms
0024   // ********************************************************************
0025   char a[4];
0026   int j = i + 100;
0027  
0028   sprintf(a, "h%d", i+100);
0029   TString idx = a;
0030 
0031   TString plot1 = "x_vs_z_"+idx+suffixps;
0032   TString plot2 = "y_vs_z_"+idx+suffixps;
0033 
0034   // 1) rechit x vs z
0035   sprintf(a, "h%d", i+100);
0036   hxvsz = (TH2F *) file->Get(a);
0037  
0038   // 2) rechit y vs z
0039   sprintf(a, "h%d", i+200);
0040   hyvsz = (TH2F *) file->Get(a);
0041 
0042   // 3) used hit on seg x vs z
0043   sprintf(a, "h%d", i+300);
0044   hxvszS = (TH2F *) file->Get(a);
0045 
0046   // 4) used hit on seg y vs z
0047   sprintf(a, "h%d", i+400);
0048   hyvszS = (TH2F *) file->Get(a);
0049 
0050   // 5) Projected segment x vs z
0051   sprintf(a, "h%d", i+500);
0052   hxvszP = (TH2F *) file->Get(a);
0053 
0054   // 6) Projected segment y vs z
0055   sprintf(a, "h%d", i+600);
0056   hyvszP = (TH2F *) file->Get(a);
0057 
0058   // 7) rechit from electron x vs z
0059   sprintf(a, "h%d", i+700);
0060   hxvszE = (TH2F *) file->Get(a);
0061  
0062   // 8) rechit from electron y vs z
0063   sprintf(a, "h%d", i+800);
0064   hyvszE = (TH2F *) file->Get(a);
0065 
0066 
0067   // Make plot
0068   gStyle->SetOptStat(kFALSE);
0069   TCanvas *c1 = new TCanvas("c1","");
0070   gStyle->SetOptStat(kFALSE);
0071   c1->SetFillColor(10);
0072   c1->SetFillColor(10);
0073   // segments
0074   hxvszP->SetMarkerSize(0.2);
0075   hxvszP->SetMarkerStyle(6);
0076   hxvszP->SetMarkerColor(kRed);
0077 //  hxvszP->SetTitle("X-local vs Z-local");
0078   hxvszP->GetXaxis()->SetTitle("local z (cm)");
0079   hxvszP->GetYaxis()->SetTitle("local x (cm)");
0080   hxvszP->Draw();
0081   // rechits
0082   hxvsz->SetMarkerSize(3);
0083   hxvsz->SetMarkerStyle(30);
0084   hxvsz->SetMarkerColor(kBlack);
0085   hxvsz->Draw("SAME");
0086   // rechits from electrons
0087   hxvszE->SetMarkerSize(3);
0088   hxvszE->SetMarkerStyle(29);
0089   hxvszE->SetMarkerColor(kGreen);
0090   hxvszE->Draw("SAME");
0091   // used hits on segments
0092   hxvszS->SetMarkerSize(2);
0093   hxvszS->SetMarkerStyle(29);
0094   hxvszS->SetMarkerColor(kBlue);
0095   hxvszS->Draw("SAME");
0096   c1->Print(plot1);  
0097 
0098 
0099   // Make plot
0100   gStyle->SetOptStat(kFALSE);
0101   TCanvas *c1 = new TCanvas("c1","");
0102   gStyle->SetOptStat(kFALSE);
0103   c1->SetFillColor(10);
0104   c1->SetFillColor(10);
0105   // segments
0106   hyvszP->SetMarkerSize(0.2);
0107   hyvszP->SetMarkerStyle(6);
0108   hyvszP->SetMarkerColor(kRed);
0109 //  hyvszP->SetTitle("Y-local vs Z-local");
0110   hyvszP->GetXaxis()->SetTitle("local z (cm)");
0111   hyvszP->GetYaxis()->SetTitle("local y (cm)");
0112   hyvszP->Draw();
0113   // rechits
0114   hyvsz->SetMarkerSize(3);
0115   hyvsz->SetMarkerStyle(30);
0116   hyvsz->SetMarkerColor(kBlack);
0117   hyvsz->Draw("SAME");
0118   // rechits from electrons
0119   hyvszE->SetMarkerSize(3);
0120   hyvszE->SetMarkerStyle(29);
0121   hyvszE->SetMarkerColor(kGreen);
0122   hyvszE->Draw("SAME");
0123   // used hits on segments
0124   hyvszS->SetMarkerSize(2);
0125   hyvszS->SetMarkerStyle(29);
0126   hyvszS->SetMarkerColor(kBlue);
0127   hyvszS->Draw("SAME");
0128   c1->Print(plot2);  
0129 
0130 }
0131  
0132 gROOT->ProcessLine(".q");
0133 
0134 }