KaliVeda
Toolkit for HIC analysis
KVLevelScheme.cpp
1 //Created by KVClassFactory on Fri Nov 20 12:23:35 2015
2 //Author: gruyer,,,
3 
4 #include "KVLevelScheme.h"
5 #include "KVFileReader.h"
6 #include "TH2F.h"
7 #include "TLine.h"
8 #include "TCanvas.h"
9 #include "TRandom.h"
10 #include "TStyle.h"
11 #include "TLatex.h"
12 using namespace std;
13 
15 
16 
17 
18 
19 
22 KVLevelScheme::KVLevelScheme(const char* symb)
23 {
24  // Default constructor
25  ncol = 1;
26  dy = 300;
27  dx = 1.3;
28  txs = 22;
29  ddx = .6;
30  hh = 0;
31  cc = 0;
32 
33  fCompNuc = new KVNucleus(symb);
34  fCompNuc->SetExcitEnergy(0);
35 
36  fFunc = new TF1("ExciEnergy", this, &KVLevelScheme::Evaluate, 0, 20000, 1, "KVLevelScheme", "Evaluate");
37  fFunc->SetNpx(5000);
38 
39  fFuncErel = new TF1("RelEnergy", this, &KVLevelScheme::EvaluateErel, 0, 20000, 1, "KVLevelScheme", "EvaluateErel");
40  fFuncErel->SetNpx(5000);
41 
42  InitStructure();
43 }
44 
45 
46 
49 
51 {
52  // Destructor
53  if (fFunc) delete fFunc;
54 }
55 
56 
57 
59 
61 {
62  TString dir = "";
63  KVBase::SearchKVFile("levels", dir, "data");
64  ReadLevels(Form("%s/%02d.dat", dir.Data(), fCompNuc->GetZ()));
65 }
66 
67 
68 
70 
72 {
73  KVFileReader fr;
74  fr.OpenFileToRead(file);
75 
76  while (fr.IsOK()) {
77  fr.ReadLine("|");
78  KVString tmp = "";
79 
80  Int_t npars = fr.GetNparRead();
81  if (npars == 0) break;
82 
83  Int_t aa, zz;
84  KVString jpi, gamma;
85  Double_t ee, de; //, gam, dgam;
86  ee = de = 0.;
87 
88  aa = fr.GetIntReadPar(0);
89  zz = fr.GetIntReadPar(1);
90 
91  if (aa != fCompNuc->GetA() || zz != fCompNuc->GetZ()) continue;
92 
93  ee = fr.GetDoubleReadPar(2);
94 
95  tmp = fr.GetReadPar(3);
97  if (!tmp.EqualTo("")) de = tmp.Atof();
98  else de = 0.;
99 
100  tmp = fr.GetReadPar(4);
102  if (!tmp.EqualTo("")) {
103  jpi = tmp.Data();
104  jpi.ReplaceAll("GE", "");
105  jpi.ReplaceAll(" & ", ",");
106  jpi.ReplaceAll("LE", "");
107  jpi.ReplaceAll("LT", "");
108  jpi.ReplaceAll("AP", "");
109  jpi.ReplaceAll(" OR ", "or");
110  jpi.ReplaceAll(" TO ", "to");
111  }
112  else jpi = "";
113 
114  gamma = "";
115  if (npars == 6) {
116  gamma = fr.GetReadPar(5);
117  }
118 
119  AddResonance(ee, jpi.Data(), gamma.Data());
120  }
121  fr.CloseFile();
122 }
123 
124 
125 
127 
128 void KVLevelScheme::AddResonance(Double_t ex, const char* jpi, const char* gam)
129 {
130  KVExcitedState* ll = new KVExcitedState;
131  ll->set(ex, jpi, gam);
132  fLevels.AddLast(ll);
133 }
134 
135 
136 
138 
140 {
141  cout << endl;
142  cout << " n." << " Jpi" << " E" << " T1/2" << endl;
143  cout << endl;
144  for (Int_t ii = (int)fLevels.GetSize() - 1; ii >= 0; ii--) {
145  KVExcitedState* ll = (KVExcitedState*)fLevels.At(ii);
146  cout << Form("%3d", ii) << ".";
147  ll->print();
148  }
149  cout << endl;
150  cout << " " << Form("%5s", fCompNuc->GetSymbol()) << endl;
151  cout << endl;
152 }
153 
154 
155 
157 
158 double KVLevelScheme::func(double xx, double tt)
159 {
160  double yy = 0.;
161  for (Int_t ii = 0; ii < fLevels.GetSize(); ii++) {
162  KVExcitedState* ll = (KVExcitedState*)fLevels.At(ii);
163  yy += ll->eval(xx);
164  }
165  yy *= TMath::Exp(-xx / tt);
166  return yy;
167 }
168 
169 
170 
172 
173 double KVLevelScheme::getWidth(double gam, TString unit)
174 {
175  unit.ToUpper();
176  if (unit.EqualTo("EV")) gam *= 1e-3;
177  else if (unit.EqualTo("KEV")) gam *= 1;
178  else if (unit.EqualTo("MEV")) gam *= 1e3;
179  else if (unit.Contains("S")) gam *= 1e-6;
180  return gam;
181 }
182 
183 
184 
186 
188 {
189  jpi.ReplaceAll("(", "");
190  jpi.ReplaceAll(")", "");
191 
192  int j = 0;
193 
194  if (jpi.EqualTo("")) return 0;
195  if (jpi.EqualTo("GE")) return 0;
196 
197  if (jpi.Contains("-")) jpi = jpi(0, jpi.Index("-"));
198  if (jpi.Contains("+")) jpi = jpi(0, jpi.Index("+"));
199  if (jpi.Contains("/2")) jpi = jpi(0, jpi.Index("/"));
200  if (jpi.Contains(",")) jpi = jpi(0, jpi.Index(","));
201 
202  j = jpi.Atoi();
203  return j;
204 }
205 
206 
207 
209 
211 {
212  if (il > (int)fLevels.GetSize()) return 0;
213  else return ((KVExcitedState*) fLevels.At(il))->fEnergy;
214 }
215 
216 
217 
219 
221 {
222  if (il > (int)fLevels.GetSize()) return 0;
223  else return ((KVExcitedState*) fLevels.At(il))->fWidth;
224 }
225 
226 
227 
229 
231 {
232  if (il > (int)fLevels.GetSize()) return 0;
233  else return ((KVExcitedState*) fLevels.At(il))->fSpin;
234 }
235 
236 
237 
239 
241 {
242  if (il > (int)fLevels.GetSize()) return 0;
243  else return ((KVExcitedState*) fLevels.At(il))->fParity;
244 }
245 
246 
247 
249 
250 const char* KVLevelScheme::GetJPiStr(int il)
251 {
252  if (il > (int)fLevels.GetSize()) return "";
253  else return ((KVExcitedState*) fLevels.At(il))->fJPi.Data();
254 }
255 
256 
258 
259 const char* KVLevelScheme::GetEGammaStr(int il)
260 {
261  if (il > (int)fLevels.GetSize()) return 0;
262  KVString gam = ((KVExcitedState*) fLevels.At(il))->fGamma.Data();
264  gam.ReplaceAll("-1.0", "-");
265  return Form("%d (%s)", TMath::Nint(((KVExcitedState*) fLevels.At(il))->fEnergy), gam.Data());
266 }
267 
268 
269 
271 
273 {
274  jpi.ReplaceAll("(", "");
275  jpi.ReplaceAll(")", "");
276 
277 // int j = 0;
278  Int_t pi = 0;
279 
280  if (jpi.EqualTo("")) return 0;
281  if (jpi.EqualTo("GE")) return 0;
282 
283  // determination of pi
284  if ((jpi.Contains("-")) && (jpi.Contains("+"))) {
285  if (jpi.Index("+") < jpi.Index("-")) pi = 1;
286  else pi = -1;
287  }
288  else if (jpi.Contains("-")) pi = -1;
289  else pi = 1;
290 
291  return pi;
292 }
293 
294 
295 
296 
298 
299 double KVLevelScheme::Evaluate(double* x, double* p)
300 {
301  double xx = x[0];
302  double tt = p[0];
303  return func(xx, tt);
304 }
305 
306 
307 
309 
310 double KVLevelScheme::EvaluateErel(double* x, double* p)
311 {
312  double xx = x[0];
313  double tt = p[0];
314  return func(xx - fQvalue * 1000, tt);
315 }
316 
317 
318 
320 
321 void KVLevelScheme::GetParticlesFromErel(KVNucleus* n1, KVNucleus* n2, double erel, bool randAngle, TVector3* vsrc)
322 {
323  KVNucleus nuc1(n1->GetZ(), n1->GetA());
324  KVNucleus nuc2(n2->GetZ(), n2->GetA());
325 
326  double m1 = nuc1.GetMass();
327  double m2 = nuc2.GetMass();
328  double e1 = (m2 / (m1 + m2)) * erel;
329  double e2 = (m1 / (m1 + m2)) * erel;
330 
331  nuc1.SetKE(e1);
332  nuc1.SetTheta(0);
333  nuc2.SetKE(e2);
334  nuc2.SetTheta(180);
335 
336  if (randAngle) {
337  TVector3 v1r = nuc1.GetVelocity();
338  TVector3 v2r = nuc2.GetVelocity();
339  TRotation rr;
340  rr.SetXEulerAngles(gRandom->Rndm() * 2.*TMath::Pi(), TMath::ACos(gRandom->Rndm() * 2. - 1.), gRandom->Rndm() * 2.*TMath::Pi());
341  v1r *= rr;
342  nuc1.SetVelocity(v1r);
343  v2r *= rr;
344  nuc2.SetVelocity(v2r);
345  }
346 
347  if (vsrc) {
348  *vsrc *= -1;
349  nuc1.SetFrame("src", *vsrc);
350  nuc2.SetFrame("src", *vsrc);
351  n1->SetVelocity(nuc1.GetFrame("src")->GetVelocity());
352  n2->SetVelocity(nuc2.GetFrame("src")->GetVelocity());
353  }
354 
355 }
356 
357 
358 
360 
361 void KVLevelScheme::GetParticlesFromExci(KVNucleus* nuc1, KVNucleus* nuc2, double erel, bool randAngle, TVector3* vsrc)
362 {
363  KVNucleus comp = *nuc1 + *nuc2;
364  double qq = -1 * (comp.GetExcitEnergy());
365  GetParticlesFromErel(nuc1, nuc2, erel + qq, randAngle, vsrc);
366 }
367 
368 
369 
371 
372 void KVLevelScheme::GetRandomParticles(KVNucleus* n1, KVNucleus* n2, double T/*keV*/, bool randAngle, TVector3* vsrc)
373 {
374  if (!fLevels.GetSize()) cout << "KVLevelScheme::GetRandomParticles: please initialize the level scheme firts..." << endl;
375  fFunc->SetParameter(0, T);
376  double excit = fFunc->GetRandom();
377 
378  GetParticlesFromExci(n1, n2, excit, randAngle, vsrc);
379 }
380 
381 
382 
384 
385 void KVLevelScheme::SetDrawStyle(double deMin, double fullWidth, double lineWidth, int textSize)
386 {
387  dy = deMin;
388  dx = fullWidth;
389  ddx = lineWidth;
390  txs = textSize;
391 }
392 
393 
394 
396 
397 void KVLevelScheme::Draw(Option_t* /*option*/)
398 {
399  ncol = 1;
400  int icol = 0;
401 
402  int cols[20];
403  for (int ic = 0; ic < 20; ic++) cols[ic] = 0;
404  cols[0] = GetLevelEnergy(0);
405 
406  for (int ii = 0; ii < GetNLevels(); ii++) {
407  for (int ic = 0; ic < TMath::Max(ii, 20); ic++) {
408  if ((GetLevelEnergy(ii) - cols[ic] > dy)) {
409  cols[ic] = GetLevelEnergy(ii);
410  icol = ic;
411  if (icol + 1 > ncol) ncol = icol + 1;
412  break;
413  }
414  }
415  }
416 
417  double max = GetLevelEnergy(GetNLevels() - 1);
418  if (hh) delete hh;
419  hh = new TH2F(Form("dumhist%s", fCompNuc->GetSymbol()), "", 1, 0, ncol * dx + 0.5 * dx, 1000, -200, max + 200);
420 
421 // TString opt = option;
422  gStyle->SetOptStat(0);
423 // if (!opt.Contains("same")) {
424  if (cc) delete cc;
425  cc = new TCanvas(Form("levels%s", fCompNuc->GetSymbol()), Form("levels%s", fCompNuc->GetSymbol()), (ncol) * (dx) * 0.4 * 400, 800);
426  cc->SetTickx(1);
427  cc->SetTicky(1);
428  cc->SetTopMargin(0.02);
429  cc->SetBottomMargin(0.02);
430  cc->SetRightMargin(0.02);
431  cc->SetLeftMargin(0.02);
432 
433  hh->GetXaxis()->SetAxisColor(0);
434  hh->GetYaxis()->SetAxisColor(0);
435 
436  hh->GetXaxis()->SetLabelSize(0);
437  hh->GetYaxis()->SetNdivisions(0);
438 
439  hh->Draw();
440 // }
441 
442  for (int ic = 0; ic < 20; ic++) cols[ic] = 0;
443  cols[0] = GetLevelEnergy(0);
444 
445  icol = 0;
446  for (int ii = 0; ii < GetNLevels(); ii++) {
447  for (int ic = 0; ic < TMath::Max(ii, 20); ic++) {
448  if ((GetLevelEnergy(ii) - cols[ic] > dy)) {
449  cols[ic] = GetLevelEnergy(ii);
450  icol = ic;
451  break;
452  }
453  }
454 
455  TLine* ll = new TLine(0.5 + dx * icol, GetLevelEnergy(ii), ddx + .5 + dx * icol, GetLevelEnergy(ii));
456  ll->Draw("same");
457  TString jpi = GetJPiStr(ii);
458  jpi.ReplaceAll("+", "^{+}");
459  jpi.ReplaceAll("-", "^{-}");
460  TLatex* tex = new TLatex(.4 + dx * icol, GetLevelEnergy(ii), jpi.Data());
461  tex->SetTextAlign(32);
462  tex->SetTextFont(133);
463  tex->SetTextSize(txs);
464  // tex->Draw();
465 
466  tex = new TLatex(ddx + .6 + dx * icol, GetLevelEnergy(ii), Form("%d", TMath::Nint(GetLevelEnergy(ii)))); //GetEGammaStr(ii));
467  tex->SetTextAlign(12);
468  tex->SetTextFont(133);
469  tex->SetTextSize(txs);
470  tex->Draw();
471  }
472 
473  TLatex* tte = new TLatex(0.5, 0.95, Form("^{%d}%s", fCompNuc->GetA(), fCompNuc->GetSymbol("EL")));
474  tte->SetNDC(1);
475  tte->SetTextAlign(23);
476  tte->SetTextFont(133);
477  tte->SetTextSize(txs);
478  tte->Draw();
479 
480  hh->GetXaxis()->SetLimits(0, dx * (ncol + 0.5));
481 }
482 
483 
484 
486 
487 void KVLevelScheme::DrawThreshold(const char* symb, Option_t* option, double ex)
488 {
489  TString opt = option;
490  ncol++;
491  KVNucleus a(symb);
492 
493  KVNucleus tmp = *fCompNuc - a;
494  double qa = tmp.GetExcitEnergy() * -1.;
495 
496  TLatex* tte = 0;
497  TLine* lq = 0;
498 
499  TString decay = Form("^{%d}%s + ^{%d}%s", a.GetA(), a.GetSymbol("EL"), tmp.GetA(), tmp.GetSymbol("EL"));
500  tte = new TLatex(0.5 + (0.5 * ddx) + dx * (ncol - 1), qa * 1000 - 100, decay.Data());
501  tte->SetTextAlign(23);
502  tte->SetTextFont(133);
503  tte->SetTextSize(txs);
504  tte->Draw();
505 
506 
507  if (opt.Contains("l")) {
508  lq = new TLine(0.5, qa * 1000, 0.5 + dx * (ncol - 1), qa * 1000);
509  lq->SetLineColor(kGray);
510  lq->SetLineStyle(7);
511  lq->Draw();
512  }
513 
514  lq = new TLine(0.5 + dx * (ncol - 1), qa * 1000, 0.5 + ddx + dx * (ncol - 1), qa * 1000);
515  lq->Draw();
516 
517  tte = new TLatex(ddx + .6 + dx * (ncol - 1), qa * 1000, Form("%d", TMath::Nint(1000 * qa)));
518  tte->SetTextAlign(12);
519  tte->SetTextFont(133);
520  tte->SetTextSize(txs);
521  tte->Draw();
522 
523  if (ex > 0) {
524  a.SetExcitEnergy(ex / 1000.);
525  tmp = *fCompNuc - a;
526  qa = tmp.GetExcitEnergy() * -1.;
527 
528  if (opt.Contains("l")) {
529  lq = new TLine(0.5, qa * 1000, 0.5 + dx * (ncol - 1), qa * 1000);
530  lq->SetLineColor(kGray);
531  lq->SetLineStyle(7);
532  lq->Draw();
533  }
534 
535  lq = new TLine(0.5 + dx * (ncol - 1), qa * 1000, 0.5 + ddx + dx * (ncol - 1), qa * 1000);
536  lq->Draw();
537  tte = new TLatex(ddx + .6 + dx * (ncol - 1), qa * 1000, Form("%d", TMath::Nint(1000 * qa)));
538  tte->SetTextAlign(12);
539  tte->SetTextFont(133);
540  tte->SetTextSize(txs);
541  tte->Draw();
542  }
543  cc->SetWindowSize((ncol) * (dx) * 0.35 * 400, 800);
544  hh->GetXaxis()->SetLimits(0, dx * (ncol + 0.5));
545 }
546 
547 
548 
550 
551 double KVLevelScheme::GetThreshold(const char* outnuc)
552 {
553  if (!strcmp(outnuc, "")) outnuc = "1H";
554  KVNucleus a(outnuc);
555  KVNucleus tmp = *fCompNuc - a;
556  return tmp.GetExcitEnergy() * -1000.; //returns the threshold in keV
557 }
558 
559 
560 
562 
563 void KVLevelScheme::SetDecayProduct(KVNucleus* nuc, double excit_energy)
564 {
565  fDecayProd = nuc;
566  fDecayProd->SetExcitEnergy(excit_energy * 0.001);
567  KVNucleus tmp = *fCompNuc - *fDecayProd;
568  fQvalue = tmp.GetExcitEnergy() * -1.;
569 }
570 
571 
int Int_t
#define e(i)
double Double_t
const char Option_t
kGray
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
int * lq
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
R__EXTERN TStyle * gStyle
static Bool_t SearchKVFile(const Char_t *name, TString &fullpath, const Char_t *kvsubdir="")
Definition: KVBase.cpp:517
Excited state of atomic nucleus.
void set(Double_t ee, Double_t ww, Double_t jj, Int_t pi)
Double_t eval(Double_t excit)
Handle reading columns of numeric data in text files.
Definition: KVFileReader.h:121
void CloseFile()
Definition: KVFileReader.h:237
ReadStatus ReadLine(const KVString &pattern="")
Definition: KVFileReader.h:243
Double_t GetDoubleReadPar(Int_t pos) const
Definition: KVFileReader.h:334
Int_t GetIntReadPar(Int_t pos) const
Definition: KVFileReader.h:338
Int_t GetNparRead() const
Definition: KVFileReader.h:325
Bool_t IsOK()
Definition: KVFileReader.h:231
KVString GetReadPar(Int_t pos) const
Definition: KVFileReader.h:342
Bool_t OpenFileToRead(const KVString &filename)
Definition: KVFileReader.h:210
Tool to simulate nucleus multi-particle decay.
Definition: KVLevelScheme.h:24
Int_t GetLevelParity(int il)
double EvaluateErel(double *x, double *p)
void SetDrawStyle(double deMin=300., double fullWidth=1.3, double lineWidth=0.6, int textSize=22)
const char * GetEGammaStr(int il)
void DrawThreshold(const char *symb, Option_t *option="", double ex=0.)
double GetThreshold(const char *outnuc="")
static void GetParticlesFromExci(KVNucleus *n1, KVNucleus *n2, double erel, bool randAngle=0, TVector3 *vsrc=0)
virtual ~KVLevelScheme()
Destructor.
int getJ(TString jpi)
int getPI(TString jpi)
static void GetParticlesFromErel(KVNucleus *n1, KVNucleus *n2, double erel, bool randAngle=0, TVector3 *vsrc=0)
Double_t GetLevelEnergy(int il)
Int_t GetLevelSpin(int il)
void SetDecayProduct(KVNucleus *nuc, double excit_energy=0)
void Draw(Option_t *option="")
void ReadLevels(const char *file)
double getWidth(double gam, TString unit)
double Evaluate(double *x, double *p)
const char * GetJPiStr(int il)
void AddResonance(Double_t ex, const char *jpi, const char *gam)
void GetRandomParticles(KVNucleus *n1, KVNucleus *n2, double T, bool randAngle=0, TVector3 *vsrc=0)
void InitStructure()
double func(double xx, double tt)
Double_t GetLevelWidth(int il)
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
const Char_t * GetSymbol(Option_t *opt="") const
Definition: KVNucleus.cpp:71
void SetExcitEnergy(Double_t e)
Definition: KVNucleus.cpp:858
Double_t GetExcitEnergy() const
Definition: KVNucleus.h:279
Int_t GetA() const
Definition: KVNucleus.cpp:792
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:763
void SetTheta(Double_t theta)
Definition: KVParticle.h:696
void SetVelocity(const TVector3 &)
Set velocity of particle (in cm/ns units)
void SetKE(Double_t ecin)
Definition: KVParticle.cpp:246
void SetFrame(const Char_t *frame, const KVFrameTransform &)
Definition: KVParticle.cpp:775
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
Definition: KVParticle.cpp:897
Double_t GetMass() const
Definition: KVParticle.h:577
TVector3 GetVelocity() const
returns velocity vector in cm/ns units
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
void RemoveAllExtraWhiteSpace()
Definition: KVString.cpp:1253
virtual void SetTextAlign(Short_t align=11)
virtual void SetTextFont(Font_t tfont=62)
virtual void SetTextSize(Float_t tsize=1)
virtual void Draw(Option_t *option="")
Double_t Rndm() override
TRotation & SetXEulerAngles(Double_t phi, Double_t theta, Double_t psi)
Int_t Atoi() const
Double_t Atof() const
const char * Data() const
Bool_t EqualTo(const char *cs, ECaseCompare cmp=kExact) const
void ToUpper()
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
TString & ReplaceAll(const char *s1, const char *s2)
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
void SetOptStat(Int_t stat=1)
virtual void SetNDC(Bool_t isNDC=kTRUE)
Double_t x[n]
Double_t ex[n]
double gamma(double x)
double T(double x)
double max(double x, double y)
Double_t ACos(Double_t)
Int_t Nint(T x)
Double_t Exp(Double_t x)
constexpr Double_t Pi()
Double_t Max(Double_t a, Double_t b)
TArc a
auto * tt
ClassImp(TPyArg)