TkN 2.1
Toolkit for Nuclei
deformation.C
1#include "TH2F.h"
2#include "TStyle.h"
3#include "TGraph2D.h"
4#include "TGraph.h"
5#include "TLatex.h"
6#include "TLegend.h"
7#include "TView.h"
8#include "TEllipse.h"
9#include "TMath.h"
10#include "TArrow.h"
11#include "TROOT.h"
12
13#include "tkmanager.h"
14#include "tknuclear_chart.h"
15
16#define COLORS 20
17
18using namespace tkn;
19
20TH2F* setup_histo(double ts = .05);
21void draw_legend(int palette[COLORS], double betamax, double ts = .05);
22
23void draw_deformation()
24{
25 gROOT->SetStyle("tkn-histo");
26 gStyle->SetOptStat(0);
27 gStyle->SetOptTitle(0);
28
29 glog.set_warnings(false);
30
31 int cs = 1300;
32 TCanvas *c = new TCanvas("c","c",cs,cs*.7);
33 c->SetLeftMargin(0.15);
34 c->SetBottomMargin(0.15);
35 c->SetRightMargin(0.05);
36 c->SetTopMargin(0.05);
37 c->SetTickx();
38 c->SetTicky();
39 c->SetBorderMode(0);
40 c->SetBorderSize(0);
41
42 TH2F *h = setup_histo();
43 h->Draw();
44
45 Double_t r[] = {0.0, 0.0, 1.0, 1.0};
46 Double_t g[] = {0.0, 1.0, 1.0, 0.0};
47 Double_t b[] = {1.0, 1.0, 0.0, 0.0};
48 Double_t stop[] = {0.0, 0.3, 0.7, 1.0};
49
50 int colu=TColor::CreateGradientColorTable(4, stop, r, g, b, COLORS);
51 int palette[COLORS];
52 for(int k=0;k<COLORS;k++) palette[k]=colu+k;
53
54 double betamax = .6;
55 int col=0;
56
58 for(const auto &nuc : dtm.get_nuclei([](auto nuc) {return (nuc->get_z()%2==0 && nuc->get_n()%2==0);})) {
59
60 auto lev_scheme = nuc->get_level_scheme();
61
62 auto lvl2 = lev_scheme->get_level("2+1",true);
63 auto lvl4 = lev_scheme->get_level("4+1",true);
64 if(!lvl2||!lvl4) continue;
65
66 if(!lvl2->has_property("energy")||!lvl4->has_property("energy")) continue;
67 double e2 = lvl2->get_energy();
68 double e4 = lvl4->get_energy();
69 double r42 = e4/e2;
70
71 auto gamma = lev_scheme->get_decay<tkgammadecay>("2+1->0+1",false);
72 if(!gamma) continue;
73 auto BE2W = gamma->get_trans_prob(true,2,true);
74 if(BE2W<=0||std::isnan(BE2W)) continue;
75
76 auto quad_def = nuc->get("quadrupoleDeformation");
77 if(!quad_def) continue;
78
79 double beta = quad_def->get_value();
80 if(beta>betamax) continue;
81 col = palette[TMath::Nint((COLORS-1)*beta/betamax)];
82
83 double radius = .008*sqrt((1+beta)*pow(nuc->get_a(),2./3.));
84 TEllipse* el = new TEllipse(r42,BE2W/nuc->get_a(),radius,radius*1.2*(1+beta));
85 el->SetLineWidth(0);
86 el->SetFillColorAlpha(col,.7);
87 el->SetLineColorAlpha(col,0);
88 el->Draw("same");
89 }
90
91 draw_legend(palette, betamax);
92}
93
94TH2F *setup_histo(double ts)
95{
96 TH2F *h = new TH2F("h","h",100,1,3.5,100,-0.1,2);
97 h->GetXaxis()->SetNdivisions(209);
98 h->GetYaxis()->SetNdivisions(209);
99 h->GetXaxis()->SetLabelSize(ts);
100 h->GetYaxis()->SetLabelSize(ts);
101 h->GetXaxis()->SetTitleSize(ts+.01);
102 h->GetYaxis()->SetTitleSize(ts+.01);
103 h->GetXaxis()->SetTickLength(.02);
104 h->GetYaxis()->SetTickLength(.02*.7);
105 h->GetXaxis()->SetTitleOffset(1.1);
106 h->GetYaxis()->SetTitle("B(E2)/A (in Weisskopf unit)");
107 h->GetXaxis()->SetTitle("E(4^{+})/E(2^{+})");
108 return h;
109}
110
111void draw_legend(int palette[], double betamax, double ts)
112{
113 TEllipse* el = 0;
114 TLatex* txt = 0;
115
116 double yy = 1.7;
117 double dyy = .02;
118 int col = 0;
119
120 double val =betamax;
121 el = new TEllipse(1.2,yy,.05,.05*1.2*(1+val));
122 col = palette[TMath::Nint((COLORS-1)*val/betamax)];
123 el->SetFillColorAlpha(col,.7);
124 el->SetLineColorAlpha(kBlack,0);
125 el->Draw("same");
126
127 txt = new TLatex(1.3,yy, "0.6");
128 txt->SetTextAlign(12);
129 txt->SetTextFont(132);
130 txt->SetTextSize(ts);
131 txt->Draw("same");
132
133 yy-= .05*1.2*(1+val)+dyy;
134 val = 0.4;
135 yy-= .05*1.2*(1+val);
136
137 el = new TEllipse(1.2,yy,.05,.05*1.2*(1+val));
138 col = palette[TMath::Nint((COLORS-1)*val/betamax)];
139 el->SetFillColorAlpha(col,.7);
140 el->SetLineColorAlpha(kBlack,0);
141 el->Draw("same");
142
143 txt = new TLatex(1.3,yy, "0.4");
144 txt->SetTextAlign(12);
145 txt->SetTextFont(132);
146 txt->SetTextSize(ts);
147 txt->Draw("same");
148
149 yy-= .05*1.2*(1+val)+dyy;
150 val = 0.2;
151 yy-= .05*1.2*(1+val);
152
153 el = new TEllipse(1.2,yy,.05,.05*1.2*(1+val));
154 col = palette[TMath::Nint((COLORS-1)*val/betamax)];
155 el->SetFillColorAlpha(col,.7);
156 el->SetLineColorAlpha(kBlack,0);
157 el->Draw("same");
158
159 txt = new TLatex(1.3,yy, "0.2");
160 txt->SetTextAlign(12);
161 txt->SetTextFont(132);
162 txt->SetTextSize(ts);
163 txt->Draw("same");
164
165 yy-= .05*1.2*(1+val)+dyy;
166 val = 0.0;
167 yy-= .05*1.2*(1+val);
168
169 el = new TEllipse(1.2,yy,.05,.05*1.2*(1+val));
170 col = palette[TMath::Nint((COLORS-1)*val/betamax)];
171 el->SetFillColorAlpha(col,.7);
172 el->SetLineColorAlpha(kBlack,0);
173 el->Draw("same");
174
175 txt = new TLatex(1.3,yy, "0.0");
176 txt->SetTextAlign(12);
177 txt->SetTextFont(132);
178 txt->SetTextSize(ts);
179 txt->Draw("same");
180}
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
vector< shared_ptr< tknucleus > > get_nuclei(std::function< bool(shared_ptr< tknucleus >)>const &_selection)
return a vector containing all the known nuclei filtered by the lambda function
Definition: tkmanager.cpp:367
Definition: tklog.cpp:39