In [1]:
%%cpp -d
// pair_timewalk_correction.C

#include <iostream>
#include <vector>
#include <cmath>

#include "TFile.h"
#include "TTree.h"
#include "TH1F.h"
#include "TH2F.h"
#include "TProfile.h"
#include "TF1.h"
#include "TCanvas.h"
#include "TROOT.h"
#include "TStyle.h"
#include "TLatex.h"
#include "TFitResultPtr.h"

const int MAXHIT = 1024;
const int NCLUSTER = 12;
const int NSEG = 7;

double E_HIGH = 600.0;
double E_REF = 600.0;
double E_MAX_OFFSET = 3000.0;

TH2F *h_tw_id[NSEG] = {0};
TProfile *hp_tw_id[NSEG] = {0};

TH1F *h_off_pair[NSEG][NSEG] = {0};
TF1  *f_off_pair[NSEG][NSEG] = {0};

double t_off[NCLUSTER][NSEG][NSEG];
int has_toff[NCLUSTER][NSEG][NSEG];

TF1 *f_tw[NCLUSTER][NSEG] = {0};
int has_fit[NCLUSTER][NSEG];

Long64_t seg_high_count[NCLUSTER][NSEG];

double seg_offset[NCLUSTER][NSEG];
int has_seg_offset[NCLUSTER][NSEG];

int ref_seg_cluster[NCLUSTER];


// ------------------------------------------------------------
// Make pair-combination time-walk histograms for one cluster
// ------------------------------------------------------------
std::vector<TH2F*> make_pair_timewalk_cluster(int clusterID = 6)
{
  TH1::AddDirectory(kFALSE);

  std::vector<TH2F*> hlist;

  if (clusterID < 0 || clusterID >= NCLUSTER) {
    std::cout << "invalid clusterID" << std::endl;
    return hlist;
  }

  TFile *fin = new TFile("eurica_event.root");
  if (!fin || fin->IsZombie()) {
    std::cout << "cannot open eurica_event.root" << std::endl;
    return hlist;
  }

  TTree *tree = (TTree*)fin->Get("tree");
  if (!tree) {
    std::cout << "cannot find tree" << std::endl;
    fin->Close();
    return hlist;
  }

  int ghit;
  int gid[MAXHIT];
  double ge[MAXHIT];
  double gt[MAXHIT];

  tree->SetBranchAddress("ghit", &ghit);
  tree->SetBranchAddress("gid", gid);
  tree->SetBranchAddress("ge", ge);
  tree->SetBranchAddress("gt", gt);

  for (int s = 0; s < NSEG; s++) {
    seg_high_count[clusterID][s] = 0;
  }

  // reset offset histograms and fit functions
  for (int i = 0; i < NSEG; i++) {
    for (int j = 0; j < NSEG; j++) {

      has_toff[clusterID][i][j] = 0;
      t_off[clusterID][i][j] = 0.0;

      if (h_off_pair[i][j]) {
        delete h_off_pair[i][j];
        h_off_pair[i][j] = 0;
      }

      if (f_off_pair[i][j]) {
        delete f_off_pair[i][j];
        f_off_pair[i][j] = 0;
      }

      if (i == j) continue;

      h_off_pair[i][j] =
        new TH1F(Form("h_off_c%d_%d_%d", clusterID, i, j),
                 Form("cluster %d: pair %d-%d; t_{%d}-t_{%d} (ns); counts",
                      clusterID, i, j, i, j),
                 120, -600, 600);

      h_off_pair[i][j]->SetDirectory(0);
    }
  }

  Long64_t nentries = tree->GetEntries();

  // first pass: fill high-energy pair offset histograms
  for (Long64_t ientry = 0; ientry < nentries; ientry++) {

    tree->GetEntry(ientry);

    for (int a = 0; a < ghit; a++) {

      int cid_a = gid[a] / 7;
      int sid_a = gid[a] % 7;

      if (cid_a != clusterID) continue;
      if (sid_a < 0 || sid_a >= NSEG) continue;
      if (ge[a] < E_HIGH || ge[a] > E_MAX_OFFSET) continue;

      seg_high_count[clusterID][sid_a]++;

      for (int b = 0; b < ghit; b++) {

        if (a == b) continue;

        int cid_b = gid[b] / 7;
        int sid_b = gid[b] % 7;

        if (cid_b != clusterID) continue;
        if (sid_b < 0 || sid_b >= NSEG) continue;
        if (sid_b == sid_a) continue;
        if (ge[b] < E_HIGH || ge[b] > E_MAX_OFFSET) continue;

        double dt = gt[a] - gt[b];

        h_off_pair[sid_a][sid_b]->Fill(dt);
      }
    }
  }

  // fit pair offsets
  for (int i = 0; i < NSEG; i++) {
    for (int j = 0; j < NSEG; j++) {

      if (i == j) continue;
      if (!h_off_pair[i][j]) continue;

      if (h_off_pair[i][j]->GetEntries() < 50) {
        continue;
      }

      int maxbin = h_off_pair[i][j]->GetMaximumBin();
      double x0 = h_off_pair[i][j]->GetBinCenter(maxbin);

      f_off_pair[i][j] =
        new TF1(Form("f_off_c%d_%d_%d", clusterID, i, j),
                "gaus",
                x0 - 120.0,
                x0 + 120.0);

      f_off_pair[i][j]->SetParameter(0, h_off_pair[i][j]->GetMaximum());
      f_off_pair[i][j]->SetParameter(1, x0);
      f_off_pair[i][j]->SetParameter(2, 60.0);

      TFitResultPtr r = h_off_pair[i][j]->Fit(f_off_pair[i][j], "RQ0");

      if ((int)r == 0) {
        t_off[clusterID][i][j] = f_off_pair[i][j]->GetParameter(1);
      } else {
        t_off[clusterID][i][j] = x0;
      }

      has_toff[clusterID][i][j] = 1;
    }
  }

  // reset time-walk histograms
  for (int s = 0; s < NSEG; s++) {

    if (h_tw_id[s]) {
      delete h_tw_id[s];
      h_tw_id[s] = 0;
    }

    h_tw_id[s] =
      new TH2F(Form("h_tw_c%d_id%d", clusterID, s),
               Form("cluster %d: segment %d, all references; #Delta t_{ij} (ns); E_{seg %d}",
                    clusterID, s, s),
               80, -600, 600,
               130, 0, 1300);

    h_tw_id[s]->SetDirectory(0);
    hlist.push_back(h_tw_id[s]);
  }

  // second pass: fill pair-offset-corrected time-walk plots
  for (Long64_t ientry = 0; ientry < nentries; ientry++) {

    tree->GetEntry(ientry);

    for (int a = 0; a < ghit; a++) {

      int cid_a = gid[a] / 7;
      int sid_a = gid[a] % 7;

      if (cid_a != clusterID) continue;
      if (sid_a < 0 || sid_a >= NSEG) continue;
      if (ge[a] < 30) continue;

      for (int b = 0; b < ghit; b++) {

        if (a == b) continue;

        int cid_b = gid[b] / 7;
        int sid_b = gid[b] % 7;

        if (cid_b != clusterID) continue;
        if (sid_b < 0 || sid_b >= NSEG) continue;
        if (sid_b == sid_a) continue;
        if (ge[b] < E_HIGH) continue;

        if (!has_toff[clusterID][sid_a][sid_b]) continue;

        double dt = gt[a] - gt[b];
        double dt_corr = dt - t_off[clusterID][sid_a][sid_b];

        h_tw_id[sid_a]->Fill(dt_corr, ge[a]);
      }
    }
  }

  fin->Close();

  std::cout << "cluster " << clusterID
            << " pair time-walk histograms filled." << std::endl;

  for (int s = 0; s < NSEG; s++) {
    std::cout << "segment " << s
              << " high-energy count = " << seg_high_count[clusterID][s]
              << ", walk entries = " << h_tw_id[s]->GetEntries()
              << std::endl;
  }

  return hlist;
}


// ------------------------------------------------------------
// Draw 7 segment time-walk plots for current cluster
// ------------------------------------------------------------
void draw_pair_timewalk_cluster(int clusterID = 6)
{
  TCanvas *c_old = (TCanvas*)gROOT->FindObject("c_pair_tw");
  if (c_old) delete c_old;

  TCanvas *c = new TCanvas("c_pair_tw",
                           Form("cluster %d: pair-combination time walk", clusterID),
                           900, 600);

  c->Divide(4, 2);
  gStyle->SetOptStat(0);

  for (int s = 0; s < NSEG; s++) {
    c->cd(s + 1);
    if (h_tw_id[s]) h_tw_id[s]->Draw("colz");
  }

  c->Draw();
}


// ------------------------------------------------------------
// Draw pair offset histograms for current cluster
// ------------------------------------------------------------
void draw_pair_offset_cluster(int clusterID = 6)
{
  TCanvas *c_old = (TCanvas*)gROOT->FindObject("c_pair_off");
  if (c_old) delete c_old;

  TCanvas *c = new TCanvas("c_pair_off",
                           Form("cluster %d: pair offsets", clusterID),
                           1050, 1050);

  c->Divide(7, 7);
  gStyle->SetOptStat(0);

  for (int i = 0; i < NSEG; i++) {
    for (int j = 0; j < NSEG; j++) {

      int ipad = i * NSEG + j + 1;
      c->cd(ipad);

      if (i == j) {
        TLatex text;
        text.SetTextAlign(22);
        text.SetTextSize(0.18);
        text.DrawLatexNDC(0.5, 0.5, Form("%d = %d", i, j));
        continue;
      }

      if (!h_off_pair[i][j]) continue;

      h_off_pair[i][j]->SetLineColor(kBlack);
      h_off_pair[i][j]->Draw();

      if (f_off_pair[i][j]) {
        f_off_pair[i][j]->SetLineColor(kRed);
        f_off_pair[i][j]->Draw("same");
      }

      TLatex label;
      label.SetTextSize(0.10);
      label.SetNDC();

      if (has_toff[clusterID][i][j]) {
        label.DrawLatex(0.18, 0.82,
                        Form("%d-%d: %.1f ns", i, j,
                             t_off[clusterID][i][j]));
      } else {
        label.DrawLatex(0.18, 0.82,
                        Form("%d-%d: no fit", i, j));
      }
    }
  }

  c->Draw();
}


// ------------------------------------------------------------
// Choose reference segment by highest high-energy count
// ------------------------------------------------------------
int ChooseReferenceSegmentByHighCount(int clusterID)
{
  int bestSeg = -1;
  Long64_t bestCount = -1;

  for (int s = 0; s < NSEG; s++) {

    if (!has_fit[clusterID][s]) continue;

    Long64_t count = seg_high_count[clusterID][s];

    std::cout << "cluster " << clusterID
              << ", reference candidate segment " << s
              << ": high-energy count = "
              << count << std::endl;

    if (count > bestCount) {
      bestCount = count;
      bestSeg = s;
    }
  }

  return bestSeg;
}


// ------------------------------------------------------------
// Absorb segment offset into p0
// ------------------------------------------------------------
void absorb_segment_offset_to_p0(int clusterID)
{
  if (clusterID < 0 || clusterID >= NCLUSTER) return;

  int refSeg = ChooseReferenceSegmentByHighCount(clusterID);
  ref_seg_cluster[clusterID] = refSeg;

  if (refSeg < 0) {
    std::cout << "cluster " << clusterID
              << ": no valid reference segment found." << std::endl;
    return;
  }

  std::cout << "cluster " << clusterID
            << ": selected reference segment = "
            << refSeg << std::endl;

  for (int s = 0; s < NSEG; s++) {

    seg_offset[clusterID][s] = 0.0;
    has_seg_offset[clusterID][s] = 0;

    if (!has_fit[clusterID][s]) continue;
    if (!f_tw[clusterID][s]) continue;

    double C = 0.0;

    if (s == refSeg) {
      C = 0.0;
      has_seg_offset[clusterID][s] = 1;
    } else if (has_toff[clusterID][s][refSeg]) {
      C = t_off[clusterID][s][refSeg];
      has_seg_offset[clusterID][s] = 1;
    } else {
      std::cout << "cluster " << clusterID
                << ", segment " << s
                << ": no direct offset to reference segment "
                << refSeg
                << ", set C = 0" << std::endl;
      C = 0.0;
    }

    seg_offset[clusterID][s] = C;

    // Shift p0 so that f(E_REF) = C.
    double fref = f_tw[clusterID][s]->Eval(E_REF);
    double p0_old = f_tw[clusterID][s]->GetParameter(0);
    double p0_new = p0_old - fref + C;

    f_tw[clusterID][s]->SetParameter(0, p0_new);

    std::cout << "cluster " << clusterID
              << ", segment " << s
              << ": C = " << C
              << " ns, p0 -> " << p0_new
              << std::endl;
  }
}


// ------------------------------------------------------------
// Fit time-walk functions for current cluster
// ------------------------------------------------------------
void fit_pair_timewalk_cluster(int clusterID = 6, bool draw = true)
{
  if (clusterID < 0 || clusterID >= NCLUSTER) {
    std::cout << "invalid clusterID" << std::endl;
    return;
  }

  double fitEmin = 30.0;
  double fitEmax = 800.0;

  double dtMin = -400.0;
  double dtMax = 400.0;

  TCanvas *c = 0;

  if (draw) {
    TCanvas *c_old = (TCanvas*)gROOT->FindObject("c_fit_tw");
    if (c_old) delete c_old;

    c = new TCanvas("c_fit_tw",
                    Form("cluster %d time-walk fits", clusterID),
                    900, 600);

    c->Divide(4, 2);
    c->SetLogy(0);
  }

  for (int s = 0; s < NSEG; s++) {

    has_fit[clusterID][s] = 0;

    if (!h_tw_id[s]) continue;
    if (h_tw_id[s]->GetEntries() < 100) continue;

    int xbin1 = h_tw_id[s]->GetXaxis()->FindBin(dtMin);
    int xbin2 = h_tw_id[s]->GetXaxis()->FindBin(dtMax);

    if (hp_tw_id[s]) {
      delete hp_tw_id[s];
      hp_tw_id[s] = 0;
    }

    hp_tw_id[s] = h_tw_id[s]->ProfileY(Form("hp_tw_c%d_s%d", clusterID, s),
                                       xbin1, xbin2);

    hp_tw_id[s]->SetDirectory(0);

    if (f_tw[clusterID][s]) {
      delete f_tw[clusterID][s];
      f_tw[clusterID][s] = 0;
    }

    f_tw[clusterID][s] =
      new TF1(Form("f_tw_c%d_s%d", clusterID, s),
              "[0]+[1]/sqrt(x)+[2]/x+[3]/(x*x)",
              30.0, 1300.0);

    f_tw[clusterID][s]->SetParameters(0.0, 1000.0, -10000.0, 100000.0);

    hp_tw_id[s]->Fit(f_tw[clusterID][s], "R0Q", "", fitEmin, fitEmax);

    has_fit[clusterID][s] = 1;

    std::cout << "cluster " << clusterID
              << ", segment " << s
              << " fit:"
              << " p0=" << f_tw[clusterID][s]->GetParameter(0)
              << ", p1=" << f_tw[clusterID][s]->GetParameter(1)
              << ", p2=" << f_tw[clusterID][s]->GetParameter(2)
              << ", p3=" << f_tw[clusterID][s]->GetParameter(3)
              << std::endl;

    if (draw) {
      c->cd(s + 1);
      gPad->SetLogy(0);

      hp_tw_id[s]->SetMarkerStyle(20);
      hp_tw_id[s]->SetMarkerSize(0.7);
      hp_tw_id[s]->SetMinimum(-300);
      hp_tw_id[s]->SetMaximum(300);

      hp_tw_id[s]->Draw();

      TF1 *f_draw = (TF1*)f_tw[clusterID][s]->Clone(
        Form("f_draw_c%d_s%d", clusterID, s)
      );

      f_draw->SetLineColor(kRed);
      f_draw->SetRange(30, 1300);
      f_draw->Draw("same");
    }
  }

  // After fitting walk shape, absorb high-energy segment offset into p0.
  absorb_segment_offset_to_p0(clusterID);

  if (draw) c->Draw();
}


// ------------------------------------------------------------
// Fit all clusters
// ------------------------------------------------------------
void fit_all_clusters_pair()
{
  for (int c = 0; c < NCLUSTER; c++) {

    std::cout << "======================================" << std::endl;
    std::cout << "Processing cluster " << c << std::endl;

    make_pair_timewalk_cluster(c);
    fit_pair_timewalk_cluster(c, false);
  }

  std::cout << "All cluster time-walk fits finished." << std::endl;
}


// ------------------------------------------------------------
// Evaluate final correction.
// p0 already includes high-energy offset.
// ------------------------------------------------------------
double GetTimeWalkCorrection(int clusterID, int segID, double E)
{
  if (clusterID < 0 || clusterID >= NCLUSTER) return 0.0;
  if (segID < 0 || segID >= NSEG) return 0.0;
  if (!has_fit[clusterID][segID]) return 0.0;
  if (!f_tw[clusterID][segID]) return 0.0;
  if (E <= 0) return 0.0;

  return f_tw[clusterID][segID]->Eval(E);
}


// ------------------------------------------------------------
// Generate corrected ROOT file
// ------------------------------------------------------------
void make_time_corrected_tree_pair(
    const char *inputFile = "eurica_event.root",
    const char *outputFile = "eurica_time_pair.root")
{
  TFile *fin = new TFile(inputFile);
  if (!fin || fin->IsZombie()) {
    std::cout << "cannot open " << inputFile << std::endl;
    return;
  }

  TTree *tin = (TTree*)fin->Get("tree");
  if (!tin) {
    std::cout << "cannot find tree in " << inputFile << std::endl;
    fin->Close();
    return;
  }

  int ghit;
  int gid_in[MAXHIT];
  double ge_in[MAXHIT];
  double gt_in[MAXHIT];

  tin->SetBranchAddress("ghit", &ghit);
  tin->SetBranchAddress("gid", gid_in);
  tin->SetBranchAddress("ge", ge_in);
  tin->SetBranchAddress("gt", gt_in);

  TFile *fout = new TFile(outputFile, "recreate");
  TTree *tout = new TTree("tree", "pair time-walk corrected tree");

  int o_ghit;
  int o_gid[MAXHIT];
  double o_ge[MAXHIT];
  double o_gt[MAXHIT];

  tout->Branch("ghit", &o_ghit, "ghit/I");
  tout->Branch("gid", o_gid, "gid[ghit]/I");
  tout->Branch("ge", o_ge, "ge[ghit]/D");
  tout->Branch("gt", o_gt, "gt[ghit]/D");

  Long64_t nentries = tin->GetEntries();

  for (Long64_t ientry = 0; ientry < nentries; ientry++) {

    tin->GetEntry(ientry);

    o_ghit = ghit;

    for (int i = 0; i < ghit; i++) {

      int clusterID = gid_in[i] / 7;
      int segID = gid_in[i] % 7;

      double corr = GetTimeWalkCorrection(clusterID, segID, ge_in[i]);

      o_gid[i] = gid_in[i];
      o_ge[i] = ge_in[i];
      o_gt[i] = gt_in[i] - corr;
    }

    tout->Fill();
  }

  fout->cd();
  tout->Write();
  fout->Close();
  fin->Close();

  std::cout << "time-corrected file saved to "
            << outputFile << std::endl;
}

同一Cluster内不同id之间的toff¶

  • 可以看出几何上相邻单元之间才有明显的关联。
In [3]:
make_pair_timewalk_cluster(6);
draw_pair_offset_cluster(6);
cluster 6 pair time-walk histograms filled.
segment 0 high-energy count = 38109, walk entries = 7236
segment 1 high-energy count = 36768, walk entries = 7083
segment 2 high-energy count = 37226, walk entries = 7120
segment 3 high-energy count = 39352, walk entries = 7345
segment 4 high-energy count = 38709, walk entries = 7518
segment 5 high-energy count = 36949, walk entries = 7314
segment 6 high-energy count = 42415, walk entries = 15350

合并后同一Cluster内不同id之的timewalk¶

In [5]:
draw_pair_timewalk_cluster(6);

timewalk拟合¶

In [7]:
fit_pair_timewalk_cluster(6, true);
cluster 6, segment 0 fit: p0=39.4166, p1=-646.399, p2=-11067.2, p3=207189
cluster 6, segment 1 fit: p0=51.1366, p1=-988.767, p2=-7951.39, p3=152485
cluster 6, segment 2 fit: p0=38.9093, p1=-560.6, p2=-12613.7, p3=268957
cluster 6, segment 3 fit: p0=48.5325, p1=-1039.83, p2=-7714.07, p3=132969
cluster 6, segment 4 fit: p0=53.0691, p1=-1076.93, p2=-6734, p3=129831
cluster 6, segment 5 fit: p0=59.4729, p1=-1265.45, p2=-5295.2, p3=130088
cluster 6, segment 6 fit: p0=32.7498, p1=-385.362, p2=-13004.3, p3=239117
cluster 6, reference candidate segment 0: high-energy count = 38109
cluster 6, reference candidate segment 1: high-energy count = 36768
cluster 6, reference candidate segment 2: high-energy count = 37226
cluster 6, reference candidate segment 3: high-energy count = 39352
cluster 6, reference candidate segment 4: high-energy count = 38709
cluster 6, reference candidate segment 5: high-energy count = 36949
cluster 6, reference candidate segment 6: high-energy count = 42415
cluster 6: selected reference segment = 6
cluster 6, segment 0: C = -1.1966 ns, p0 -> 43.0623
cluster 6, segment 1: C = 2.05581 ns, p0 -> 55.2508
cluster 6, segment 2: C = -3.75621 ns, p0 -> 39.4059
cluster 6, segment 3: C = 0.458135 ns, p0 -> 55.3964
cluster 6, segment 4: C = 2.94166 ns, p0 -> 57.7697
cluster 6, segment 5: C = -2.86643 ns, p0 -> 57.2594
cluster 6, segment 6: C = 0 ns, p0 -> 36.742

将timewalk修正写入root文件¶

In [9]:
fit_all_clusters_pair();
make_time_corrected_tree_pair("eurica_event.root", "eurica_time_pair.root");
======================================
Processing cluster 0
cluster 0 pair time-walk histograms filled.
segment 0 high-energy count = 23477, walk entries = 4070
segment 1 high-energy count = 21592, walk entries = 4094
segment 2 high-energy count = 24359, walk entries = 4604
segment 3 high-energy count = 27349, walk entries = 5046
segment 4 high-energy count = 27472, walk entries = 4820
segment 5 high-energy count = 23838, walk entries = 4330
segment 6 high-energy count = 25520, walk entries = 9093
cluster 0, segment 0 fit: p0=46.9143, p1=-908.998, p2=-9770.05, p3=199874
cluster 0, segment 1 fit: p0=39.1479, p1=-545.049, p2=-11566.5, p3=220851
cluster 0, segment 2 fit: p0=18.8904, p1=5.97335, p2=-16258.6, p3=285243
cluster 0, segment 3 fit: p0=17.6447, p1=174.833, p2=-17614.4, p3=319376
cluster 0, segment 4 fit: p0=31.1526, p1=-209.942, p2=-16754.1, p3=328047
cluster 0, segment 5 fit: p0=61.401, p1=-1571.08, p2=-814.486, p3=29498.1
cluster 0, segment 6 fit: p0=48.6658, p1=-1018.22, p2=-6539.74, p3=111728
cluster 0, reference candidate segment 0: high-energy count = 23477
cluster 0, reference candidate segment 1: high-energy count = 21592
cluster 0, reference candidate segment 2: high-energy count = 24359
cluster 0, reference candidate segment 3: high-energy count = 27349
cluster 0, reference candidate segment 4: high-energy count = 27472
cluster 0, reference candidate segment 5: high-energy count = 23838
cluster 0, reference candidate segment 6: high-energy count = 25520
cluster 0: selected reference segment = 4
cluster 0, segment 0: no direct offset to reference segment 4, set C = 0
cluster 0, segment 0: C = 0 ns, p0 -> 52.8379
cluster 0, segment 1: no direct offset to reference segment 4, set C = 0
cluster 0, segment 1: C = 0 ns, p0 -> 40.9155
cluster 0, segment 2: no direct offset to reference segment 4, set C = 0
cluster 0, segment 2: C = 0 ns, p0 -> 26.0615
cluster 0, segment 3: C = 4.39136 ns, p0 -> 25.724
cluster 0, segment 4: C = 0 ns, p0 -> 35.5832
cluster 0, segment 5: C = 2.29096 ns, p0 -> 67.7054
cluster 0, segment 6: C = 2.8805 ns, p0 -> 55.0382
======================================
Processing cluster 1
cluster 1 pair time-walk histograms filled.
segment 0 high-energy count = 19302, walk entries = 2376
segment 1 high-energy count = 18293, walk entries = 2226
segment 2 high-energy count = 18380, walk entries = 2370
segment 3 high-energy count = 20309, walk entries = 2501
segment 4 high-energy count = 20688, walk entries = 2342
segment 5 high-energy count = 19840, walk entries = 2239
segment 6 high-energy count = 0, walk entries = 0
cluster 1, segment 0 fit: p0=2.72733, p1=453.987, p2=-19531.8, p3=342211
cluster 1, segment 1 fit: p0=42.4434, p1=-662.267, p2=-9465.33, p3=174666
cluster 1, segment 2 fit: p0=22.0485, p1=-391.402, p2=-10882.6, p3=186893
cluster 1, segment 3 fit: p0=29.613, p1=-302.946, p2=-15302.7, p3=310241
cluster 1, segment 4 fit: p0=18.2845, p1=240.706, p2=-16612.5, p3=265068
cluster 1, segment 5 fit: p0=18.4379, p1=-297.994, p2=-12455.5, p3=222050
cluster 1, reference candidate segment 0: high-energy count = 19302
cluster 1, reference candidate segment 1: high-energy count = 18293
cluster 1, reference candidate segment 2: high-energy count = 18380
cluster 1, reference candidate segment 3: high-energy count = 20309
cluster 1, reference candidate segment 4: high-energy count = 20688
cluster 1, reference candidate segment 5: high-energy count = 19840
cluster 1: selected reference segment = 4
cluster 1, segment 0: no direct offset to reference segment 4, set C = 0
cluster 1, segment 0: C = 0 ns, p0 -> 13.0684
cluster 1, segment 1: no direct offset to reference segment 4, set C = 0
cluster 1, segment 1: C = 0 ns, p0 -> 42.3273
cluster 1, segment 2: no direct offset to reference segment 4, set C = 0
cluster 1, segment 2: C = 0 ns, p0 -> 33.5974
cluster 1, segment 3: C = -13.6848 ns, p0 -> 23.3257
cluster 1, segment 4: C = 0 ns, p0 -> 17.1244
cluster 1, segment 5: C = 12.0502 ns, p0 -> 44.3581
======================================
Processing cluster 2
cluster 2 pair time-walk histograms filled.
segment 0 high-energy count = 26164, walk entries = 4582
segment 1 high-energy count = 23175, walk entries = 4498
segment 2 high-energy count = 22341, walk entries = 4284
segment 3 high-energy count = 23482, walk entries = 4199
segment 4 high-energy count = 22835, walk entries = 3944
segment 5 high-energy count = 24352, walk entries = 4295
segment 6 high-energy count = 25164, walk entries = 8763
cluster 2, segment 0 fit: p0=61.7783, p1=-1539.6, p2=-1014.77, p3=56923
cluster 2, segment 1 fit: p0=47.0527, p1=-1110.43, p2=-5379.93, p3=136598
cluster 2, segment 2 fit: p0=55.0697, p1=-1247.85, p2=-4614.96, p3=109363
cluster 2, segment 3 fit: p0=51.9785, p1=-904.966, p2=-9566.43, p3=215939
cluster 2, segment 4 fit: p0=53.2853, p1=-1158.24, p2=-5787.81, p3=149887
cluster 2, segment 5 fit: p0=51.0463, p1=-1226.16, p2=-4157.79, p3=118586
cluster 2, segment 6 fit: p0=38.9778, p1=-650.119, p2=-11682.5, p3=261651
cluster 2, reference candidate segment 0: high-energy count = 26164
cluster 2, reference candidate segment 1: high-energy count = 23175
cluster 2, reference candidate segment 2: high-energy count = 22341
cluster 2, reference candidate segment 3: high-energy count = 23482
cluster 2, reference candidate segment 4: high-energy count = 22835
cluster 2, reference candidate segment 5: high-energy count = 24352
cluster 2, reference candidate segment 6: high-energy count = 25164
cluster 2: selected reference segment = 0
cluster 2, segment 0: C = 0 ns, p0 -> 64.3872
cluster 2, segment 1: C = 17.8343 ns, p0 -> 71.7544
cluster 2, segment 2: no direct offset to reference segment 0, set C = 0
cluster 2, segment 2: C = 0 ns, p0 -> 58.3312
cluster 2, segment 3: no direct offset to reference segment 0, set C = 0
cluster 2, segment 3: C = 0 ns, p0 -> 52.2893
cluster 2, segment 4: no direct offset to reference segment 0, set C = 0
cluster 2, segment 4: C = 0 ns, p0 -> 56.5149
cluster 2, segment 5: C = 18.9946 ns, p0 -> 75.6526
cluster 2, segment 6: C = 18.1755 ns, p0 -> 63.4604
======================================
Processing cluster 3
cluster 3 pair time-walk histograms filled.
segment 0 high-energy count = 0, walk entries = 0
segment 1 high-energy count = 28798, walk entries = 3661
segment 2 high-energy count = 28299, walk entries = 5180
segment 3 high-energy count = 25926, walk entries = 4693
segment 4 high-energy count = 24756, walk entries = 4466
segment 5 high-energy count = 25560, walk entries = 3217
segment 6 high-energy count = 27795, walk entries = 8733
cluster 3, segment 1 fit: p0=31.2095, p1=-462.621, p2=-11266.1, p3=200300
cluster 3, segment 2 fit: p0=18.7066, p1=280.982, p2=-21087.9, p3=495165
cluster 3, segment 3 fit: p0=44.9472, p1=-884.959, p2=-8289.06, p3=164535
cluster 3, segment 4 fit: p0=72.8989, p1=-1881.61, p2=-1001.38, p3=96289.6
cluster 3, segment 5 fit: p0=32.9805, p1=-355.419, p2=-14102.2, p3=282203
cluster 3, segment 6 fit: p0=28.6165, p1=-166.268, p2=-15706.7, p3=284249
cluster 3, reference candidate segment 1: high-energy count = 28798
cluster 3, reference candidate segment 2: high-energy count = 28299
cluster 3, reference candidate segment 3: high-energy count = 25926
cluster 3, reference candidate segment 4: high-energy count = 24756
cluster 3, reference candidate segment 5: high-energy count = 25560
cluster 3, reference candidate segment 6: high-energy count = 27795
cluster 3: selected reference segment = 1
cluster 3, segment 1: C = 0 ns, p0 -> 37.1069
cluster 3, segment 2: C = -2.99823 ns, p0 -> 19.3018
cluster 3, segment 3: no direct offset to reference segment 1, set C = 0
cluster 3, segment 3: C = 0 ns, p0 -> 49.4863
cluster 3, segment 4: no direct offset to reference segment 1, set C = 0
cluster 3, segment 4: C = 0 ns, p0 -> 78.218
cluster 3, segment 5: no direct offset to reference segment 1, set C = 0
cluster 3, segment 5: C = 0 ns, p0 -> 37.2296
cluster 3, segment 6: C = -7.61034 ns, p0 -> 24.5658
======================================
Processing cluster 4
cluster 4 pair time-walk histograms filled.
segment 0 high-energy count = 22972, walk entries = 4309
segment 1 high-energy count = 25337, walk entries = 4507
segment 2 high-energy count = 26014, walk entries = 4747
segment 3 high-energy count = 24935, walk entries = 4591
segment 4 high-energy count = 23178, walk entries = 4211
segment 5 high-energy count = 23043, walk entries = 4206
segment 6 high-energy count = 25167, walk entries = 9251
cluster 4, segment 0 fit: p0=43.6768, p1=-830.266, p2=-7831.74, p3=176549
cluster 4, segment 1 fit: p0=46.9207, p1=-970.319, p2=-6882.34, p3=140910
cluster 4, segment 2 fit: p0=31.5695, p1=-501.841, p2=-10980.3, p3=219712
cluster 4, segment 3 fit: p0=39.3321, p1=-843.919, p2=-10875.2, p3=247024
cluster 4, segment 4 fit: p0=67.1792, p1=-1453.32, p2=-4484.63, p3=123739
cluster 4, segment 5 fit: p0=41.8103, p1=-764.65, p2=-9550.24, p3=205897
cluster 4, segment 6 fit: p0=56.8719, p1=-1252.92, p2=-5575.73, p3=147492
cluster 4, reference candidate segment 0: high-energy count = 22972
cluster 4, reference candidate segment 1: high-energy count = 25337
cluster 4, reference candidate segment 2: high-energy count = 26014
cluster 4, reference candidate segment 3: high-energy count = 24935
cluster 4, reference candidate segment 4: high-energy count = 23178
cluster 4, reference candidate segment 5: high-energy count = 23043
cluster 4, reference candidate segment 6: high-energy count = 25167
cluster 4: selected reference segment = 2
cluster 4, segment 0: no direct offset to reference segment 2, set C = 0
cluster 4, segment 0: C = 0 ns, p0 -> 46.458
cluster 4, segment 1: C = -4.44724 ns, p0 -> 46.245
cluster 4, segment 2: C = 0 ns, p0 -> 38.1778
cluster 4, segment 3: C = 3.8488 ns, p0 -> 55.7407
cluster 4, segment 4: no direct offset to reference segment 2, set C = 0
cluster 4, segment 4: C = 0 ns, p0 -> 66.4623
cluster 4, segment 5: no direct offset to reference segment 2, set C = 0
cluster 4, segment 5: C = 0 ns, p0 -> 46.5618
cluster 4, segment 6: C = -3.43383 ns, p0 -> 56.5995
======================================
Processing cluster 5
cluster 5 pair time-walk histograms filled.
segment 0 high-energy count = 24976, walk entries = 4635
segment 1 high-energy count = 25802, walk entries = 4869
segment 2 high-energy count = 27221, walk entries = 5216
segment 3 high-energy count = 28442, walk entries = 5246
segment 4 high-energy count = 27764, walk entries = 5090
segment 5 high-energy count = 25846, walk entries = 4649
segment 6 high-energy count = 27581, walk entries = 10079
cluster 5, segment 0 fit: p0=34.6657, p1=-530.603, p2=-10413.6, p3=166687
cluster 5, segment 1 fit: p0=41.0514, p1=-795.632, p2=-6779.22, p3=112648
cluster 5, segment 2 fit: p0=43.1535, p1=-609.049, p2=-12555.1, p3=239242
cluster 5, segment 3 fit: p0=37.1823, p1=-455.975, p2=-13598.6, p3=256779
cluster 5, segment 4 fit: p0=48.8997, p1=-940.716, p2=-8884.55, p3=195052
cluster 5, segment 5 fit: p0=41.3897, p1=-544.429, p2=-13345.8, p3=288446
cluster 5, segment 6 fit: p0=32.5894, p1=-443.573, p2=-14100.3, p3=281080
cluster 5, reference candidate segment 0: high-energy count = 24976
cluster 5, reference candidate segment 1: high-energy count = 25802
cluster 5, reference candidate segment 2: high-energy count = 27221
cluster 5, reference candidate segment 3: high-energy count = 28442
cluster 5, reference candidate segment 4: high-energy count = 27764
cluster 5, reference candidate segment 5: high-energy count = 25846
cluster 5, reference candidate segment 6: high-energy count = 27581
cluster 5: selected reference segment = 3
cluster 5, segment 0: no direct offset to reference segment 3, set C = 0
cluster 5, segment 0: C = 0 ns, p0 -> 38.5548
cluster 5, segment 1: no direct offset to reference segment 3, set C = 0
cluster 5, segment 1: C = 0 ns, p0 -> 43.4673
cluster 5, segment 2: C = 1.05977 ns, p0 -> 46.1847
cluster 5, segment 3: C = 0 ns, p0 -> 40.5662
cluster 5, segment 4: C = 5.77046 ns, p0 -> 58.4408
cluster 5, segment 5: no direct offset to reference segment 3, set C = 0
cluster 5, segment 5: C = 0 ns, p0 -> 43.6679
cluster 5, segment 6: C = -3.26797 ns, p0 -> 37.5605
======================================
Processing cluster 6
cluster 6 pair time-walk histograms filled.
segment 0 high-energy count = 38109, walk entries = 7236
segment 1 high-energy count = 36768, walk entries = 7083
segment 2 high-energy count = 37226, walk entries = 7120
segment 3 high-energy count = 39352, walk entries = 7345
segment 4 high-energy count = 38709, walk entries = 7518
segment 5 high-energy count = 36949, walk entries = 7314
segment 6 high-energy count = 42415, walk entries = 15350
cluster 6, segment 0 fit: p0=39.4166, p1=-646.399, p2=-11067.2, p3=207189
cluster 6, segment 1 fit: p0=51.1366, p1=-988.767, p2=-7951.39, p3=152485
cluster 6, segment 2 fit: p0=38.9093, p1=-560.6, p2=-12613.7, p3=268957
cluster 6, segment 3 fit: p0=48.5325, p1=-1039.83, p2=-7714.07, p3=132969
cluster 6, segment 4 fit: p0=53.0691, p1=-1076.93, p2=-6734, p3=129831
cluster 6, segment 5 fit: p0=59.4729, p1=-1265.45, p2=-5295.2, p3=130088
cluster 6, segment 6 fit: p0=32.7498, p1=-385.362, p2=-13004.3, p3=239117
cluster 6, reference candidate segment 0: high-energy count = 38109
cluster 6, reference candidate segment 1: high-energy count = 36768
cluster 6, reference candidate segment 2: high-energy count = 37226
cluster 6, reference candidate segment 3: high-energy count = 39352
cluster 6, reference candidate segment 4: high-energy count = 38709
cluster 6, reference candidate segment 5: high-energy count = 36949
cluster 6, reference candidate segment 6: high-energy count = 42415
cluster 6: selected reference segment = 6
cluster 6, segment 0: C = -1.1966 ns, p0 -> 43.0623
cluster 6, segment 1: C = 2.05581 ns, p0 -> 55.2508
cluster 6, segment 2: C = -3.75621 ns, p0 -> 39.4059
cluster 6, segment 3: C = 0.458135 ns, p0 -> 55.3964
cluster 6, segment 4: C = 2.94166 ns, p0 -> 57.7697
cluster 6, segment 5: C = -2.86643 ns, p0 -> 57.2594
cluster 6, segment 6: C = 0 ns, p0 -> 36.742
======================================
Processing cluster 7
cluster 7 pair time-walk histograms filled.
segment 0 high-energy count = 0, walk entries = 0
segment 1 high-energy count = 0, walk entries = 0
segment 2 high-energy count = 0, walk entries = 0
segment 3 high-energy count = 0, walk entries = 0
segment 4 high-energy count = 0, walk entries = 0
segment 5 high-energy count = 0, walk entries = 0
segment 6 high-energy count = 0, walk entries = 0
cluster 7: no valid reference segment found.
======================================
Processing cluster 8
cluster 8 pair time-walk histograms filled.
segment 0 high-energy count = 30849, walk entries = 4008
segment 1 high-energy count = 31579, walk entries = 0
segment 2 high-energy count = 0, walk entries = 0
segment 3 high-energy count = 31161, walk entries = 4405
segment 4 high-energy count = 0, walk entries = 0
segment 5 high-energy count = 30951, walk entries = 3948
segment 6 high-energy count = 0, walk entries = 0
cluster 8, segment 0 fit: p0=42.1219, p1=-770.11, p2=-9801.26, p3=202896
cluster 8, segment 3 fit: p0=39.4478, p1=-757.884, p2=-8054.64, p3=133820
cluster 8, segment 5 fit: p0=10.986, p1=445.772, p2=-18840.4, p3=321026
cluster 8, reference candidate segment 0: high-energy count = 30849
cluster 8, reference candidate segment 3: high-energy count = 31161
cluster 8, reference candidate segment 5: high-energy count = 30951
cluster 8: selected reference segment = 3
cluster 8, segment 0: C = 2.29613 ns, p0 -> 49.5076
cluster 8, segment 3: C = 0 ns, p0 -> 43.9931
cluster 8, segment 5: C = 0.291655 ns, p0 -> 12.602
======================================
Processing cluster 9
cluster 9 pair time-walk histograms filled.
segment 0 high-energy count = 31634, walk entries = 5709
segment 1 high-energy count = 30677, walk entries = 3895
segment 2 high-energy count = 0, walk entries = 0
segment 3 high-energy count = 29651, walk entries = 3678
segment 4 high-energy count = 30298, walk entries = 5439
segment 5 high-energy count = 29608, walk entries = 5649
segment 6 high-energy count = 30863, walk entries = 9685
cluster 9, segment 0 fit: p0=29.3227, p1=-301.69, p2=-12189.3, p3=222431
cluster 9, segment 1 fit: p0=37.7305, p1=-727.206, p2=-7802.84, p3=120859
cluster 9, segment 3 fit: p0=33.5071, p1=-693.217, p2=-9943.4, p3=187050
cluster 9, segment 4 fit: p0=46.413, p1=-826.257, p2=-9506.35, p3=174143
cluster 9, segment 5 fit: p0=23.4127, p1=-104.019, p2=-16549.8, p3=314715
cluster 9, segment 6 fit: p0=36.8393, p1=-561.554, p2=-10264.3, p3=173887
cluster 9, reference candidate segment 0: high-energy count = 31634
cluster 9, reference candidate segment 1: high-energy count = 30677
cluster 9, reference candidate segment 3: high-energy count = 29651
cluster 9, reference candidate segment 4: high-energy count = 30298
cluster 9, reference candidate segment 5: high-energy count = 29608
cluster 9, reference candidate segment 6: high-energy count = 30863
cluster 9: selected reference segment = 0
cluster 9, segment 0: C = 0 ns, p0 -> 32.0141
cluster 9, segment 1: C = -0.623385 ns, p0 -> 41.7337
cluster 9, segment 3: no direct offset to reference segment 0, set C = 0
cluster 9, segment 3: C = 0 ns, p0 -> 44.3532
cluster 9, segment 4: no direct offset to reference segment 0, set C = 0
cluster 9, segment 4: C = 0 ns, p0 -> 49.092
cluster 9, segment 5: C = -2.88797 ns, p0 -> 28.0674
cluster 9, segment 6: C = -2.25608 ns, p0 -> 37.2934
======================================
Processing cluster 10
cluster 10 pair time-walk histograms filled.
segment 0 high-energy count = 30611, walk entries = 5755
segment 1 high-energy count = 29777, walk entries = 3666
segment 2 high-energy count = 0, walk entries = 0
segment 3 high-energy count = 31417, walk entries = 3700
segment 4 high-energy count = 31165, walk entries = 5762
segment 5 high-energy count = 31149, walk entries = 5866
segment 6 high-energy count = 32097, walk entries = 9725
cluster 10, segment 0 fit: p0=34.8898, p1=-572.249, p2=-10549.5, p3=200767
cluster 10, segment 1 fit: p0=32.7966, p1=-496.914, p2=-10679.2, p3=189795
cluster 10, segment 3 fit: p0=63.3226, p1=-1571.29, p2=-3129.45, p3=72410.1
cluster 10, segment 4 fit: p0=36.0404, p1=-416.213, p2=-11889.7, p3=224314
cluster 10, segment 5 fit: p0=52.7488, p1=-1281.27, p2=-4887.8, p3=102030
cluster 10, segment 6 fit: p0=10.8916, p1=409.656, p2=-19952.9, p3=373266
cluster 10, reference candidate segment 0: high-energy count = 30611
cluster 10, reference candidate segment 1: high-energy count = 29777
cluster 10, reference candidate segment 3: high-energy count = 31417
cluster 10, reference candidate segment 4: high-energy count = 31165
cluster 10, reference candidate segment 5: high-energy count = 31149
cluster 10, reference candidate segment 6: high-energy count = 32097
cluster 10: selected reference segment = 6
cluster 10, segment 0: C = -2.90838 ns, p0 -> 37.4783
cluster 10, segment 1: C = 2.44311 ns, p0 -> 40.001
cluster 10, segment 3: C = 10.334 ns, p0 -> 79.4962
cluster 10, segment 4: C = -0.418961 ns, p0 -> 35.7659
cluster 10, segment 5: C = 1.83962 ns, p0 -> 62.0102
cluster 10, segment 6: C = 0 ns, p0 -> 15.4938
======================================
Processing cluster 11
cluster 11 pair time-walk histograms filled.
segment 0 high-energy count = 29802, walk entries = 5518
segment 1 high-energy count = 29940, walk entries = 5403
segment 2 high-energy count = 30507, walk entries = 5653
segment 3 high-energy count = 31249, walk entries = 5749
segment 4 high-energy count = 31492, walk entries = 5889
segment 5 high-energy count = 30202, walk entries = 5574
segment 6 high-energy count = 30510, walk entries = 11467
cluster 11, segment 0 fit: p0=25.6912, p1=-161.685, p2=-15973, p3=293904
cluster 11, segment 1 fit: p0=40.3897, p1=-625.603, p2=-13500.6, p3=285462
cluster 11, segment 2 fit: p0=37.4116, p1=-324.749, p2=-16863.1, p3=342162
cluster 11, segment 3 fit: p0=40.5205, p1=-476.462, p2=-16280.4, p3=332150
cluster 11, segment 4 fit: p0=36.8603, p1=-276.541, p2=-16777.6, p3=337295
cluster 11, segment 5 fit: p0=32.9735, p1=-405.425, p2=-13416.4, p3=248888
cluster 11, segment 6 fit: p0=22.7251, p1=108.401, p2=-18156.7, p3=321132
cluster 11, reference candidate segment 0: high-energy count = 29802
cluster 11, reference candidate segment 1: high-energy count = 29940
cluster 11, reference candidate segment 2: high-energy count = 30507
cluster 11, reference candidate segment 3: high-energy count = 31249
cluster 11, reference candidate segment 4: high-energy count = 31492
cluster 11, reference candidate segment 5: high-energy count = 30202
cluster 11, reference candidate segment 6: high-energy count = 30510
cluster 11: selected reference segment = 4
cluster 11, segment 0: no direct offset to reference segment 4, set C = 0
cluster 11, segment 0: C = 0 ns, p0 -> 32.4059
cluster 11, segment 1: no direct offset to reference segment 4, set C = 0
cluster 11, segment 1: C = 0 ns, p0 -> 47.2482
cluster 11, segment 2: no direct offset to reference segment 4, set C = 0
cluster 11, segment 2: C = 0 ns, p0 -> 40.4124
cluster 11, segment 3: C = -4.42229 ns, p0 -> 41.2406
cluster 11, segment 4: C = 0 ns, p0 -> 38.3155
cluster 11, segment 5: C = 7.89635 ns, p0 -> 46.1171
cluster 11, segment 6: C = 2.27839 ns, p0 -> 27.2221
All cluster time-walk fits finished.
time-corrected file saved to eurica_time_pair.root
In [ ]: