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