DSSD front-back Correlation(I) - dssd1¶
Front-Back Energy Correlation¶
When a single charged particle hits the intersection of front strip $X_i$ and back strip $Y_j$ in a DSSD, the charges collected on the two sides originate from the same deposited energy. For an ideal single-pixel event,
$$ E_{\mathrm{true}} = E_{x,i} = E_{y,j}. $$
Therefore, the front-side and back-side signals should show a one-to-one energy correlation.
In the DAQ system, however, the directly recorded quantities are the raw ADC amplitudes $A_{x,i}$ and $A_{y,j}$, rather than the physical energies themselves. Because the electronic response differs from strip to strip, the raw amplitudes on the two sides are correlated, but not yet directly comparable on a common scale.
The purpose of the present tutorial is to use this front-back correlation of purified single-pixel events to establish a detector-internal relative energy normalization.

Relative Energy Calibration¶
To compare signals from different strips, we introduce a unified relative energy scale inside the detector.
Choose one front strip $X_r$ as the baseline reference and define
$$ E^{(\mathrm{rel})} \equiv A_{x,r}, \qquad k_{x,r}=1, \qquad b_{x,r}=0. $$
Then each strip is mapped onto this common relative scale by a linear transformation:
$$ E_{x,i}^{(\mathrm{rel})} = k_{x,i} A_{x,i} + b_{x,i}, \qquad E_{y,j}^{(\mathrm{rel})} = k_{y,j} A_{y,j} + b_{y,j}. $$
Here $A_{x,i}$ and $A_{y,j}$ are the raw amplitudes measured on strips $X_i$ and $Y_j$, while $(k,b)$ are the strip-wise calibration parameters.
For a correctly paired single-pixel event, the normalized front-side and back-side energies should satisfy
$$ E_{x,i}^{(\mathrm{rel})} \approx E_{y,j}^{(\mathrm{rel})}. $$
The task of relative calibration is therefore to determine the parameters $(k,b)$ so that this front-back consistency is achieved on purified correlation data.
Multi-hit position matching¶
The relative energy normalization established here is the basis for later front-back pairing in multi-hit events. Once all strip signals are mapped onto a common scale, the consistency between normalized front-side and back-side energies can be used as an additional constraint for resolving pairing ambiguity.
Assume two particles strike the detector at the same time, producing true hit points at $(X_1, Y_1)$ and $(X_2, Y_2)$ with distinct energies ($E_1 \neq E_2$). The readout system records signals from front strips $X_1, X_2$ and back strips $Y_1, Y_2$. Relying solely on geometric topology, it is impossible to distinguish whether the genuine physical combination is $\{(X_1, Y_1), (X_2, Y_2)\}$ or the spurious "ghost" combination $\{(X_1, Y_2), (X_2, Y_1)\}$.
By normalizing all strip signals to the unified scale $E^{(\mathrm{rel})}$, we can construct a residual criterion based on the physical principle of front-back energy consistency:
$$ \Delta E = |E^{(\mathrm{rel})}_{x} - E^{(\mathrm{rel})}_{y}| $$
Physical Matching Criterion: For a correctly paired front-back strip combination, their normalized energy difference $\Delta E$ should be minimized (approaching zero). Therefore, relative normalization serves as an essential prerequisite for precise multi-particle identification and position reconstruction.

In the present tutorial, however, we focus only on purified single-pixel events and the extraction of the relative calibration coefficients.
Pixel-Pixel Calibration Method¶
In this tutorial, we consider the simplest calibration strategy: pixel-pixel correlation.
Fix one front strip $X_{i^*}$ as the reference. For each back strip $Y_j$, select events belonging to the pixel $(i^*,j)$ and fit
$$ A_{x,i^*} = k_{y,j} A_{y,j} + b_{y,j}. $$
This determines the relative calibration of the back strip $Y_j$:
$$ E_{y,j}^{(\mathrm{rel})} = k_{y,j} A_{y,j} + b_{y,j}. $$
Then choose one calibrated back strip $Y_{j^*}$ as the new reference and fit each front strip $X_k$,
$$ E_{y,j^*}^{(\mathrm{rel})} = k_{x,k} A_{x,k} + b_{x,k}, $$
so that
$$ E_{x,k}^{(\mathrm{rel})} = k_{x,k} A_{x,k} + b_{x,k}. $$
This method directly uses the front-back correlation at the level of individual pixels. It is conceptually straightforward, but its quality depends on whether the selected pixels contain sufficient statistics.
In the present tutorial, we restrict ourselves to this pixel-pixel method. A more global strip-strip normalization strategy will be introduced in a later tutorial.
Detector Setup¶
ROOT File: data_16C.root
- Reaction System: 25 MeV/A $^{16}$C + $^9$Be, measuring large-angle multiple charged particles produced from the fragmentation of $^{16}$C.
- Detector Array: A telescope array consisting of 3 DSSDs + CsI, positioned at zero degrees behind the target.
- D1, D2, and D3 are all double-sided silicon strip detectors with 32 (X) $\times$ 32 (Y) strips.
- Strip width: 2 mm, inter-strip spacing: 0.1 mm.
- Hardware Trigger: Requires at least two hits on the X-side of both D1 and D2 (Multiplicity $\ge$ 2).
- Data Preprocessing: Pedestal events have been removed.

Branch:
// xenergy, yenergy, xtime (Array of 32 for each detector)
d1x[32], d1y[32], d1t[32];
d2x[32], d2y[32], d2t[32];
d3x[32], d3y[32], d3t[32];
// hit multiplicity, energy, time, strip ID (Compact arrays based on hits)
d1xhit, d1xe[d1xhit], d1xt[d1xhit], d1xs[d1xhit];
d2xhit, d2xe[d2xhit], d2xt[d2xhit], d2xs[d2xhit];
d3xhit, d3xe[d3xhit], d3xt[d3xhit], d3xs[d3xhit];
d1yhit, d1ye[d1yhit], d1yt[d1yhit], d1ys[d1yhit];
d2yhit, d2ye[d2yhit], d2yt[d2yhit], d2ys[d2yhit];
d3yhit, d3ye[d3yhit], d3yt[d3yhit], d3ys[d3yhit];
TFile *ipf = new TFile("./data/data_16C.root");
TTree *tree = (TTree*)ipf->Get("tree");
TCanvas *c1 = new TCanvas("c1","c1");
Multiplicity for X and Y Sides
The hardware trigger condition ($M \ge 2$) ensures that the vast majority of events have a front and back multiplicity of 2 or higher. Events with energy below the pedestal threshold have been filtered out.
tree->Draw("d1xhit: d1yhit>>(15, 0, 15, 15, 0, 15)", "", "colz");
gPad->SetLogz();
c1->Draw();
- Most of these $M \ge 2$ events are Charge Sharing events between adjacent strips, while a small fraction consists of genuine multi-particle incident events.
tree->Scan("d1xe: d1xs: d1ye: d1ys", "", "", 10, 1);
*********************************************************************** * Row * Instance * d1xe * d1xs * d1ye * d1ys * *********************************************************************** * 1 * 0 * 4800 * 17 * 5775 * 18 * * 1 * 1 * 982 * 18 * * * * 2 * 0 * 3069 * 23 * 3352 * 11 * * 2 * 1 * 341 * 24 * 111 * 10 * * 3 * 0 * 5325 * 23 * 4522 * 15 * * 3 * 1 * * * 922 * 16 * * 4 * 0 * 5822 * 21 * 5688 * 15 * * 4 * 1 * * * 271 * 14 * * 5 * 0 * 4048 * 11 * 5213 * 12 * * 5 * 1 * 1149 * 12 * * * * 6 * 0 * 3719 * 14 * 4255 * 16 * * 6 * 1 * 374 * 13 * * * * 7 * 0 * 4546 * 12 * 4961 * 9 * * 7 * 1 * 297 * 11 * * * * 8 * 0 * 4051 * 22 * 3877 * 9 * * 8 * 1 * * * 329 * 10 * * 9 * 0 * 6720 * 13 * 6513 * 12 * * 9 * 1 * 142 * 10 * 418 * 13 * * 9 * 2 * 109 * 14 * * * * 10 * 0 * 3873 * 20 * 3007 * 8 * * 10 * 1 * * * 974 * 9 * ***********************************************************************
Front-Back Correlation
A clear linear correlation band is observed in the central region. However, prominent scattered points outside this band result from Charge Sharing—where the deposited energy splits between adjacent strips, causing the energy recorded in a single strip to be less than the total deposited energy.
gPad->SetLogy(0);
tree->Draw("d1x[12]: d1y[13]>>(1000, 0, 8000, 1000, 0, 8000)", "", "colz");
c1->Draw();//front-back correlation
X-Y Correlations Unaffected by Adjacent Strip Energy Sharing
To mitigate the impact of energy sharing between adjacent strips, two approaches are traditionally considered:
1. Global Single Multiplicity Requirement: Require exactly one hit on both X and Y sides (xhit==1 && yhit==1). This effectively eliminates scattered points but results in severely low event statistics due to our specific hardware trigger condition ($M \ge 2$).
c1->Clear();
TCut chit1 = "d1xs==12 && d1ys==13 && d1xhit==1 && d1yhit==1";
tree->Draw("d1xe: d1ye>>(1000, 0, 8000, 1000, 0, 8000)", chit1, "colz");
c1->Draw();
2. Local Isolation Cut (Anti-coincidence Condition): Instead of restricting global multiplicity, strictly require that the selected pixel has no signal in its immediate spatially adjacent strips.
- Use a threshold of
< 50channels (representing the hardware pedestal/noise limit). Ensuring neighboring strips are below this threshold topologically guarantees no charge sharing. - For
d1x[12]: Strips 11 and 13 must have no signal (< 50). - For
d1y[13]: Strips 12 and 14 must have no signal (< 50).
Applying this local isolation condition significantly recovers usable event statistics while maintaining high purity of the linear response.
c1->Clear();
TCut cveto = "d1x[11]<50 && d1x[13]<50 && d1y[12]<50 && d1y[14]<50"; // no sharing
TCut c1213 = "d1x[12]>200 && d1y[13]>200 && d1y[13]<8000" && cveto;
tree->Draw("d1x[12]: d1y[13]>>h2(1000, 0, 8000, 1000, 0, 8000)", c1213, "colz");
c1->Draw();
1. Extracting Unbinned Data for Regression
To extract the exact mathematical correlation, we must perform an unbinned regression. We bypass the binned TH2 histogram and directly map the purified TTree events into a TGraph.
TGraph *gr = new TGraph(tree->GetSelectedRows(), tree->GetV2(), tree->GetV1());
gr->SetMarkerSize(0.2);
gr->Draw("A*");
c1->Draw();
2. Standard vs. Robust Regression
Standard Least Squares (pol1): Severely biased by residual, uncleaned cross-talk outliers. These sparse extreme points act as "levers," pulling the fit away from the true diagonal.
Robust Fitting ("ROB"): ROOT's robust fitting utilizes Least Trimmed Squares to automatically down-weight extreme outliers not originating from the main data-generating process.
ROBUST REGRESSION¶
In the presence of outliers that do not come from the same data-generating process as the rest of the data, least squares estimation is inefficient and can be biased. Robust regression methods are designed to be not overly affected by violations of assumptions by the underlying data-generating process. ROOT supports robust fitting for polynomial functions using the "ROB" option.
https://root.cern/doc/master/fitLinearRobust_8C.html

ROOT only supports robust fitting for polynomial functions.
// 1. Standard Fit (Least Squares)
TF1 *fStd = new TF1("fStd", "pol1", 200, 8000);
fStd->SetLineColor(kRed);
gr->Fit(fStd, "RQ");
// 2. Robust Fit (Least Trimmed Squares)
TF1 *fRob = new TF1("fRob", "pol1", 200, 8000);
fRob->SetLineColor(kBlue);
fRob->SetLineWidth(1);
// The "+" option adds the new function to the graph's list, preserving the previous one
gr->Fit(fRob, "R+ ROB Q");
// Draw the graph with points
gr->Draw("A*");
// 3. Add Legend for visual comparison
// (Positioned at top-left to avoid overlapping the y=x diagonal data)
TLegend *leg = new TLegend(0.15, 0.75, 0.45, 0.88);
leg->SetBorderSize(0);
leg->AddEntry(fStd, "Standard Fit (Pulled by outliers)", "l");
leg->AddEntry(fRob, "Robust Fit (Resists outliers)", "l");
leg->Draw();
c1->Draw();
3. Residual Distribution & The "Density Trap"
Let's examine the robust fit's residual distribution: $Y - (p_1 \cdot X + p_0)$.
Even with robust fitting, a slight curvature at high energies remains. This is not a physical nonlinearity of the amplifier. It is a statistical artifact caused by the Density Weight: the overwhelmingly high statistics of low-energy fragments deeply anchor the cost function, slightly compromising the mathematically correct high-energy slope.
// Extract the parameters from the Robust Fit
double p0 = fRob->GetParameter(0);
double p1 = fRob->GetParameter(1);
// Construct the draw command for the residual: Y vs (X - (p1*Y + p0))
TString stree;
stree.Form("d1y[13]:d1x[12]-(%f*d1y[13]+%f)>>ha(100,-40,40,1000,0,8000)", p1, p0);
// Ensure the pad is cleared and ready for a 2D color plot
c1->Clear();
tree->Draw(stree.Data(), c1213, "colz");
c1->Draw();
Assignment¶
Save the highly purified single-pixel events from DSSD 1/2/3—strictly applying the local isolation cut to exclude adjacent charge sharing—into structured d1/2/3xy.root files. TThese files will serve as the data foundation for the subsequent front-back relative calibration studies.
1. Global 2D Correlation Filter (cut1/2/3.C)
Create a rough geometric boundary to enclose the main $Y \approx X$ diagonal band. Draw a 2D scatter plot under ideal clean conditions (xhit==1 && yhit==1), use the ROOT TBrowser Graphical Cut to draw a polygon, and save it as a macro (e.g., cut1.C). This topologically rejects massive far-field noises.
2. Event Extraction Code Structure
// Load the global 2D graphical filter
gROOT->Macro("./data/cut1.C");
TCutG *cut1 = (TCutG*)gROOT->GetListOfSpecials()->FindObject("cut1");
// Helper function: Ensure strip is isolated (neighbors < pedestal threshold)
bool IsStripIsolated(Int_t* energyArray, int stripId) {
if (stripId > 0 && energyArray[stripId - 1] > 50) return false;
if (stripId < 31 && energyArray[stripId + 1] > 50) return false;
return true;
}
// ... [Initialize TFile, TTree, branches] ...
TFile *fout = new TFile("d1xy.root", "RECREATE");
TTree *tout = new TTree("tree", "Purified XY Correlation");
Int_t ix, iy, xe, ye;
// ... [Branch bindings for ix, iy, xe, ye] ...
for(Long64_t jentry=0; jentry<nentries; jentry++) {
tree->GetEntry(jentry);
for(int i=0; i<32; i++) {
if(d1x[i] < 200 || !IsStripIsolated(d1x, i)) continue; // X criteria
for(int j=0; j<32; j++) {
if(d1y[j] < 200 || !IsStripIsolated(d1y, j)) continue; // Y criteria
if(!cut1->IsInside(d1x[i], d1y[j])) continue; // Global 2D limit
ix = i; iy = j; xe = d1x[i]; ye = d1y[j];
tout->Fill();
}
}
}
fout->Write(); fout->Close();
Expected results for DSSD1¶
TCanvas *c1 = new TCanvas("c1","c1");
TFile *fin = new TFile("./data/d1xy.root");
TTree *tree =(TTree *)fin->Get("tree");
tree->Draw("ye:xe>>(4000,0,8000,4000,0,8000)","","colz");
c1->SetLogz();
c1->Draw();
Warning in <TCanvas::Constructor>: Deleting canvas with same name: c1