File indexing completed on 2025-08-05 08:15:07
0001
0002 #include <TFile.h>
0003 #include <TGraphAsymmErrors.h>
0004 #include <TGraphErrors.h>
0005 #include <TH2.h>
0006 #include <TLatex.h>
0007 #include <TLegend.h>
0008 #include <TString.h>
0009 #include <TTree.h>
0010 #include <cassert>
0011 #include <cmath>
0012 #include "SaveCanvas.C"
0013 #include "sPhenixStyle.C"
0014
0015 TFile *_file0 = NULL;
0016 TString description;
0017
0018 void DrawTPCIntegratedCharge(
0019
0020
0021
0022 const TString infile = "/sphenix/user/jinhuang/TPC/Multiplicity/AuAu200MB_170kHz_Iter2/AuAu200MB_170kHz_Iter2_SUM.xml_TPCIntegratedCharge.root",
0023 const TString disc = "Au+Au MB Triggered + 170 kHz collision"
0024 )
0025 {
0026 SetsPhenixStyle();
0027 gStyle->SetOptStat(0);
0028 gStyle->SetOptFit(1111);
0029 TVirtualFitter::SetDefaultFitter("Minuit2");
0030
0031 _file0 = new TFile(infile);
0032 assert(_file0->IsOpen());
0033 description = disc;
0034
0035 vector<int> layers;
0036 layers.push_back(1);
0037 layers.push_back(16);
0038 layers.push_back(17);
0039 layers.push_back(32);
0040 layers.push_back(33);
0041 layers.push_back(48);
0042
0043 Check();
0044
0045 ChargePerLayer(layers);
0046 }
0047
0048 void ChargePerLayer(vector<int> layers)
0049 {
0050 assert(_file0);
0051
0052 TH2 *hLayerCellCharge = (TH2 *) _file0->GetObjectChecked("hLayerCellCharge", "TH2");
0053 assert(hLayerCellCharge);
0054
0055 TH1 *chargePDFs[100000];
0056 TH1 *chargeCDFs[100000];
0057
0058 Color_t colors[] = {kPink + 2, kRed + 2, kSpring + 2, kGreen + 2, kAzure + 2, kBlue + 2};
0059
0060 for (int i = 0; i < layers.size(); ++i)
0061 {
0062 Color_t color = colors[i];
0063
0064 const int layer = layers[i];
0065 TH1 *chargePDF = hLayerCellCharge->ProjectionY(Form("hLayerCellCharge_PDF_Layer%d", layer), layer, layer);
0066 chargePDF->Sumw2();
0067
0068 chargePDF->SetMarkerStyle(kFullCircle);
0069 chargePDF->SetMarkerColor(color);
0070 chargePDF->SetLineColor(color);
0071
0072 chargePDF->Scale(1 / chargePDF->Integral(0, -1));
0073
0074 TH1 *chargeCDF = chargePDF->Clone(Form("hLayerCellCharge_CDF_Layer%d", layer));
0075 for (int bin = chargeCDF->GetNbinsX() + 1; bin >= 0; --bin)
0076 {
0077 double cdf = 0;
0078 double cdf_err = 0;
0079 cdf = chargePDF->IntegralAndError(bin, -1, cdf_err);
0080
0081 chargeCDF->SetBinContent(bin, cdf);
0082 chargeCDF->SetBinError(bin, cdf_err);
0083 }
0084
0085 chargePDF->Rebin(10);
0086
0087 chargePDFs[i] = (chargePDF);
0088 chargeCDFs[i] = chargeCDF;
0089 }
0090
0091 TCanvas *c1 = new TCanvas("ChargePerLayer", "ChargePerLayer", 1800, 960);
0092 c1->Divide(2, 1);
0093 int idx = 1;
0094 TPad *p;
0095
0096 p = (TPad *) c1->cd(idx++);
0097 c1->Update();
0098 p->SetLogy();
0099
0100 p->DrawFrame(0, 1e-4, 610, 1, ";Charge [fC];Probability/bin");
0101
0102 TLegend *leg = new TLegend(.17, .8, .93, .93);
0103 leg->AddEntry("", "#it{#bf{sPHENIX}} Simulation, CD-1 configuration", "");
0104 leg->AddEntry("", description, "");
0105 leg->AddEntry("", "TPC charge PDF per FEE chan. over 13 us drift window", "");
0106 leg->Draw();
0107
0108 leg = new TLegend(.55, .5, .93, .75);
0109 for (int i = 0; i < layers.size(); ++i)
0110 {
0111 chargePDFs[i]->Draw("same");
0112
0113 TString performance = Form("Layer #%d: <Q> = %.0f fC", layers[i], chargePDFs[i]->GetMean());
0114 cout << "ChargePerLayer : " << performance << endl;
0115 leg->AddEntry(chargePDFs[i], performance, "lp");
0116 }
0117
0118 leg->Draw();
0119
0120 p = (TPad *) c1->cd(idx++);
0121 c1->Update();
0122 p->SetLogy();
0123 p->DrawFrame(0, 1e-3, 610, 1, ";Charge Threshold [fC];Probability[Q > Charge Threshold]");
0124
0125 TLegend *leg = new TLegend(.17, .8, .93, .93);
0126 leg->AddEntry("", "#it{#bf{sPHENIX}} Simulation, CD-1 configuration", "");
0127 leg->AddEntry("", description, "");
0128 leg->AddEntry("", "TPC charge CCDF per FEE chan. over 13 us drift window", "");
0129 leg->Draw();
0130
0131 leg = new TLegend(.2, .2, .65, .4);
0132 for (int i = 0; i < layers.size(); ++i)
0133 {
0134 chargeCDFs[i]->Draw("same");
0135
0136 TString cdf_desk = Form("Layer #%d: P[Q>300fC] = %.1f%%",
0137 layers[i],
0138 100 * chargeCDFs[i]->GetBinContent(chargeCDFs[i]->FindBin(300)));
0139 cout << "ChargePerLayer : " << cdf_desk << endl;
0140 leg->AddEntry(chargeCDFs[i], cdf_desk, "lp");
0141 }
0142 leg->Draw();
0143
0144 SaveCanvas(c1, TString(_file0->GetName()) + TString(c1->GetName()), kTRUE);
0145 }
0146
0147 void Check()
0148 {
0149 assert(_file0);
0150
0151 TCanvas *c1 = new TCanvas("Check", "Check", 1800, 960);
0152 c1->Divide(3, 2);
0153 int idx = 1;
0154 TPad *p;
0155
0156 p = (TPad *) c1->cd(idx++);
0157 c1->Update();
0158 p->SetLogy();
0159
0160 TH1 *hNormalization = (TH1 *) _file0->GetObjectChecked("hNormalization", "TH1");
0161 assert(hNormalization);
0162
0163 hNormalization->Draw();
0164
0165 p = (TPad *) c1->cd(idx++);
0166 c1->Update();
0167 p->SetLogz();
0168
0169 TH2 *hLayerCellHit = (TH2 *) _file0->GetObjectChecked("hLayerCellHit", "TH2");
0170 assert(hLayerCellHit);
0171
0172 hLayerCellHit->Draw("colz");
0173
0174 TProfile *hLayerCellHit_px = hLayerCellHit->ProfileX();
0175 hLayerCellHit_px->Draw("same");
0176
0177 p = (TPad *) c1->cd(idx++);
0178 c1->Update();
0179 p->SetLogz();
0180
0181 TH2 *hLayerCellCharge = (TH2 *) _file0->GetObjectChecked("hLayerCellCharge", "TH2");
0182 assert(hLayerCellCharge);
0183
0184 hLayerCellCharge->Draw("colz");
0185
0186 p = (TPad *) c1->cd(idx++);
0187 c1->Update();
0188 p->SetLogz();
0189
0190 TH2 *hLayerSumCellHit = (TH2 *) _file0->GetObjectChecked("hLayerSumCellHit", "TH2");
0191 assert(hLayerSumCellHit);
0192
0193 hLayerSumCellHit->Draw("colz");
0194
0195 p = (TPad *) c1->cd(idx++);
0196 c1->Update();
0197 p->SetLogz();
0198
0199 TH2 *hLayerSumCellCharge = (TH2 *) _file0->GetObjectChecked("hLayerSumCellCharge", "TH2");
0200 assert(hLayerSumCellCharge);
0201
0202 hLayerSumCellCharge->Draw("colz");
0203
0204 SaveCanvas(c1, TString(_file0->GetName()) + TString(c1->GetName()), kTRUE);
0205 }