TkN 2.1
Toolkit for Nuclei

Explore nuclear deformation island

This example produces four different figures illustrating, with four different physical observations, the sudden offset of deformation at N=60 in the A~100 region:

  • Two neutron separation energy (MeV)
  • charge radius (fm)
  • Reduced electric transition probability of the first E2 transition. (Weisskopf per mass unit)
  • Energy of the first 2+ state (MeM)

Start tkn-root, compile and execute this macro :

tkn-root
_____ _ _ | Documentation: https://tkn.in2p3.fr/
(_ _) | \ | | | Source: https://gitlab.in2p3.fr/tkn/tkn-lib
| |_ _| \| | |
| | |/ / | | Version 1.0
| | <| |\ | |
|_|_|\_\_| \_| | Database: TkN_ensdf_221101_xundl_210701_v1.0.db
tkn [0] .L deformation_island.C+
tkn [1] draw_deformation_island()

It produces the following plot:

Source code :

#include "TMultiGraph.h"
#include "TGraph.h"
#include "tkmanager.h"
#include "TROOT.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TColor.h"
#include "TLegend.h"
#include "TLegendEntry.h"
#include "TLine.h"
using namespace tkn;
const char* ccs[7] = {"#1ddbdc","#1aabbf","#1a83a5","#1d6190","#1e427a","#182562","#120952"};
int cols[7];
int marks[7] = {20,33,21,47,22,34,23};
double masize[7] = {1.3,1.8,1,1.4,1.4,1.4,1.4};
TGraph* get_graph(int zz);
TLegend* get_legend();
TLine* get_line(TMultiGraph* mg);
void draw_deformation_island(){
gROOT->SetStyle("tkn-histo");
for(int ii=0; ii<7; ii++) cols[ii] = TColor::GetColor(ccs[ii]);
tkmanager dtm;
TMultiGraph* mgr = new TMultiGraph;
TMultiGraph* mgs = new TMultiGraph;
TMultiGraph* mgb2 = new TMultiGraph;
TMultiGraph* mgr42 = new TMultiGraph;
for(int zz=32; zz<=45; zz+=2)
{
TGraph* gr = get_graph(zz);
TGraph* gs = get_graph(zz);
TGraph* gb2 = get_graph(zz);
TGraph* gr42 = get_graph(zz);
for(int nn=48; nn<=68; nn+=2)
{
auto nuc = dtm.get_nucleus(zz,nn+zz);
if(!nuc) continue;
if(nuc->has_property("radius"))
gr->SetPoint(gr->GetN(),nuc->get_n(), nuc->get("radius")->get_value());
if(nuc->has_property("twoNeutronSeparationEnergy")&&nuc->get("twoNeutronSeparationEnergy")->get_info()==tkn::tkproperty::kKnown)
gs->SetPoint(gs->GetN(),nuc->get_n(), nuc->get("twoNeutronSeparationEnergy")->get_value()*0.001);
auto lev_scheme = nuc->get_level_scheme();
auto gamma = lev_scheme->get_decay<tkgammadecay>("2+1->0+1",false);
if(!gamma){}
else{
auto BE2W = gamma->get_trans_prob(true,2,true);
if(BE2W<=0||std::isnan(BE2W)) {}
else gb2->SetPoint(gb2->GetN(),nuc->get_n(),BE2W/(nuc->get_a()));
}
auto lvl2 = lev_scheme->get_level("2+1",false);
if(!lvl2||!lvl2->has_property("energy")) {}
else{
double e2 = lvl2->get_energy();
gr42->SetPoint(gr42->GetN(),nuc->get_n(),e2*0.001);
}
}
if(gr->GetN()>1) mgr->Add(gr,"PL");
else delete gr;
if(gs->GetN()>1) mgs->Add(gs,"PL");
else delete gs;
if(gb2->GetN()>1) mgb2->Add(gb2,"PL");
else delete gb2;
if(gr42->GetN()>1) mgr42->Add(gr42,"PL");
else delete gr42;
}
TCanvas* cc = new TCanvas("Deformation","Deformation island",1200,1000);
cc->Divide(2,2);
cc->cd(1);
mgs->Draw("PAL");
mgs->GetYaxis()->SetTitle("S_{2n} (MeV)");
mgs->GetXaxis()->SetTitle("N");
get_line(mgs)->Draw();
cc->cd(2);
mgr->Draw("PAL");
mgr->GetYaxis()->SetTitle("r_{ch} (fm)");
mgr->GetXaxis()->SetTitle("N");
get_line(mgr)->Draw();
cc->cd(3);
mgb2->Draw("PAL");
mgb2->GetYaxis()->SetTitle("B(E2)/A (W.u./A)");
mgb2->GetXaxis()->SetTitle("N");
get_legend()->Draw();
get_line(mgb2)->Draw();
cc->cd(4);
mgr42->Draw("PAL");
mgr42->GetYaxis()->SetTitle("E(2^{+}) (MeV)");
mgr42->GetXaxis()->SetTitle("N");
get_line(mgr42)->Draw();
return;
}
TGraph* get_graph(int zz)
{
TGraph* gg = new TGraph;
gg->SetMarkerStyle(marks[zz/2-16]);
gg->SetMarkerSize(masize[zz/2-16]);
gg->SetMarkerColor(cols[zz/2-16]);
gg->SetLineColor(cols[zz/2-16]);
gg->SetName(Form("Z=%d",zz));
return gg;
}
TLine* get_line(TMultiGraph* mg)
{
TLine *ll = new TLine(60,mg->GetYaxis()->GetXmin(),60,mg->GetYaxis()->GetXmax());
ll->SetLineStyle(9);
return ll;
}
TLegend* get_legend()
{
TLegend *leg = new TLegend(0.203499,0.35485,0.542606,0.937585,NULL,"brNDC");
leg->SetBorderSize(0);
leg->SetTextFont(132);
leg->SetLineColor(1);
leg->SetLineStyle(1);
leg->SetLineWidth(1);
leg->SetFillColor(0);
leg->SetFillStyle(1001);
TLegendEntry *entry=0;
Int_t ci;
entry=leg->AddEntry("Z=44","Z=44","lpf");
entry->SetFillStyle(1000);
ci = TColor::GetColor("#120952");
entry->SetLineColor(ci);
entry->SetLineStyle(1);
entry->SetLineWidth(1);
entry->SetMarkerColor(ci);
entry->SetMarkerStyle(20);
entry->SetMarkerSize(1.3);
entry->SetTextFont(132);
entry=leg->AddEntry("Z=42","Z=42","lpf");
entry->SetFillStyle(1000);
ci = TColor::GetColor("#182562");
entry->SetLineColor(ci);
entry->SetLineStyle(1);
entry->SetLineWidth(1);
entry->SetMarkerColor(ci);
entry->SetMarkerStyle(20);
entry->SetMarkerSize(1.3);
entry->SetTextFont(132);
entry=leg->AddEntry("Z=40","Z=40","lpf");
entry->SetFillStyle(1000);
ci = TColor::GetColor("#1e427a");
entry->SetLineColor(ci);
entry->SetLineStyle(1);
entry->SetLineWidth(1);
entry->SetMarkerColor(ci);
entry->SetMarkerStyle(20);
entry->SetMarkerSize(1.3);
entry->SetTextFont(132);
entry=leg->AddEntry("Z=38","Z=38","lpf");
entry->SetFillStyle(1000);
ci = TColor::GetColor("#1d6190");
entry->SetLineColor(ci);
entry->SetLineStyle(1);
entry->SetLineWidth(1);
entry->SetMarkerColor(ci);
entry->SetMarkerStyle(20);
entry->SetMarkerSize(1.3);
entry->SetTextFont(132);
entry=leg->AddEntry("Z=36","Z=36","lpf");
entry->SetFillStyle(1000);
ci = TColor::GetColor("#1a83a5");
entry->SetLineColor(ci);
entry->SetLineStyle(1);
entry->SetLineWidth(1);
entry->SetMarkerColor(ci);
entry->SetMarkerStyle(20);
entry->SetMarkerSize(1.3);
entry->SetTextFont(132);
entry=leg->AddEntry("Z=34","Z=34","lpf");
entry->SetFillStyle(1000);
ci = TColor::GetColor("#1aabbf");
entry->SetLineColor(ci);
entry->SetLineStyle(1);
entry->SetLineWidth(1);
entry->SetMarkerColor(ci);
entry->SetMarkerStyle(20);
entry->SetMarkerSize(1.3);
entry->SetTextFont(132);
entry=leg->AddEntry("Z=32","Z=32","lpf");
entry->SetFillStyle(1000);
ci = TColor::GetColor("#1ddbdc");
entry->SetLineColor(ci);
entry->SetLineStyle(1);
entry->SetLineWidth(1);
entry->SetMarkerColor(ci);
entry->SetMarkerStyle(20);
entry->SetMarkerSize(1.3);
entry->SetTextFont(132);
return leg;
}
Stores information on a gamma-ray decay.
Definition: tkdecay.h:131
double get_trans_prob(bool _elec, int _L, bool _WU)
returns the gamma gamma transition probability value (ex: BE2 in Weiksopf units, get_trans_prob(1,...
Definition: tkdecay.cpp:262
Manages the database loading and provides access to the physics properties.
Definition: tkmanager.h:55
shared_ptr< tknucleus > get_nucleus(const tkstring &_nuc)
return a shared pointer to a nucleus from its name
Definition: tkmanager.h:126
Definition: tklog.cpp:39