基于望远镜法的带电粒子鉴别¶

实验设置:¶

  • 束流:能量为 25MeVA 的 $^{20}O$ 束流轰击 $^9Be$ 靶,在靶上产生碎裂反应。

  • 探测器:反应产物由三片厚度为 $1000 \mu m (\Delta E_1)$, $500 \mu m (\Delta E_2)$ 和 $1000 \mu m (\Delta E_3)$的 Si 探测器组成的望远镜探测。每个Si探测器为位置灵敏探测器,可以给出入射粒子的x,y的位置。

所有粒子不穿出 $\Delta E_3$探测器

In [1]:
//%jsroot on
TCanvas *c1 = new TCanvas;
TFile *f = new TFile("telenew.root");
tree->SetAlias("de1","d1.e[0]");
tree->SetAlias("de2","d2.e[0]");
tree->SetAlias("de3","d3.e[0]");
tree->SetAlias("eall","de1+de2+de3");
TCut cut = "d1.hit==2 && ssd1<1";
Warning in <TClass::Init>: no dictionary for class ppac is available
Warning in <TClass::Init>: no dictionary for class dssd is available

入射粒子在$\Delta E_1$探测器上的x-y 位置分布¶

In [2]:
tree->Draw("d1.x:d1.y>>hxy(32,0,32,32,0,32)",cut,"colz");
hxy->SetTitle("x ~ y; y;x");
c1->Draw();

$\Delta E_1$ ~ $\Delta E_2$ 二维关联图¶

In [3]:
c1->Clear();
tree->Draw("de1:de2>>hd1d2(500,0,250,1000,0,500)",cut,"colz");
hd1d2->SetTitle("#Delta E_{1} ~ #Delta E_{2}; #Delta E_{2}; #Delta E_{1}");
gPad->SetLogz();
c1->Draw();

$\Delta E_1$ ~ $\Delta E_2$¶

  • 蓝色:入射粒子未穿出 $\Delta E_2$, 即 $\Delta E_3=0$, 随着 $\Delta E_2$ 的增加 $\Delta E_1$ 减小。

  • 红色:入射粒子穿出 $\Delta E_2$, 即 $\Delta E_3>0$,随着 $\Delta E_2$ 的增加 $\Delta E_1$ 增加。

  • 每个同位素的红色和蓝色的交界处,称为穿透点

In [4]:
c1->Clear();
tree->Draw("de1:de2>>h2a(500,0,250,1000,0,500)",cut && "de3<1","colz");
tree->Draw("de1:de2>>h2b(500,0,250,1000,0,500)",cut && "de3>2","colz");
auto h2a = (TH2*) gROOT->FindObject("h2a");
auto h2b = (TH2*) gROOT->FindObject("h2b");
h2a->SetMarkerColor(kBlue);
h2b->SetMarkerColor(kRed);
h2a->SetTitle("#Delta E_{1} ~ #Delta E_{2} ; #Delta E_{2}; #Delta E_{1}");
h2a->Draw();
h2b->Draw("same");
c1->Draw();

$\Delta E_1$ ~ $\Delta E_2$¶

  • 左图:入射粒子未穿出 $\Delta E_2$,即$\Delta E_3=0$
  • 右图: 入射粒子穿出 $\Delta E_2$, 即$\Delta E_3>0$
  • 粒子未穿出第二片探测器: 可进行同位素鉴别
  • 粒子穿出第二片探测器: 无法进行同位素鉴别
In [5]:
tree->Draw("de1:de2>>hd1d2_1(500,0,250,1000,0,500)",cut && "de3<1","colz");
tree->Draw("de1:de2>>hd1d2_2(500,0,250,1000,0,500)",cut && "de3>2","colz");
TCanvas *c2=new TCanvas("c2","c2",1200,600);
c2->Clear();
c2->SetLogz();
c2->Divide(2,1);
c2->cd(1);
hd1d2_1->SetTitle("#Delta E_{1} ~ #Delta E_{2} : punch through #Delta E_{2}   ; #Delta E_{2}; #Delta E_{1}");
hd1d2_1->Draw("colz");
c2->cd(2);
hd1d2_2->SetTitle("#Delta E_{1} ~ #Delta E_{2} : stop at #Delta E_{2}   ; #Delta E_{2}; #Delta E_{1}");
hd1d2_2->Draw("colz");
gPad->SetLogz();
c2->Draw();

$\Delta E_2$ ~ $\Delta E_3$ 二维关联图¶

In [6]:
c1->Clear();
tree->Draw("de2:de3>>hd2d3(500,0,250,500,0,250)",cut,"colz");
hd2d3->SetTitle("#Delta E_{2} ~ #Delta E_{3} ; #Delta E_{3}; #Delta E_{2}");
gPad->SetLogz();
c1->Draw();

粒子鉴别- $E_f=Δ𝐸·𝐸$¶

给定原子核在一个极薄硅层中的能量损失与粒子能量成反比,即 Δ𝐸·𝐸 = 常数。如果探测器的厚度足够小,则对于相同原子核的不同能量, $𝐸_𝑓$=Δ𝐸·𝐸 的平方根是常数。

  • 下图显示上述结论并不成立,其原因是$\Delta E_1, \Delta E_2$的探测器不能看成极薄的探测器。
In [7]:
c2->Clear();
c2->Clear();
c2->Divide(2,1);
c2->cd(1);
tree->Draw("sqrt(de1*de2):de2 >>hpid2da(300,0,300,300,0,300)",cut&&"de3<1","colz");
hpid2da->SetTitle("E_{f} ~ E : #Delta E_{1}-#Delta E_{2}; E_{f}; E");
gPad->SetLogz();

c2->cd(2);
tree->Draw("sqrt(de2*de3):de3 >>hpid2db(300,0,300,200,0,200)",cut,"colz");
hpid2db->SetTitle("E_{f} ~ E : #Delta E_{2}-#Delta E_{3}; E_{f}; E");
gPad->SetLogz();
c2->Draw();

粒子鉴别- $E_f = \sqrt{ΔE · E + αΔE^2} + βE$¶

上一个表达式的主要偏差来源于第一片探测器的厚度不可忽略。

如果从 E 到 ΔE+E 的区间对Δ𝐸·𝐸值积分,即在探测器的能量区间上进行积分,可得到改进的表达式:

$$ E_f^2 = \int_E^{E+\Delta E}{E dE} = \frac{{(E+\Delta E)}^2}{2} - \frac{E^2}{2} = ΔE · E + 0.5ΔE^2 $$

经验表明,0.7 比通过积分得到的 0.5 系数更为合适,并且存在一个$E$的轻微线性依赖关系,通过引入线性因子 βE(β = 0.05)来进行补偿。

$$ E_f = \sqrt{ΔE · E + αΔE^2} + βE $$

无量纲常数的初值选为 α = 0.7 和 β = 0.05,通过调节(拟合)取值, 使不同能量的相同原子核的 $E_f$ 值尽可能恒定。式中$E,ΔE 和 E_f$单位为MeV。

  • [Eur. Phys. J. A (2015) 51: 93]
In [8]:
c2->Clear();
c2->Clear();
c2->Divide(2,1);
c2->cd(1);
tree->Draw("sqrt(de1*de2+0.7*de1*de1)+0.05*de2:de2 >>hpid2dc(300,0,300,400,0,400)",cut&&"de3<1","colz");
hpid2dc->SetTitle("E_{f} ~ E : #Delta E_{1}-#Delta E_{2}; E; E_{f}");
gPad->SetLogz();

c2->cd(2);
tree->Draw("sqrt(de2*de3+0.62*de2*de2)-0.01*de3:de3  >>hpid2dd(300,0,300,200,0,200)",cut&&"de3>2","colz");
hpid2dd->SetTitle("E_{f} ~ E : #Delta E_{2}-#Delta E_{3}; E; E_{f}");
gPad->SetLogz();
c2->Draw();

粒子鉴别-$E_f$投影图¶

  • 下图中每个高斯峰,代表一种同位素。
In [9]:
c2->Clear();
c2->Divide(2,1);
c2->cd(1);
tree->Draw("sqrt(de1*de2+0.7*de1*de1)+0.05*de2 >>hpid1da(1000,0,400)",cut&&"de3<1","colz");
hpid1da->SetTitle("E_{f} : #Delta E_{1}-#Delta E_{2}");
gPad->SetLogy();

c2->cd(2);
tree->Draw("sqrt(de2*de3+0.6*de2*de2)-0.01*de3 >>hpid1db(500,0,200)",cut&&"de3>2","colz");
hpid1db->SetTitle("E_{f} : #Delta E_{2}-#Delta E_{3}");
gPad->SetLogy();
c2->Draw();
In [ ]: