KaliVeda
Toolkit for HIC analysis
KVIsoscaling.cpp
1 #include "KVIsoscaling.h"
2 #include "KVString.h"
3 
4 #include "TGraphErrors.h"
5 #include "TString.h"
6 #include "TObjArray.h"
7 #include "TObjString.h"
8 #include "TMath.h"
9 #include "TF1.h"
10 #include "TFile.h"
11 #include "TAxis.h"
12 
13 #include <fstream>
14 
16 
17 // BEGIN_HTML <!--
19 /* -->
20 <h2>KVIsoscaling</h2>
21 <h4>KVClass</h4>
22 <!-- */
23 // --> END_HTML
25 
26 //____________________________________________________________________________//
27 
28 
36 void KVIsoscaling::ReadYieldsFile(const std::string& system_name, const Char_t* file_path)
37 {
38  // Read and save yields from file provided by the user and fill the corresponding structs and maps
39  // Input ASCII file format is :
40  // Z A Yield Yield_err
41  //
42  // \param[in] system name
43  // \param[in] path for the input ASCII file
44 
45  Info("ReadYieldsFile", "===== Starting reading system %s =====\n", system_name.c_str());
46 
47  // --- Open the file for a new system ---
48  std::ifstream my_file;
49  my_file.open(file_path);
50 
51  if (my_file.is_open()) {
52 
53  // - Save file path -
54  fSystemList_[system_name].file_path = std::string(file_path);
55 
56  // --- Read the file ---
57  Int_t nline = 0;
58  std::string line;
59  while (getline(my_file, line)) {
60  TString my_str(line);
61 
62  // --- Ignore comments ---
63  if (my_str.BeginsWith("#")) continue;
64 
65  // --- Extract infos ---
66  TObjArray* my_str_array = my_str.Tokenize(" \t,;");
67 
68  Int_t zz = ((TObjString*) my_str_array->At(0))->String().Atoi();
69  Int_t aa = ((TObjString*) my_str_array->At(1))->String().Atoi();
70  Float_t yy = ((TObjString*) my_str_array->At(2))->String().Atof();
71  Float_t yy_err = ((TObjString*) my_str_array->At(3))->String().Atof();
72 
73  // save new data
74  fSystemList_[system_name].yields[zz][aa].yield = yy;
75  fSystemList_[system_name].yields[zz][aa].yield_err = yy_err;
76 
77  // --- Debug ---
78  if (fdebug_) printf("(line=%d) (zz=%d, aa=%d, yy=%lf, yy_err=%lf)\n", nline, zz, aa, yy, yy_err);
79 
80  nline++;
81  }//end getline
82 
83  // --- Close file ---
84  my_file.close();
85 
86  //debug
87  Info("ReadYieldsFile", "[system_name=%s] file %s closed, masses and yields saved\n", system_name.c_str(), file_path);
88 
89  }//end file.is_opened
90 
91  // --- Handle the errors ---
92  else {
93  Error("ReadYieldsFile", "!!! Can't open file '%s' !!!", file_path);
94  }
95 
96  // --- Apply a gaussian fit to the provided yields ---
97  BuildGaussianPlots(system_name);
98 }
99 
100 
101 
107 
108 void KVIsoscaling::BuildGaussianPlots(const std::string& system_name)
109 {
110  // Method to build the plots used to check the gaussian approximation
111  // Also compute the average masses for each charge \f$\langle A(Z)\rangle\f$
112  //
113  // \param[in] system name
114 
115  if (fdebug_) Info("BuildGaussianPlots", "====== STARTING BUILDING GAUSSIAN PLOTS FOR SYSTEM %s ======", system_name.c_str());
116 
117  // - Find the yields for the given system -
118  YieldData_t my_yields = fSystemList_[system_name].yields;
119 
120  // - Iterate over the yields for the given element -
121  for (auto it = my_yields.begin(); it != my_yields.end(); it++){
122  Int_t zz = it->first;
123  Element_t my_element = it->second;
124 
125  // - Create then the graph for the given element -
126  TGraphErrors* gr = new TGraphErrors();
127  gr->SetName(Form("gr_Yield_vs_N_Z%d_%s", zz, system_name.c_str()));
128  gr->GetXaxis()->SetTitle("N");
129  gr->GetYaxis()->SetTitle("Yield");
130  gr->SetMarkerStyle(20);
131 
132  TGraphErrors* grln = new TGraphErrors();
133  grln->SetName(Form("gr_LnYield_vs_N_Z%d_%s", zz, system_name.c_str()));
134  grln->GetXaxis()->SetTitle("N");
135  grln->GetYaxis()->SetTitle("Ln(Y)");
136  grln->SetMarkerStyle(20);
137 
138  if (fdebug_) Info("BuildGaussianPlots", "Graph %s created", gr->GetName());
139 
140  // - Iterate over isotopes of the given element and fill the graph -
141  Int_t npp = 0;
142 
143  for (auto itt = my_element.begin(); itt != my_element.end(); itt++){
144  Int_t aa = itt->first;
145  IsotopeYield_t my_isoyield = itt->second;
146 
147  Int_t nn = aa-zz;
148  Float_t yield = my_isoyield.yield;
149  Float_t yielderr = my_isoyield.yield_err;
150 
151  gr->SetPoint(npp, nn, yield);
152  gr->SetPointError(npp, 0, yielderr);
153 
154  grln->SetPoint(npp, nn, TMath::Log(yield));
155  grln->SetPointError(npp, 0, TMath::Abs(yielderr / yield));
156 
157  if (fdebug_) Info("BuildGaussianPlots", "(Z=%d, N=%d, A=%d) Yield=%lf +/- %lf", zz, nn, aa, yield, yielderr);
158 
159  npp++;
160  }//end isotope iteration
161 
162  // --- Apply individual gaussian fit and save the results ---
163  TF1* gaus = new TF1(Form("gausfit_yield_vs_N_Z%d_%s", zz, system_name.c_str()), "gaus(0)", 0, 60);
164  gr->Fit(Form("gausfit_yield_vs_N_Z%d_%s", zz, system_name.c_str()), "Q");
165 
166  // --- Save results ---
167  fSystemList_[system_name].graphs[zz] = gr;
168  fSystemList_[system_name].lngraphs[zz] = grln;
169  fSystemList_[system_name].gauss[zz] = gaus;
170 
171  }//end yields iteration
172 
173  if (fdebug_) Info("BuildGaussianPlots", "====== END BUILDING GAUSSIAN PLOTS FOR SYSTEM %s ======", system_name.c_str());
174 }
175 
176 
177 
188 
189 void KVIsoscaling::TestGaussianApprox(const std::string& system_name1, const std::string& system_name2, Int_t zz, Double_t tol)
190 {
191  // Method to test the tolerance threshold of the gaussian approximation for a given charge (element) \f$Z\f$.
192  // Assuming a tolerance \f$tol\f$ set by the user and an average mass \f$\langle A(Z)\rangle\f$ with a standard deviation \f$ rms(Z) \f$,
193  // the region of experimental points will be limited to \f$ \langle A(Z)\rangle \pm tol \cdot rms(Z) \f$.
194  // This method allows to test and plot the results of various tolerance thresholds.
195  //
196  // \param[in] system_name1 name of the n-deficient system
197  // \param[in] system_name2 name of the n-rich system
198  // \param[in] zz fixed charge for the isotopic yields to be tested
199  // \param[in] tol threshold to be tested
200 
201  TGraphErrors* gr_sys1_all = fSystemList_[system_name1].graphs[zz];
202  TGraphErrors* gr_sys2_all = fSystemList_[system_name2].graphs[zz];
203  TF1* gaus1 = fSystemList_[system_name1].gauss[zz];
204  TF1* gaus2 = fSystemList_[system_name2].gauss[zz];
205 
206  // --- Build the graphs with thresholds and ratios of the yields ---
207  // Use only the A (N) shared between the 2 systems
208  TGraphErrors* gr_sys1_tol = new TGraphErrors();
209  gr_sys1_tol->SetName(Form("gr_yields_vs_N_Z%d_%s_tol", zz, system_name1.c_str()));
210  TGraphErrors* gr_sys2_tol = new TGraphErrors();
211  gr_sys2_tol->SetName(Form("gr_yields_vs_N_Z%d_%s_tol", zz, system_name1.c_str()));
212 
213  TGraphErrors* gr_ratio_all = new TGraphErrors();
214  gr_ratio_all->SetName(Form("gr_ratio_all_Z%d_%s_%s", zz, system_name1.c_str(), system_name2.c_str()));
215  TGraphErrors* gr_ratio_tol = new TGraphErrors();
216  gr_ratio_tol->SetName(Form("gr_ratio_tol_Z%d_%s_%s", zz, system_name1.c_str(), system_name2.c_str()));
217 
218  // - Extract all isotopes shared between the 2 systems -
219  KVNumberList nl_inter_a = GetSharedANumberList(system_name1.c_str(), system_name2.c_str(), zz);
220 
221  if (fdebug_) Info("TestGaussianApproxe", "=== Starting iteration over isotope list='%s' for Z=%d ===", nl_inter_a.GetList(), zz);
222 
223  // - Build graphs with tolerance applied -
224  Int_t npoint_all = 0;
225  Int_t np1 = 0;
226  Int_t np2 = 0;
227  Int_t npoint_tol = 0;
228  nl_inter_a.Begin();
229  while (!nl_inter_a.End()) {
230  Int_t next_aa = nl_inter_a.Next();
231  Int_t nn = next_aa - zz;
232 
233  Float_t yield1 = fSystemList_[system_name1].yields[zz][next_aa].yield;
234  Float_t yield2 = fSystemList_[system_name2].yields[zz][next_aa].yield;
235 
236  Float_t yielderr1 = fSystemList_[system_name1].yields[zz][next_aa].yield_err;
237  Float_t yielderr2 = fSystemList_[system_name2].yields[zz][next_aa].yield_err;;
238 
239  Float_t ratio = yield2 / yield1;
240  Float_t ratio_err = TMath::Sqrt(TMath::Power(yielderr2 / yield1, 2.) + TMath::Power(yielderr1 * yield2 / yield1 / yield1, 2.));
241 
242  // - Fill the graph of the ratios with all points -
243  gr_ratio_all->SetPoint(npoint_all, nn, ratio);
244  gr_ratio_all->SetPointError(npoint_all, 0, ratio_err);
245  npoint_all++;
246 
247  // - Apply tolerance cut for the gaussian approximation -
248  Double_t diff1 = TMath::Abs(gaus1->GetParameter(1) - nn);
249  Double_t diff2 = TMath::Abs(gaus2->GetParameter(1) - nn);
250 
251  if (diff1 < tol * gaus1->GetParameter(2)) {
252  gr_sys1_tol->SetPoint(np1, nn, yield1);
253  gr_sys1_tol->SetPointError(np1, 0., yielderr1);
254  np1++;
255  }
256 
257  if (diff2 < tol * gaus2->GetParameter(2)) {
258  gr_sys2_tol->SetPoint(np2, nn, yield2);
259  gr_sys2_tol->SetPointError(np2, 0., yielderr2);
260  np2++;
261  }
262 
263  if ((diff1 < tol * gaus1->GetParameter(2)) && (diff2 < tol * gaus2->GetParameter(2))) {
264  gr_ratio_tol->SetPoint(npoint_tol, nn, ratio);
265  gr_ratio_tol->SetPointError(npoint_tol, 0, ratio_err);
266  npoint_tol++;
267  }
268  }//end shared A iter
269 
270  // --- Build ratio of the gaussian fits ---
271  TF1* gaus_ratio = new TF1(Form("gaus_ratio_Z%d_%ss_%s", zz, system_name1.c_str(), system_name2.c_str()), "gaus(0)/(gaus(3))", 0, 60);
272  gaus_ratio->SetLineColor(kBlack);
273  gaus_ratio->SetParameters(gaus2->GetParameter(0), gaus2->GetParameter(1), gaus2->GetParameter(2), gaus1->GetParameter(0), gaus1->GetParameter(1), gaus1->GetParameter(2));
274 
275  // --- Set some draw opt ---
276  gr_sys1_all->SetMarkerColor(kRed);
277  gr_sys1_all->SetMarkerStyle(24);
278  gr_sys2_all->SetMarkerColor(kBlue);
279  gr_sys2_all->SetMarkerStyle(26);
280  gr_ratio_all->SetMarkerColor(kBlack);
281  gr_ratio_all->SetMarkerStyle(25);
282 
283  gr_sys1_tol->SetMarkerColor(kRed);
284  gr_sys1_tol->SetMarkerStyle(20);
285  gr_sys2_tol->SetMarkerColor(kBlue);
286  gr_sys2_tol->SetMarkerStyle(22);
287  gr_ratio_tol->SetMarkerColor(kBlack);
288  gr_ratio_tol->SetMarkerStyle(21);
289 
290  gaus1->SetLineColor(kRed);
291  gaus2->SetLineColor(kBlue);
292  gaus_ratio->SetLineColor(kBlack);
293 
294  // --- Draw ---
295  gr_ratio_all->Draw("ap");
296  //gPad->SetLogy();
297  gr_sys1_all->Draw("p");
298  gr_sys2_all->Draw("p");
299  gaus1->Draw("same");
300  gaus2->Draw("same");
301  gaus_ratio->Draw("same");
302  gr_sys1_tol->Draw("p");
303  gr_sys2_tol->Draw("p");
304  gr_ratio_tol->Draw("p");
305 }
306 
307 
308 
318 
319 void KVIsoscaling::BuildLnR21vsNPlots(const std::string& system_name1, const std::string& system_name2)
320 {
321  // Method to compute the ratio \f$ln(R_{21})\f$ of the yields \f$Y_{i}(A,Z)\f$ of 2 asymetric reactions \f$(2)\f$ and \f$(1)\f$
322  // (projectiles and targets of reactions \f$(2)\f$ and \f$(1)\f$ only differ in their respective number of neutrons).
323  //
324  // We are interested in extracting the alpha isoscaling parameter, thus we want to draw \f$ln(R_{21})\f$ as a function of \f$N\f$ for a given \f$Z\f$.
325  // The tolerance parameter \f$tol\f$ is used to limit the number of points used to build the plots (see GetRMSTolerance() method).
326  //
327  // [in] system_name1 name of the n-deficient system
328  // [in] system_name2 name of the n-rich system index
329 
330  if (fdebug_) Info("BuildLnR21vsNPlots", "====== STARTING BUILDING PLOTS FOR %s/%s ======", system_name1.c_str(), system_name2.c_str());
331 
332  // --- Check the systems ---
333  if(system_name1!=system_name2){
334  std::string my_str = system_name1+"_"+system_name2;
335 
336  // --- Create a vector of the ratios for each Z (element) available ---
337  KVNumberList nl_inter = GetSharedZNumberList(system_name1, system_name2);
338  if (fdebug_) Info("BuildLnR21vsNPlots", "%s/%s - list of shared Z ='%s'", system_name1.c_str(), system_name2.c_str(), nl_inter.GetList());
339 
340  // --- Find the corresponding masses, yields and yields errors ---
341  Int_t ip = 0;
342  nl_inter.Begin();
343  while (!nl_inter.End()) {
344  Int_t next_zz = nl_inter.Next();
345 
346  // - Create the associated graph for ratios -
347  TGraphErrors* gr_all = new TGraphErrors();
348  gr_all->SetName(Form("gr_lnR21_vs_N_Z%d_%s_%s_all", next_zz, system_name1.c_str(), system_name2.c_str()));
349  gr_all->GetXaxis()->SetTitle("N");
350  gr_all->GetYaxis()->SetTitle("ln(R21)");
351  gr_all->SetMarkerStyle(20);
352  gr_all->SetMarkerColor(next_zz % 9 + 1);
353 
354  TGraphErrors* gr = new TGraphErrors();
355  gr->SetName(Form("gr_lnR21_vs_N_Z%d_%s_%s", next_zz, system_name1.c_str(), system_name2.c_str()));
356  gr->GetXaxis()->SetTitle("N");
357  gr->GetYaxis()->SetTitle("ln(R21)");
358  gr->SetMarkerStyle(20);
359  gr->SetMarkerColor(next_zz % 9 + 1);
360 
361  if(fdebug_) Info("BuildLnR21vsNPlots", "Graph %s created", gr->GetName());
362 
363  // - Find individual fits of the yields -
364  TF1* gaus_np = fSystemList_[system_name1].gauss[next_zz];
365  TF1* gaus_nr = fSystemList_[system_name2].gauss[next_zz];
366  Double_t par0 = gaus_nr->GetParameter(0);
367  Double_t par0_err = gaus_nr->GetParError(0);
368  Double_t par1 = gaus_nr->GetParameter(1);
369  Double_t par1_err = gaus_nr->GetParError(1);
370  Double_t par2 = gaus_nr->GetParameter(2);
371  Double_t par2_err = gaus_nr->GetParError(2);
372  Double_t par3 = gaus_np->GetParameter(0);
373  Double_t par3_err = gaus_np->GetParError(0);
374  Double_t par4 = gaus_np->GetParameter(1);
375  Double_t par4_err = gaus_np->GetParError(1);
376  Double_t par5 = gaus_np->GetParameter(2);
377  Double_t par5_err = gaus_np->GetParError(2);
378 
379  // - Compute gauss2/gauss1 function -
380  TF1* gaus_ratio = new TF1(Form("gaus_ratio_Z%d_%s_%s", next_zz, system_name1.c_str(), system_name2.c_str()), "gaus(0)/(gaus(3))", 0, 60);
381  gaus_ratio->SetLineColor(kBlack);
382  gaus_ratio->SetParameters(par0, par1, par2, par3, par4, par5);
383  fIsoscalingList_[my_str].gauss_ratios[next_zz] = gaus_ratio;
384 
385  TF1* gaus_ratio_ln = new TF1(Form("gaus_ratio_ln_Z%d_%s_%s", next_zz, system_name1.c_str(), system_name2.c_str()), "TMath::Log(gaus(0)/(gaus(3)))", 0, 60);
386  gaus_ratio_ln->SetLineColor(kBlack);
387  gaus_ratio_ln->SetParameters(par0, par1, par2, par3, par4, par5);
388  fIsoscalingList_[my_str].gauss_ratiosln[next_zz] = gaus_ratio_ln;
389 
390  // - Compute ln(R21) -
391  // - Extract all A shared between the 2 systems -
392  KVNumberList nl_inter_a = GetSharedANumberList(system_name1.c_str(), system_name2.c_str(), next_zz);
393 
394  if (fdebug_) Info("BuildLnR21vsNPlots", "=== Starting iteration over shared A list='%s' for Z=%d ===", nl_inter_a.GetList(), next_zz);
395 
396  Int_t npp = 0;
397  Int_t npp_all = 0;
398  nl_inter_a.Begin();
399  while (!nl_inter_a.End()) {
400  Int_t next_aa = nl_inter_a.Next();
401  Int_t nn = next_aa - next_zz;
402 
403  Float_t yield1 = fSystemList_[system_name1].yields[next_zz][next_aa].yield;
404  Float_t yield2 = fSystemList_[system_name2].yields[next_zz][next_aa].yield;
405 
406  Float_t lnr21 = TMath::Log(yield2 / yield1);
407 
408  Float_t yielderr1 = fSystemList_[system_name1].yields[next_zz][next_aa].yield_err;
409  Float_t yielderr2 = fSystemList_[system_name1].yields[next_zz][next_aa].yield_err;
410  Float_t err = TMath::Sqrt(TMath::Power(yielderr2 / yield2, 2.) + TMath::Power(yielderr1 / yield1, 2.));
411  //Float_t err = TMath::Abs(yielderr2/yield2) + TMath::Abs(yielderr1/yield1);
412 
413  // --- Fill the graph with all points ---
414  gr_all->SetPoint(npp_all, nn, lnr21);
415  gr_all->SetPointError(npp_all, 0, err);
416  npp_all++;
417 
418  // --- Apply tolerance cut for the gaussian approximation ---
419  Double_t diff1 = TMath::Abs(par1 - nn);
420  Double_t diff2 = TMath::Abs(par4 - nn);
421  if (diff1 < ftol_ * par2 && diff2 < ftol_ * par5) {
422  gr->SetPoint(npp, nn, lnr21);
423  gr->SetPointError(npp, 0, err);
424  npp++;
425  }
426  else if (fdebug_) {
427  printf("[Z=%d, N=%d, A=%d] Point refused : (tol1=%lf, diff1=%lf), (tol2=%lf, diff2=%lf) - (mean1=%lf, sigma1=%lf), (mean2=%lf, sigma2=%lf)",
428  next_zz, next_aa - next_zz, next_aa, ftol_ * par2, diff1, ftol_ * par5, diff2, par1, par2, par4, par5);
429  }
430  if (fdebug_) Info("BuildLnR21vsNPlots", "[Z=%d, N=%d, A=%d] Yield1(%s)=%lf, Yield2(%s)=%lf", next_zz, next_aa - next_zz, next_aa, system_name1.c_str(), yield1, system_name2.c_str(), yield2);
431  }
432 
433  ip++;
434 
435  // --- Save graph for the element ---
436  fIsoscalingList_[my_str].graphs[next_zz] = gr;
437 
438  }//end z list iter
439 
440  if (fdebug_) Info("BuildLnR21vsNPlots", "====== END OF BUILDING PLOTS FOR %s/%s ======", system_name1.c_str(), system_name2.c_str());
441 
442  }//end of systems OK
443  else {
444  Error("BuildLnR21Plots", "!!! Something is wrong in the input systems !!!");
445  }
446 
447  return;
448 }
449 
450 
451 
459 
460 void KVIsoscaling::FitLnR21vsNPlots(const std::string& system_name1, const std::string& system_name2, Option_t* option, Option_t* gooption)
461 {
462  // This method applies a linear fit to the \f$ln(R_{21})\f$ ratios to extract \f$\alpha\f$ isoscaling parameters (for each element shared between the two systems).
463  //
464  // \param[in] system_name1 name of the n-deficient system
465  // \param[in] system_name2 name of the n-rich system
466  // \param[in] option fit options
467  // \param[in] gooption graphic fit options
468 
469  if (fdebug_) Info("FitLnR21vsNPlots", "====== STARTING FITS FOR %s/%s ======", system_name1.c_str(), system_name2.c_str());
470 
471  // --- Init ---
472  TGraphErrors* gr = nullptr;
473  TF1* func = nullptr;
474 
475  std::string my_str = system_name1+"_"+system_name2;
476  IsoGraph_t my_isograph = fIsoscalingList_[my_str].graphs;
477 
478  // --- Iterate over all graphs and apply the fit ---
479  for(auto it = my_isograph.begin(); it != my_isograph.end(); it++){
480  Int_t zz = it->first;
481  gr = it->second;
482 
483  KVNumberList nl = GetSharedANumberList(system_name1, system_name2, zz);
484 
485  Int_t nmin = nl.First() - zz - 2;
486  Int_t nmax = nl.Last() - zz + 2;
487  func = new TF1(Form("fit_lnR21_vs_N_Z%d_%s_%s", zz, system_name1.c_str(), system_name2.c_str()), "pol1", float(nmin), float(nmax));
488  func->SetLineColor(zz % 9 + 1);
489 
490  for (Int_t ii = 0; ii < 10; ii++) gr->Fit(func, option, gooption, float(nmin + 2), float(nmax - 2));
491 
492  // --- Save fit results ---
493  fIsoscalingList_[my_str].fits[zz] = func;
494 
495  }//end graph iteration
496 
497  if (fdebug_) Info("BuildLnR21vsNPlots", "====== END OF FITS FOR %s/%s ======", system_name1.c_str(), system_name2.c_str());
498 }
499 
500 
501 
507 
508 void KVIsoscaling::DrawLnR21vsNFits(const std::string& system_name1, const std::string& system_name2)
509 {
510  // This method draws the results of \f$ln(R_{21})\f$ vs \f$N\f$ linear fits for the given system pair.
511  //
512  // \param[in] system_name1 name of the n-deficient system
513  // \param[in] system_name2 name of the n-rich system
514 
515  std::string my_str = system_name1+"_"+system_name2;
516 
517  TMultiGraph* mgr = new TMultiGraph();
518  TGraphErrors* gr = nullptr;
519  TF1* func = nullptr;
520 
521  // --- Iterate over all graphs ---
522  IsoGraph_t my_isograph = fIsoscalingList_[my_str].graphs;
523  for(auto it = my_isograph.begin(); it != my_isograph.end(); it++){
524  Int_t zz = it->first;
525  gr = it->second;
526 
527  mgr->Add(gr);
528  }//end graph iter
529 
530  mgr->Draw("ap");
531  mgr->GetXaxis()->SetTitle("N");
532  mgr->GetYaxis()->SetTitle("ln(R21)");
533 
534  // --- Iterate over all fits ---
535  IsoFit_t my_isofit = fIsoscalingList_[my_str].fits;
536  for(auto itt = my_isofit.begin(); itt != my_isofit.end(); itt++){
537  Int_t zz = itt->first;
538  func = itt->second;
539 
540  func->Draw("same");
541  }//end fit iter
542 }
543 
544 
545 
569 
570 void KVIsoscaling::BuildIsoscalingPlots(const std::string& system_name1, const std::string& system_name2, Int_t mcolor, Int_t mstyle, Bool_t draw)
571 {
572  // This method allows to create isoscaling results plots for a given system pair :
573  //
574  // \f$\alpha\f$ vs \f$Z\f$
575  // \f$\Delta\f$ vs \f$Z\f$
576  // \f$\alpha\f$ vs \f$\Delta\f$
577  // \f$Csym/T\f$ vs \f$Z\f$
578  //
579  // To compute \f$C_{sym}/T\f$, we use the usual isoscaling formula (see for ex. Raduta & Gulminelli, PRC 75 044605 bis (2007), Eq.(16)) :
580  //\f[
581  // C_{sym}(\langle A(Z) \rangle)/T = \alpha(Z)/\Delta(Z)
582  //\f]
583  //
584  // where
585  //\f[
586  //\Delta = (Z/\langle A_{1}(Z) \rangle)^{2}-(Z/\langle A_{2}(Z) \rangle)^{2}
587  //\f]
588  //
589  // \param[in] system_name1 name of the n-deficient system
590  // \param[in] system_name2 name of the n-rich system
591  // \param[in] mcolor color to use for drawing
592  // \param[in] mstyle style to use for drawing
593  // \param[in] draw =kTRUE if we want to plot the result
594 
595  if (fdebug_) Info("BuildIsoscalingPlots", "====== STARTING BUILDING PLOTS FOR %s/%s ======", system_name1.c_str(), system_name2.c_str());
596  std::string my_str = system_name1+"_"+system_name2;
597 
598  // --- Create the graphs ---
599  TGraphErrors* gr_alpha = new TGraphErrors();
600  gr_alpha->SetName(Form("Alpha_vs_Z_%s", my_str.c_str()));
601  gr_alpha->SetMarkerStyle(20);
602  gr_alpha->GetXaxis()->SetTitle("Z");
603  gr_alpha->GetYaxis()->SetTitle("#alpha");
604  gr_alpha->SetMarkerStyle(mstyle);
605  gr_alpha->SetMarkerColor(mcolor);
606 
607  TGraphErrors* gr_delta = new TGraphErrors();
608  gr_delta->SetName(Form("Delta_vs_Z_%s", my_str.c_str()));
609  gr_delta->SetMarkerStyle(20);
610  gr_delta->GetXaxis()->SetTitle("Z");
611  gr_delta->GetYaxis()->SetTitle("#Delta");
612  gr_delta->SetMarkerStyle(mstyle);
613  gr_delta->SetMarkerColor(mcolor);
614 
615  TGraphErrors* gr_alpha_delta = new TGraphErrors();
616  gr_alpha_delta->SetName(Form("Alpha_vs_Delta_%s", my_str.c_str()));
617  gr_alpha_delta->SetMarkerStyle(20);
618  gr_alpha_delta->GetYaxis()->SetTitle("#alpha");
619  gr_alpha_delta->GetXaxis()->SetTitle("#Delta");
620  gr_alpha_delta->SetMarkerStyle(mstyle);
621  gr_alpha_delta->SetMarkerColor(mcolor);
622 
623  TGraphErrors* gr_csymT = new TGraphErrors();
624  gr_csymT->SetName(Form("CsymOverT_vs_Z_%s", my_str.c_str()));
625  gr_csymT->SetMarkerStyle(mstyle);
626  gr_csymT->SetMarkerColor(mcolor);
627  gr_csymT->GetXaxis()->SetTitle("Z");
628  gr_csymT->GetYaxis()->SetTitle("#alpha/4#Delta");
629 
630  // --- Iterate over shared elements ---
631  KVNumberList nl = GetSharedZNumberList(system_name1, system_name2);
632 
633  Int_t np = 0;
634  nl.Begin();
635  while (!nl.End()) {
636  Int_t zz = nl.Next();
637 
638  Float_t alpha, alpha_err, delta, delta_err, csymT, csymT_err;
639  Bool_t is_ok_alpha = GetAlpha(system_name1, system_name2, zz, alpha, alpha_err);
640  Bool_t is_ok_delta = GetDeltaZA2(system_name1, system_name2, zz, delta, delta_err, fdebug_);
641  Bool_t is_ok_csymT = GetCsymOverT(system_name1, system_name2, zz, csymT, csymT_err, fdebug_);
642 
643  if (is_ok_delta && is_ok_alpha && is_ok_csymT) {
644 
645  gr_alpha->SetPoint(np, zz, alpha);
646  gr_alpha->SetPointError(np, 0., alpha_err);
647 
648  gr_delta->SetPoint(np, zz, delta);
649  gr_delta->SetPointError(np, 0., delta_err);
650 
651  gr_alpha_delta->SetPoint(np, delta, alpha);
652  gr_alpha_delta->SetPointError(np, delta_err, alpha_err);
653 
654  gr_csymT->SetPoint(np, zz, csymT);
655  gr_csymT->SetPointError(np, 0., csymT_err);
656 
657  np++;
658  }//end ok
659  }//end element iteration
660 
661  // --- Save results ---
662  fIsoscalingList_[my_str].gr_alpha = gr_alpha;
663  fIsoscalingList_[my_str].gr_delta = gr_delta;
664  fIsoscalingList_[my_str].gr_alpha_delta = gr_alpha_delta;
665  fIsoscalingList_[my_str].gr_csymT = gr_csymT;
666 
667  // --- Draw ---
668  //if (draw) gr->Draw("ap");
669 
670  if (fdebug_) Info("BuildIsoscalingPlots", "====== END BUILDING PLOTS FOR %s/%s ======", system_name1.c_str(), system_name2.c_str());
671 }
672 
673 
674 
679 
681 {
682  // This method creates a TMultiGraph containing all the \f$C_{sym}/T\f$ isoscaling results (for all combination of systems pairs already tested).
683  //
684  // \param[in] mgr pointer to the TMultiGraph provided by the user
685  TGraphErrors* gr = nullptr;
686 
687  // --- Iterate over all systems ---
688  for(auto it = fIsoscalingList_.begin(); it != fIsoscalingList_.end(); it++){
689  std::string str_sys = it->first;
690  Isoscaling_t my_iso = it->second;
691 
692  gr = my_iso.gr_csymT;
693  mgr->Add(gr);
694 
695  Info("CreateCsymOverTMultiGraph", "%s added to the MultiGraph", gr->GetName());
696  }
697 }
698 
699 
700 
705 
706 void KVIsoscaling::SaveResultsROOT(const Char_t* file_name)
707 {
708  // This method allows to save the yields and the isoscaling fits and graphs (\f$ln(R_{21})\f$, \f$\alpha\f$, \f$\Delta\f$, \f$C_{sym}/T\f$ graphs) in a *.root file.
709  //
710  // \param[in] file_name output file path
711 
712  TFile my_file(Form("%s", file_name), "RECREATE");
713 
714  // --- Save yields ---
715  TGraphErrors *gr_yield = nullptr;
716  TGraphErrors *gr_lnyield = nullptr;
717  TF1 *gr_gauss = nullptr;
718 
719  for(auto it0 = fSystemList_.begin(); it0 != fSystemList_.end(); it0++){
720  std::string str_sys = it0->first;
721  System_t my_system = it0->second;
722 
723  // - yields -
724  YieldGraph_t my_graphs = my_system.graphs;
725  for(auto itt0 = my_graphs.begin(); itt0 != my_graphs.end(); itt0++){
726  Int_t zz = itt0->first;
727  gr_yield = itt0->second;
728 
729  gr_yield->Write();
730  }//end yields graph iteration
731 
732  // - ln(yield) -
733  YieldGraph_t my_lngraphs = my_system.lngraphs;
734  for(auto ittt0 = my_lngraphs.begin(); ittt0 != my_lngraphs.end(); ittt0++){
735  Int_t zz = ittt0->first;
736  gr_lnyield = ittt0->second;
737 
738  gr_lnyield->Write();
739  }//end Ln(yield) graphs iteration
740 
741  // - gaussian fits -
742  YieldGauss_t my_gauss = my_system.gauss;
743  for(auto itttt0 = my_gauss.begin(); itttt0 != my_gauss.end(); itttt0++){
744  Int_t zz = itttt0->first;
745  gr_gauss= itttt0->second;
746 
747  gr_gauss->Write();
748  }//end fits iteration
749  }//end yield iteration
750 
751  // --- Save isoscaling results ---
752  TGraphErrors *gr_lnR21 = nullptr;
753  TF1 *fit_lnR21= nullptr;
754  TGraphErrors *gr_alpha = nullptr;
755  TGraphErrors *gr_delta = nullptr;
756  TGraphErrors *gr_alpha_delta = nullptr;
757  TGraphErrors *gr_csymT = nullptr;
758 
759  for(auto it = fIsoscalingList_.begin(); it != fIsoscalingList_.end(); it++){
760  std::string str_sys = it->first;
761  Isoscaling_t my_iso = it->second;
762 
763  IsoGraph_t my_graphs = my_iso.graphs;
764  for(auto itt = my_graphs.begin(); itt != my_graphs.end(); itt++){
765  Int_t zz = itt->first;
766  gr_lnR21 = itt->second;
767 
768  gr_lnR21->Write();
769  }//end lnR21 graph iteration
770 
771  IsoFit_t my_fits = my_iso.fits;
772  for(auto ittt = my_fits.begin(); ittt != my_fits.end(); ittt++){
773  Int_t zz = ittt->first;
774  fit_lnR21 = ittt->second;
775 
776  fit_lnR21->Write();
777  }//end lnR21 fits iteration
778 
779  gr_alpha = my_iso.gr_alpha;
780  gr_delta = my_iso.gr_delta;
781  gr_alpha_delta = my_iso.gr_alpha_delta;
782  gr_csymT = my_iso.gr_csymT;
783 
784  gr_alpha->Write();
785  gr_delta->Write();
786  gr_csymT->Write();
787  gr_alpha_delta->Write();
788  }//end isoscaling list iter
789 
790  if (fdebug_) Info("SaveResultsROOT", "File '%s' saved", file_name);
791 }
792 
793 
794 
799 
801 {
802  // This method allows to save results in a *.txt file.
803  //
804  // \param[in] file_name output file path
805  KVString fname(file_name);
806 
807  std::ofstream my_file;
808  my_file.open(fname.Data(), std::ofstream::out);
809  my_file << "### system_pair Z Alpha d(Alpha) <A_1> d<A_1> <A_2> d<A_2> denum d(denum) Csym/T d(Csym/T) ###" << std::endl;
810 
811  // --- Iterate over all isoscaling results ---
812  for(auto it0 = fIsoscalingList_.begin(); it0 != fIsoscalingList_.end(); it0++){
813  std::string str_sys = it0->first;
814  Isoscaling_t my_iso = it0->second;
815 
816  // !!! only works if system_name does not contain _ !!!
817  KVString mystr(str_sys);
818  TObjArray *arr = mystr.Tokenize("_");
819  std::string sys1(((TObjString*) arr->At(0))->GetString().Data());
820  std::string sys2(((TObjString*) arr->At(1))->GetString().Data());
821 
822  my_file << "### " << sys1 << " " << sys2 << " ###" << std::endl;
823 
824  IsoFit_t my_isofit = my_iso.fits;
825 
826  // --- Iterate over all elements for the given isoscaling fit result ---
827  for(auto it1 = my_isofit.begin(); it1 != my_isofit.end(); it1++){
828  Int_t zz = it1->first;
829  TF1 *fit = it1->second;
830 
831  // Alpha
832  Float_t alpha = -666.;
833  Float_t alpha_err = -666.;
834  GetAlpha(sys1, sys2, zz, alpha, alpha_err);
835 
836  // Denum
837  Float_t denum = -666.;
838  Float_t denum_err = -666.;
839  GetDeltaZA2(sys1, sys2, zz, denum, denum_err, kFALSE);
840 
841  //Csym/T
842  Float_t csymT = -666.;
843  Float_t csymT_err = -666.;
844  GetCsymOverT(sys1, sys2, zz, csymT, csymT_err, kFALSE);
845 
846  //<A1> and <A2>
847  Float_t meanA1 = -666.;
848  Float_t meanA1_err = -666.;
849  Float_t meanA2 = -666.;
850  Float_t meanA2_err = -666.;
851  GetAMean(sys1, zz, meanA1, meanA1_err);
852  GetAMean(sys2, zz, meanA2, meanA2_err);
853 
854 
855  KVString line(Form("%s %s %d %lf %lf %lf %lf %lf %lf %lf %lf %lf %lf", sys1.c_str(), sys2.c_str(), zz, alpha, alpha_err, meanA1, meanA1_err, meanA2, meanA2_err, denum, denum_err, csymT, csymT_err));
856  my_file << line.Data() << std::endl;
857  }//end element iteration
858 
859  if (fdebug_) Info("SaveResultsASCII", "[%s/%s] results saved\n", sys1.c_str(), sys2.c_str());
860 
861  }//end isoscaling results iteration
862 
863  my_file.close();
864 
865  if (fdebug_) Info("SaveResultsASCII", "File '%s' saved", fname.Data());
866 }
867 
868 
869 
870 //===================================
871 //== GETTERS YIELDS ==
872 //===================================
873 
882 
883 Bool_t KVIsoscaling::GetAMean(const std::string& system_name, Int_t zz, Float_t& meanA, Float_t& meanA_err)
884 {
885  // This method returns the average mass \f$ \langle A(Z) \rangle \f$ of the given element (charge \f$Z\f$) for the given system.
886  // The average mass and associated uncertainty are obtained from the result of a gaussian fit applied to the isotopic distribution.
887  //
888  // \param[in] system_name name of the system
889  // \param[in] zz charge of the element
890  // \param[out] meanA average mass
891  // \param[out] meanA_err error on the average mass
892  Bool_t is_ok = kTRUE;
893 
894  // find the corresponding graph
895  KVNumberList nl = GetZNumberList(system_name);
896 
897  if(nl.Contains(zz)){
898  TF1* my_gaus = fSystemList_[system_name].gauss[zz];
899  meanA = my_gaus->GetParameter(1);
900  meanA_err= my_gaus->GetParError(1);
901  }
902  else{
903  Error("GetAMean", "[sys=%s] Z=%d is out of range", system_name.c_str(), zz);
904  is_ok = kFALSE;
905  }
906 
907  return is_ok;
908 }
909 
910 
911 
917 
918 Int_t KVIsoscaling::FindZFromAmean(const std::string& system_name, Int_t aa)
919 {
920  // Returns the charge \f$Z\f$ (element) associated to the provided mass \f$A\f$ such as \f$\langle A(Z) \rangle = A \f$.
921  //
922  // \param[in] system_name name of the system
923  // \param[in] aa mass
924 
925  Int_t zz = -666;
926  Double_t min_val = 10.;
927 
928  KVNumberList nl = GetZNumberList(system_name);
929  nl.Begin();
930  while (!nl.End()) {
931  Int_t next_zz = nl.Next();
932 
933  // extract mean aa
934  Float_t aa_real, aa_real_err;
935  GetAMean(system_name, next_zz, aa_real, aa_real_err);
936 
937  Double_t delta_aa = TMath::Abs(aa - aa_real);
938  if (delta_aa < min_val) {
939  zz = next_zz;
940  min_val = delta_aa;
941 
942  if (fdebug_) Info("FindZFromAmean", "[%s - A=%d] : Z=%d -> <A(Z)>=%lf +/- %lf [FOUND] delta=%lf, min_val=%lf", system_name.c_str(), aa, next_zz, aa_real, aa_real_err, delta_aa, min_val);
943  }
944  else if (fdebug_) Info("FindZFromAmean", "[%s - A=%d] : Z=%d -> <A(Z)>=%lf +/- %lf [NOT FOUND] delta=%lf, min_val=%lf", system_name.c_str(), aa, next_zz, aa_real, aa_real_err, delta_aa, min_val);
945  }
946 
947  return zz;
948 }
949 
950 
951 
952 //=======================================
953 //== GETTERS ISOSCALING ==
954 //=======================================
955 
964 
965 Bool_t KVIsoscaling::GetAlpha(const std::string& system_name1, const std::string& system_name2, Int_t zz, Float_t& alpha, Float_t& alpha_err)
966 {
967  // This method extracts the isoscaling \f$\alpha\f$ parameter.
968  //
969  // \param[in] system_name1 name of the n-deficient system
970  // \param[in] system_name2 name of the n-rich system
971  // \param[in] zz charge for the element the isoscaling fit was applied
972  // \param[out] alpha isoscaling fit parameter \f$\alpha\f$
973  // \param[out] alpha_err error on the isoscaling fit parameter \f$\alpha\f$
974 
975  Bool_t is_ok = kTRUE;
976 
977  // --- Find the corresponding alpha parameter in fit function list ---
978  std::string my_str = system_name1+"_"+system_name2;
979  TF1* func = fIsoscalingList_[my_str].fits[zz];
980 
981  if(func!=nullptr){
982  alpha = func->GetParameter(1);
983  alpha_err = func->GetParError(1);
984  }
985  else{
986  Error("GetAlpha", "!!! (%s / %s) [Z=%d] : Fit function not found !!!", system_name1.c_str(), system_name2.c_str(), zz);
987  is_ok=kFALSE;
988  }
989 
990  return is_ok;
991 }
992 
993 
994 
1006 
1007 Bool_t KVIsoscaling::GetDeltaZA2(const std::string& system_name1, const std::string& system_name2, Int_t zz, Float_t& denum, Float_t& denum_err, Bool_t debug)
1008 {
1009  // This method returns the value of the asymmetry parameter \f$ \Delta \f$ (denumerator for \f$ C_{sym}/T \f$ estimation), such as
1010  //\f[
1011  //\Delta = (Z/\langle A_{1} \rangle)^{2}-(Z/ \langle A_{2} \rangle)^{2}
1012  //\f]
1013  //
1014  // \param[in] system_name1 name of the n-deficient system
1015  // \param[in] system_name2 name of the n-rich system
1016  // \param[in] zz charge for the element the isoscaling fit was applied
1017  // \param[out] denum asymmetry parameter \f$\Delta\f$
1018  // \param[out] denum_err error on the asymmetry parameter \f$\Delta\f$
1019 
1020  Bool_t is_ok = kFALSE;
1021 
1022  Float_t meanA1 = -666.;
1023  Float_t meanA2 = -666.;
1024  Float_t meanA1_err = -666.;
1025  Float_t meanA2_err = -666.;
1026 
1027  Bool_t ok1 = GetAMean(system_name1, zz, meanA1, meanA1_err);
1028  Bool_t ok2 = GetAMean(system_name2, zz, meanA2, meanA2_err);
1029 
1030  if (ok1 && ok2) {
1031  is_ok = kTRUE;
1032 
1033  denum = TMath::Power(1.*zz / meanA1, 2.) - TMath::Power(1.*zz / meanA2, 2.);
1034  denum_err = TMath::Sqrt(TMath::Power(zz * zz * 2.*meanA1_err / meanA1 / meanA1 / meanA1, 2.) + TMath::Power(zz * zz * 2.*meanA2_err / meanA2 / meanA2 / meanA2, 2.));
1035 
1036  if (debug) Info("GetDeltaZA2", "(%s / %s) [Z=%d] : denum=%lf +/- %lf", system_name1.c_str(), system_name2.c_str(), zz, denum, denum_err);
1037  }
1038  else {
1039  if (!ok1) Warning("GetDeltaZA2", "<A1> was not found for Z=%d of system %s/%s", zz, system_name1.c_str(), system_name2.c_str());
1040  if (!ok2) Warning("GetDeltaZA2", "<A2> was not found for Z=%d of system %s/%s", zz, system_name1.c_str(), system_name2.c_str());
1041  }
1042 
1043  return is_ok;
1044 }
1045 
1046 
1047 
1059 
1060 Bool_t KVIsoscaling::GetCsymOverT(const std::string& system_name1, const std::string& system_name2, Int_t zz, Float_t& csymT, Float_t& csymT_err, Bool_t debug)
1061 {
1062  // This method returns the value of \f$C_{sym}/T\f$ from isoscaling method such as (see for ex. Raduta & Gulminelli, PRC 75 044605 (2007) Eq.(16))
1063  //\f[
1064  // C_{sym}/T = \alpha/4\Delta
1065  //\f]
1066  //
1067  // \param[in] system_name1 name of the n-deficient system
1068  // \param[in] system_name2 name of the n-rich system
1069  // \param[in] zz charge for the element the isoscaling fit was applied
1070  // \param[out] csymT value of the ratio \f$C_{sym}/T = \alpha/4\Delta\f$
1071  // \param[out] csymT_err error on the ratio \f$C_{sym}/T = \alpha/4\Delta\f$
1072 
1073  Bool_t is_ok = kFALSE;
1074  Bool_t is_oka = kFALSE;
1075  Bool_t is_okd = kFALSE;
1076 
1077  // --- Get corresponding alpha isoscaling parameter ---
1078  Float_t alpha = -666.;
1079  Float_t alpha_err = -666.;
1080  is_oka = GetAlpha(system_name1, system_name2, zz, alpha, alpha_err);
1081 
1082  // --- Get the denum ---
1083  Float_t denum = -666.;
1084  Float_t denum_err = -666.;
1085  is_okd = GetDeltaZA2(system_name1, system_name2, zz, denum, denum_err, debug);
1086 
1087  // --- Compute Csym/T ---
1088  if (is_oka && is_okd) {
1089  is_ok = kTRUE;
1090  csymT = alpha / 4. / denum;
1091  csymT_err = TMath::Sqrt(TMath::Power(alpha_err / 4. / denum, 2.) + TMath::Power(alpha * denum_err / 4. / denum / denum, 2.));
1092 
1093  if (debug) {
1094  Info("GetCsymOverT", "(%s / %s) [Z=%d] : alpha=%lf +/- %lf, denum=%lf +/- %lf, Csym/T=%lf +/- %lf", system_name1.c_str(), system_name2.c_str(), zz, alpha, alpha_err, denum, denum_err, csymT, csymT_err);
1095  }
1096  }
1097 
1098  return is_ok;
1099 }
1100 
1101 
1102 
1109 
1110 KVNumberList KVIsoscaling::GetZNumberList(const std::string& system_name)
1111 {
1112  // This method returns the list of elements of a given system
1113  //
1114  // \param[in] system_name name of the system
1115  //
1116  // \returns KVNumberList of the all the charges \f$Z\f$
1117 
1118  KVNumberList nl;
1119  YieldData_t my_yields = fSystemList_[system_name].yields;
1120 
1121  for(auto it = my_yields.begin(); it != my_yields.end(); it++){
1122  Int_t zz = it->first;
1123  nl.Add(zz);
1124  }
1125 
1126  return nl;
1127 }
1128 
1129 
1130 
1138 
1139 KVNumberList KVIsoscaling::GetSharedZNumberList(const std::string& system_name1, const std::string& system_name2)
1140 {
1141  // This method returns the list of all elements shared between the two provided systems.
1142  //
1143  // \param[in] system_name1 name of the first system
1144  // \param[in] system_name2 name of the second system
1145  //
1146  // \returns KVNumberList of the shared charges \f$Z\f$
1147 
1148  KVNumberList nl_z1 = GetZNumberList(system_name1);
1149  KVNumberList nl_z2 = GetZNumberList(system_name2);
1150  KVNumberList nl_inter;
1151  nl_z1.Copy(nl_inter);
1152  nl_inter.Inter(nl_z2);
1153 
1154  if (fdebug_) {
1155  Info("GetSharedZNumberList", "Zlist[%s]= %s, Zlist[%s]= %s, shared Zlist= %s",
1156  system_name1.c_str(), nl_z1.GetList(), system_name2.c_str(), nl_z2.GetList(), nl_inter.GetList()) ;
1157  }
1158 
1159  return nl_inter;
1160 }
1161 
1162 
1163 
1171 
1172 KVNumberList KVIsoscaling::GetANumberList(const std::string& system_name, Int_t zz)
1173 {
1174  // This method returns the list of all isotopes for the given element.
1175  //
1176  // \param[in] system_name name of the system
1177  // \param[in] zz charge of the element
1178  //
1179  // \returns KVNumberList of the masses \f$A\f$
1180 
1181  KVNumberList nl;
1182  Element_t my_element = fSystemList_[system_name].yields[zz];
1183 
1184  for(auto it = my_element.begin(); it != my_element.end(); it++){
1185  Int_t aa = it->first;
1186  nl.Add(aa);
1187  }
1188 
1189  return nl;
1190 }
1191 
1192 
1193 
1202 
1203 KVNumberList KVIsoscaling::GetSharedANumberList(const std::string& system_name1, const std::string& system_name2, Int_t zz)
1204 {
1205  // This method returns the list of all isotopes shared between the two provided systems for the given element.
1206  //
1207  // \param[in] system_name1 name of the first system
1208  // \param[in] system_name2 name of the second system
1209  // \param[in] zz charge of the element
1210  //
1211  // \returns KVNumberList of the shared masses \f$A\f$
1212 
1213  KVNumberList nl_a1 = GetANumberList(system_name1, zz);
1214  KVNumberList nl_a2 = GetANumberList(system_name2, zz);
1215  KVNumberList nl_inter_a;
1216  nl_a1.Copy(nl_inter_a);
1217  nl_inter_a.Inter(nl_a2);
1218 
1219  if (fdebug_){
1220  Info("GetSharedANumberList", "Z=%d : Alist[%s]= %s, Alist[%s]= %s, shared Alist= %s",
1221  zz, system_name1.c_str(), nl_a1.GetList(), system_name2.c_str(), nl_a2.GetList(), nl_inter_a.GetList()) ;
1222  }
1223 
1224  return nl_inter_a;
1225 }
1226 
1227 
1228 //============================
1229 //== PRINTERS ==
1230 //============================
1231 
1233 
1235 {
1236  Info("PrintYields", "Printing the list of yields : ");
1237 
1238  // --- Iterate over systems ---
1239  for(auto it0 = fSystemList_.begin(); it0 != fSystemList_.end(); it0++){
1240 
1241  std::string sys_name = it0->first;
1242  System_t my_system = it0->second;
1243 
1244  std::string my_path = my_system.file_path;
1245  YieldData_t my_yields = my_system.yields;
1246 
1247  printf("\n====== %s =====\n", sys_name.c_str());
1248  printf("extracted from path = '%s'\n", my_path.c_str());
1249 
1250  // --- Iterate over elements for the given system ---
1251  for(auto it1 = my_yields.begin(); it1 != my_yields.end(); it1++){
1252  Int_t zz = it1->first;
1253  Element_t my_ele = it1->second;
1254 
1255  // --- Iterate over isostopes of the given element ---
1256  for(auto it2 = my_ele.begin(); it2 != my_ele.end(); it2++){
1257  Int_t aa = it2->first;
1258  IsotopeYield_t yield = it2->second;
1259 
1260  printf(" Y(Z=%d, A=%d) = %lf +/- %lf\n", zz, aa, yield.yield, yield.yield_err);
1261  }//end isostope iteration
1262  }//end element iteration
1263  }//end system iteration
1264 
1265  Info("PrintYields", "That's all folks !");
1266 }
1267 
1268 
1269 
1271 
1273 {
1274  Info("PrintSystemsList", "Printing the list of available systems : ");
1275 
1276  // --- Iterate over systems ---
1277  for(auto it0 = fSystemList_.begin(); it0 != fSystemList_.end(); it0++){
1278 
1279  std::string sys_name = it0->first;
1280  System_t my_system = it0->second;
1281 
1282  std::string my_path = my_system.file_path;
1283  YieldData_t my_yields = my_system.yields;
1284 
1285  printf("\n====== %s =====\n", sys_name.c_str());
1286  printf("extracted from path = '%s'\n", my_path.c_str());
1287  }//end system iteration
1288 }
1289 
1290 
1291 
int Int_t
bool Bool_t
char Char_t
float Float_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
kRed
kBlack
kBlue
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
char * Form(const char *fmt,...)
Isoscaling class.
Definition: KVIsoscaling.h:198
KVNumberList GetANumberList(const std::string &system_name, Int_t zz)
std::unordered_map< std::string, System_t > fSystemList_
(hash) map by name
Definition: KVIsoscaling.h:254
void PrintYieldsList()
— Printers —
Bool_t GetAMean(const std::string &system_name, Int_t zz, Float_t &meanA, Float_t &meanA_err)
— Getters —
void BuildIsoscalingPlots(const std::string &system_name1, const std::string &system_name2, Int_t mcolor, Int_t mstyle, Bool_t draw=kFALSE)
Bool_t fdebug_
verbose mode for debugging
Definition: KVIsoscaling.h:257
void FitLnR21vsNPlots(const std::string &system_name1, const std::string &system_name2, Option_t *option="MNVR", Option_t *gooption="goff")
Double_t ftol_
tolerance for the gaussian approximation (in sigma)
Definition: KVIsoscaling.h:258
KVNumberList GetSharedZNumberList(const std::string &system_name1, const std::string &system_name2)
void PrintSystemsList()
void BuildLnR21vsNPlots(const std::string &system_name1, const std::string &system_name2)
void DrawLnR21vsNFits(const std::string &system_name1, const std::string &system_name2)
Bool_t GetCsymOverT(const std::string &system_name1, const std::string &system_name2, Int_t zz, Float_t &csymT, Float_t &csymT_err, Bool_t debug)
Bool_t GetAlpha(const std::string &system_name1, const std::string &system_name2, Int_t zz, Float_t &alpha, Float_t &alpha_err)
void TestGaussianApprox(const std::string &system_name1, const std::string &system_name2, Int_t zz, Double_t tol)
std::unordered_map< std::string, Isoscaling_t > fIsoscalingList_
(hash) map by name
Definition: KVIsoscaling.h:255
void BuildGaussianPlots(const std::string &system_name)
void SaveResultsROOT(const Char_t *file_name="./isoscaling_output_file.root")
KVNumberList GetSharedANumberList(const std::string &system_name1, const std::string &system_name2, Int_t zz)
void CreateCsymOverTMultiGraph(TMultiGraph *mgr)
KVNumberList GetZNumberList(const std::string &system_name)
void SaveResultsASCII(const Char_t *file_name="./isoscaling_output_file.txt")
Bool_t GetDeltaZA2(const std::string &system_name1, const std::string &system_name2, Int_t zz, Float_t &denum, Float_t &denum_err, Bool_t debug)
Int_t FindZFromAmean(const std::string &system_name, Int_t aa)
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:85
void Inter(const KVNumberList &list)
void Copy(TObject &) const
Copy content of this number list into 'o'.
const Char_t * GetList() const
Bool_t Contains(Int_t val) const
returns kTRUE if the value 'val' is contained in the ranges defined by the number list
Int_t First() const
Returns smallest number included in list.
Bool_t End(void) const
Definition: KVNumberList.h:199
void Begin(void) const
void Add(Int_t)
Add value 'n' to the list.
Int_t Last() const
Returns largest number included in list.
Int_t Next(void) const
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
virtual void SetLineColor(Color_t lcolor)
virtual void SetMarkerColor(Color_t mcolor=1)
virtual void SetMarkerStyle(Style_t mstyle=1)
virtual Double_t GetParameter(const TString &name) const
virtual Double_t GetParError(Int_t ipar) const
void Draw(Option_t *option="") override
virtual void SetParameters(const Double_t *params)
virtual void SetPointError(Double_t ex, Double_t ey)
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Axis_t xmin=0, Axis_t xmax=0)
void SetName(const char *name="") override
void Draw(Option_t *chopt="") override
TAxis * GetXaxis() const
TAxis * GetYaxis() const
virtual void Add(TGraph *graph, Option_t *chopt="")
void Draw(Option_t *chopt="") override
TAxis * GetYaxis()
TAxis * GetXaxis()
virtual void SetTitle(const char *title="")
const char * GetName() const override
TObject * At(Int_t idx) const override
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
const char * Data() const
TObjArray * Tokenize(const TString &delim) const
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
TLine * line
TGraphErrors * gr
fit(model, train_loader, val_loader, num_epochs, batch_size, optimizer, criterion, save_best, scheduler)
tuple file_name
void Error(const char *location, const char *fmt,...)
void Info(const char *location, const char *fmt,...)
Double_t Power(Double_t x, Double_t y)
Double_t Log(Double_t x)
Double_t Sqrt(Double_t x)
Double_t Abs(Double_t d)
const char * String
TGraphErrors * gr_alpha_delta
Definition: KVIsoscaling.h:193
TGraphErrors * gr_csymT
Definition: KVIsoscaling.h:194
TGraphErrors * gr_alpha
Definition: KVIsoscaling.h:191
IsoGraph_t graphs
Definition: KVIsoscaling.h:187
TGraphErrors * gr_delta
Definition: KVIsoscaling.h:192
IsoFit_t fits
Definition: KVIsoscaling.h:188
Float_t yield_err
Definition: KVIsoscaling.h:163
YieldGauss_t gauss
Definition: KVIsoscaling.h:174
YieldGraph_t lngraphs
Definition: KVIsoscaling.h:173
YieldData_t yields
Definition: KVIsoscaling.h:171
YieldGraph_t graphs
Definition: KVIsoscaling.h:172
std::string file_path
Definition: KVIsoscaling.h:175
ClassImp(TPyArg)