KaliVeda
Toolkit for HIC analysis
KVIDZAFromZGrid.cpp
1 //Created by KVClassFactory on Tue Mar 8 10:00:16 2016
2 //Author: Diego Gruyer
3 
4 #include "KVIDZAFromZGrid.h"
5 #include "TMultiGraph.h"
6 #include "KVIDZALine.h"
7 #include "TCanvas.h"
8 
12 
13 
14 
15 
20 {
21  // Default constructor
22  // Grid is declared as a 'ZOnlyGrid' by default (this is internal mechanics)
23 
24  init();
25  fFits.SetOwner();
26  SetOnlyZId();
27  fIgnoreMassID = false;
28 }
29 
30 
31 
32 
40 
42 {
43  // This method copies the current state of 'this' object into 'obj'
44  // You should add here any member variables, for example:
45  // (supposing a member variable KVIDZAFromZGrid::fToto)
46  // CastedObj.fToto = fToto;
47  // or
48  // CastedObj.SetToto( GetToto() );
49 
50  KVIDZAGrid::Copy(obj);
52  i.fZmaxInt = fZmaxInt;
53  i.fPIDRange = fPIDRange;
54  i.fZminInt = fZminInt;
55  fTables.Copy(i.fTables);
57 }
58 
59 
60 
61 
70 
71 void KVIDZAFromZGrid::ReadFromAsciiFile(std::ifstream& gridfile)
72 {
73 // fPIDRange = kFALSE;
74 // KVIDGraph::ReadFromAsciiFile(gridfile);
75 // if (GetParameters()->HasParameter("PIDRANGE")) {
76 // fPIDRange = kTRUE;
77 // fZmaxInt = GetParameters()->GetIntValue("PIDRANGE");
78 // LoadPIDRanges();
79 // }
80 
81  fPIDRange = kFALSE;
82 // fHasMassCut = kFALSE;
84 
85 // if (GetIdentifier("MassID")) fHasMassCut = kTRUE;
86 
87  if (GetParameters()->HasParameter("PIDRANGE")) {
88  fPIDRange = kTRUE;
89  TString pidrange = GetParameters()->GetStringValue("PIDRANGE");
90  if (pidrange.Contains("-")) {
91  TString min = (pidrange(0, pidrange.Index("-")));
92  fZminInt = min.Atoi();
93  min = (pidrange(pidrange.Index("-") + 1, pidrange.Length()));
94  fZmaxInt = min.Atoi();
95  }
96  else {
97  fZminInt = 1;
98  fZmaxInt = pidrange.Atoi();
99  }
100  LoadPIDRanges();
101  }
102  // if <PARAMETER> IgnoreMassID=1 appears in file, we are only using the PID intervals to clean
103  // up a messy de-e plot, not to give mass identification. particles will only be identified in Z.
104  if (GetParameters()->HasParameter("IgnoreMassID") && GetParameters()->GetIntValue("IgnoreMassID") == 1)
105  fIgnoreMassID = true;
106  else
107  fIgnoreMassID = false;
108 }
109 
110 
111 
113 
114 void KVIDZAFromZGrid::WriteToAsciiFile(std::ofstream& gridfile)
115 {
116  ExportToGrid();
117  KVIDGraph::WriteToAsciiFile(gridfile);
118 }
119 
120 
121 
123 
125 {
126  fZminInt = 100000;
127  fZmaxInt = 0;
128 
129  KVIDentifier* id = 0;
130  TIter it(GetIdentifiers());
131  while ((id = (KVIDentifier*)it())) {
132  int zz = id->GetZ();
133  if (!GetParameters()->HasParameter(Form("PIDRANGE%d", zz))) continue;
134  KVString mes = GetParameters()->GetStringValue(Form("PIDRANGE%d", zz));
135  if (mes.IsWhitespace()) continue;
136  int type = (mes.Contains(",") ? 2 : 1);
137  interval_set* itv = new interval_set(zz, type);
138  itv->SetName(GetName());
139  mes.Begin("|");
140  while (!mes.End()) {
141  KVString tmp = mes.Next();
142  tmp.Begin(":");
143  int aa = tmp.Next().Atoi();
144  KVString val = tmp.Next();
145  double pidmin, pidmax, pid;
146  if (type == 1) itv->add(aa, val.Atof());
147  else if (type == 2) {
148  val.Begin(",");
149  pidmin = val.Next().Atof();
150  pid = val.Next().Atof();
151  pidmax = val.Next().Atof();
152  itv->add(aa, pid, pidmin, pidmax);
153 // itv->add(aa, pid, pid-0.02, pid+0.02);
154  }
155  }
156  if (zz < fZminInt) fZminInt = zz;
157  if (zz > fZmaxInt) fZmaxInt = zz;
158  fTables.Add(itv);
159  }
160  fPIDRange = kTRUE;
161  // PrintPIDLimits();
162 }
163 
164 
165 
167 
169 {
170  fTables.Clear("all");
171  fPIDRange = kFALSE;
172 }
173 
174 
175 
177 
179 {
180  fTables.Clear("all");
181  LoadPIDRanges();
182 }
183 
184 
185 
187 
189 {
190  interval_set* itv = 0;
191  TIter it(&fTables);
192  while ((itv = (interval_set*)it())) if (itv->GetZ() == zint) return itv;
193  return 0;
194 }
195 
196 
197 
204 
206 {
207 // ((interval_set*)fTables.At(12))->fIntervals.ls();
208 
209 // for (int zz = fZminInt; zz <= fZmaxInt; zz++) {
210 // Info("PrintPIDLimits", "Z=%2d [%.4lf %.4lf]", zz, ((interval_set*)fTables.At(zz - fZminInt))->fPIDmins.at(0),
211 // ((interval_set*)fTables.At(zz - fZminInt))->fPIDmaxs.at(((interval_set*)fTables.At(zz - fZminInt))->fNPIDs - 1));
212 // }
213 }
214 
215 
216 
218 
220 {
221  if (GetParameters()->HasParameter("PIDRANGE")) GetParameters()->RemoveParameter("PIDRANGE");
222  for (int ii = 1; ii < 30; ii++) {
223  if (GetParameters()->HasParameter(Form("PIDRANGE%d", ii))) GetParameters()->RemoveParameter(Form("PIDRANGE%d", ii));
224  }
225 }
226 
227 
228 
235 
236 int KVIDZAFromZGrid::is_inside(double pid) const
237 {
238  // Look for a set of mass-interval definitions in which the given PID
239  // falls (PID from linearisation of Z identification).
240  // In principle this should be the set corresponding to Z=nint(PID),
241  // but if not Z+/-1 are also tried.
242  // Returns the value of Z for the set found (or 0 if no set found)
243 
244  int zint = TMath::Nint(pid);
245  interval_set* it = GetIntervalSet(zint);
246  if (it) {
247  if (it->is_inside(pid)) return zint;
248  else if (it->is_above(pid)) {
249 
250  it = GetIntervalSet(zint + 1);
251  if (it && it->is_inside(pid)) return zint + 1;
252  else return 0;
253  }
254  else {
255  it = GetIntervalSet(zint - 1);
256  if (it && it->is_inside(pid)) return zint - 1;
257  else return 0;
258  }
259  }
260  else return 0;
261 }
262 
263 
264 
273 
275 {
276  // General initialisation method for identification grid.
277  //
278  // This method MUST be called once before using the grid for identifications.
279  //
280  // + The ID lines are sorted.
281  // + The natural line widths of all ID lines are calculated.
282  // + The line with the largest Z (Zmax line) is found.
283 
284  SetOnlyZId();
286 
287  // set up mass fits (if any)
288  fFits.Clear();
289  if (GetParameters()->HasStringParameter("MASSFITS")) {
290  KVNumberList zlist(GetParameters()->GetStringValue("MASSFITS"));
291  for (auto z : zlist) {
292  auto massfit = GetParameters()->GetTStringValue(Form("MASSFIT_%d", z));
293  massfit.ReplaceAll(":", "=");
294  KVNameValueList fitparams;
295  fitparams.Set(massfit);
296  fFits.Add(new KVMultiGaussIsotopeFit(z, fitparams));
297  }
298  }
299 }
300 
301 
302 
304 
306 {
307  double P;
308  auto A = fitfunc->GetA(idr->PID, P);
309  if (A) {
310  idr->A = A;
311  idr->PID = fitfunc->GetInterpolatedA(idr->PID);
312  if (P > 0.5) idr->IDquality = KVIDZAGrid::kICODE0; // probability of A is >50%
313  else idr->IDquality = KVIDZAGrid::kICODE3;// OK, slight ambiguity of A
314  }
315  else {
316  // returned A=0 => background noise
318  return false;
319  }
320  return true;
321 }
322 
323 
324 
336 
338 {
339  // Fill the KVIdentificationResult object with the results of identification for point (x,y)
340  // corresponding to some physically measured quantities related to a reconstructed nucleus.
341  //
342  // If identification is successful, idr->IDOK = true.
343  // In this case, idr->Aident and idr->Zident indicate whether isotopic or only Z identification
344  // was acheived.
345  //
346  // In case of unsuccessful identification, idr->IDOK = false,
347  // BUT idr->Zident and/or idr->Aident may be true: this is to indicate which kind of
348  // identification was attempted but failed (this changes the meaning of the quality code)
349 
350  idr->Aident = idr->Zident = kFALSE;
351 
352  KVIDZAGrid::Identify(x, y, idr);
353  idr->Zident = kTRUE; // meaning Z identification was attempted, even if it failed
354  if (!idr->IDOK) return;
355 
356  bool have_pid_range_for_Z = fPIDRange && (idr->Z <= fZmaxInt) && (idr->Z > fZminInt - 1);
357  auto mass_fit_for_Z = GetMultiGaussFit(idr->Z);
358  bool have_mass_fit_for_Z = (mass_fit_for_Z != nullptr);
359  bool mass_id_success = false;
360 
361  if ((have_mass_fit_for_Z || have_pid_range_for_Z)
362  && (!fHasMassIDRegion || idr->HasFlag(GetName(), "MassID"))) { // if a mass ID region is defined, we must be inside it
363  // try mass identification
364  if (have_mass_fit_for_Z)
365  mass_id_success = MassIdentificationFromMultiGaussFit(mass_fit_for_Z, idr);
366  else
367  mass_id_success = (DeduceAfromPID(idr) > 0);
368  if (mass_id_success) {
369  // mass identification was at least attempted
370  // make sure grid's quality code is consistent with KVIdentificationResult
371  const_cast<KVIDZAFromZGrid*>(this)->fICode = idr->IDquality;
372  idr->Aident = kTRUE; // meaning A identification was attempted, even if it failed
373  }
374  else {
375  // the pid falls outside of any mass ranges for a Z which has assigned isotopes
376  // therefore although the Z identification was good, we cannot consider this
377  // particle to be identified
378  const_cast<KVIDZAFromZGrid*>(this)->fICode = kICODE4;
379  idr->IDquality = fICode; // otherwise identfication result quality code is not coherent with comment (see below)
380  }
381  idr->IDOK = (fICode < kICODE4);
382  }
383 
384  // ignore isotopic successful isotopic identification if fIgnoreMassID=true
385  if (fIgnoreMassID && idr->IDOK && idr->Aident) idr->Aident = false;
386 
387  // set comments in identification result
388  switch (fICode) {
389  case kICODE0:
390  idr->SetComment("ok");
391  break;
392  case kICODE1:
393  if (mass_id_success) idr->SetComment("slight ambiguity of A, which could be larger");
394  else idr->SetComment("slight ambiguity of Z, which could be larger");
395  break;
396  case kICODE2:
397  if (mass_id_success) idr->SetComment("slight ambiguity of A, which could be smaller");
398  else idr->SetComment("slight ambiguity of Z, which could be smaller");
399  break;
400  case kICODE3:
401  if (mass_id_success) idr->SetComment("slight ambiguity of A, which could be larger or smaller");
402  else idr->SetComment("slight ambiguity of Z, which could be larger or smaller");
403  break;
404  case kICODE4:
405  if (mass_id_success) idr->SetComment("point is outside of mass identification range");
406  else idr->SetComment("Z identification correct but no mass identification");
407  break;
408  case kICODE5:
409  if (mass_id_success) idr->SetComment("point is in between two isotopes A & A+2 (e.g. 5He, 8Be, 9B)");
410  else idr->SetComment("point is in between two lines of different Z, too far from either to be considered well-identified");
411  break;
412  case kICODE6:
413  idr->SetComment("(x,y) is below first line in grid");
414  break;
415  case kICODE7:
416  idr->SetComment("(x,y) is above last line in grid");
417  break;
418  default:
419  idr->SetComment("no identification: (x,y) out of range covered by grid");
420  }
421 }
422 
423 
424 
425 
431 
433 {
434  // First look for a set of mass intervals in which the PID of the identification result falls,
435  // if there is one (see KVIDZAFromZGrid::is_inside).
436  // If an interval set is found for a Z different to the original identification, idr->Z is changed.
437  // Then call interval_set::eval for the mass interval for this Z.
438 
439  int zint = is_inside(idr->PID);
440  if (!zint) return -1;
441  if (zint != idr->Z) idr->Z = zint;
442 
443  double res = 0.;
444  interval_set* it = GetIntervalSet(zint);
445  if (it) res = it->eval(idr);
446  return res;
447 }
448 
449 
450 
451 
453 
455 {
457  KVNumberList pids;
458  interval_set* itvs = 0;
459  TIter npid(GetIntervalSets());
460  while ((itvs = (interval_set*)npid())) {
461  if (!itvs->GetNPID()) continue;
462  pids.Add(itvs->GetZ());
463  }
464  GetParameters()->SetValue("PIDRANGE", pids.AsString());
465 
466  itvs = 0;
467  TIter next(GetIntervalSets());
468  while ((itvs = (interval_set*)next())) {
469  if (!itvs->GetNPID()) continue;
470  KVString par = Form("PIDRANGE%d", itvs->GetZ());
471  KVString val = "";
472  interval* itv = 0;
473  TIter ni(itvs->GetIntervals());
474  while ((itv = (interval*)ni())) {
475  val += Form("%d:%lf,%lf,%lf|", itv->GetA(), itv->GetPIDmin(), itv->GetPID(), itv->GetPIDmax());
476  }
477  val.Remove(val.Length() - 1);
478  GetParameters()->SetValue(par.Data(), val.Data());
479  }
480 }
481 
482 
483 
484 
485 
487 
489 {
490  double pid = idr->PID;
491  if (pid < 0.5) return 0.;
492  // calculate interpolated mass from PID
493  double res = fPIDs.Eval(pid);
494  int ares = 0;
495 
497 
498  // look for mass interval PID is in
499  // in case it falls between two intervals remember also the interval
500  // immediately to the left & right of the PID
501  interval* left_int(nullptr), *right_int(nullptr);
502  interval* inter;
503  TIter it(&fIntervals);
504  while ((inter = (interval*)it())) {
505  if (inter->is_inside(pid)) {
506  ares = inter->GetA();
507  break;
508  }
509  else if (inter->is_left_of(pid)) {
510  left_int = inter;
511  }
512  else if (!right_int && inter->is_right_of(pid)) {
513  right_int = inter;
514  }
515  }
516  if (ares != 0) {
517  // the PID is inside a defined mass interval
518  idr->A = ares;
519  idr->PID = res;
521  }
522  else {
523  // the PID is not inside a defined mass interval
524  //
525  // * if it is in between two consecutive masses i.e. A and A+1 then it is
526  // Z- and A-identified with a slight ambiguity of A
527  // * if it is in between two non-consecutive masses i.e. A and A+2 then it
528  // is not identified (e.g. 5He, 8Be, 9B)
529  if (!right_int || !left_int) {
530  // case where no left or right interval were found
531  // to prevent from crashes but should not appen
532  idr->A = ares;
533  idr->PID = res;
535  }
536  else {
537  int dA = right_int->GetA() - left_int->GetA();
538  if (dA == 1) {
539  // OK, slight ambiguity of A
540  ares = TMath::Nint(res);
541  idr->A = ares;
542  idr->PID = res;
544  }
545  else {
546  // in a hole where no isotopes should be (e.g. 5He, 8Be, 9B)
547  idr->A = ares;
548  idr->PID = res;
550  }
551  }
552  }
553  }
554  else {
555  ares = TMath::Nint(res);
556  idr->A = ares;
557  idr->PID = res;
558  if (ares > fPIDs.GetX()[0] && ares < fPIDs.GetX()[fNPIDs - 1]) {
560  }
561  else {
563  }
564  }
565  return res;
566 }
567 
568 
569 
571 
572 bool interval_set::is_inside(double pid)
573 {
574  if (fType != KVIDZAFromZGrid::kIntType) return kTRUE;
575 
576 // Info("is_inside","min: %d max:%d npids:%d", ((interval*)fIntervals.At(0))->GetA(), ((interval*)fIntervals.At(fNPIDs-1))->GetA(), fNPIDs);
577 
578  if (pid > ((interval*)fIntervals.At(0))->GetPIDmin() && pid < ((interval*)fIntervals.At(fNPIDs - 1))->GetPIDmax()) return kTRUE;
579  else return kFALSE;
580 }
581 
582 
583 
585 
586 bool interval_set::is_above(double pid)
587 {
588  if (fType != KVIDZAFromZGrid::kIntType) return kTRUE;
589 
590  if (pid > ((interval*)fIntervals.At(fNPIDs - 1))->GetPIDmax()) return kTRUE;
591  else return kFALSE;
592 }
593 
594 
595 
596 
598 
600 {
601  if (!GetNPID()) return "-";
602  KVNumberList alist;
603  for (int ii = 0; ii < GetNPID(); ii++) alist.Add(((interval*)fIntervals.At(ii))->GetA());
604  return alist.AsString();
605 }
606 
607 
608 
610 
611 interval_set::interval_set(int zz, int type)
612 {
613  fType = type;
614  fZ = zz;
615  fNPIDs = 0;
616 }
617 
618 
619 
625 
626 void interval_set::add(int aa, double pid, double pidmin, double pidmax)
627 {
628 // if (fNPIDs && pid < fPIDs.GetX()[fNPIDs - 1]) {
629 // Error("add", "Please give me peaks in the right order for Z=%d and A=%d...", fZ, aa);
630 // return;
631 // }
632  if (fType == KVIDZAFromZGrid::kIntType && !(pid > pidmin && pid < pidmax)) {
633  return;
634  }
635 
636  fPIDs.SetPoint(fNPIDs, pid, aa);
638  if (pid) fIntervals.AddLast(new interval(fZ, aa, pid, pidmin, pidmax));
639  }
640  fNPIDs++;
641 }
642 
643 
644 
645 
646 
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
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 winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char * Form(const char *fmt,...)
virtual void WriteToAsciiFile(std::ofstream &gridfile)
Definition: KVIDGraph.cpp:442
virtual void ReadFromAsciiFile(std::ifstream &gridfile)
Definition: KVIDGraph.cpp:601
const Char_t * GetName() const override
Definition: KVIDGraph.cpp:1339
const KVNameValueList * GetParameters() const
Definition: KVIDGraph.h:324
const KVList * GetIdentifiers() const
Definition: KVIDGraph.h:334
Hybrid charge & mass identification grid.
bool MassIdentificationFromMultiGaussFit(KVMultiGaussIsotopeFit *, KVIdentificationResult *) const
void Initialize() override
KVList fTables
intervals for mass id
interval_set * GetIntervalSet(int zint) const
void WriteToAsciiFile(std::ofstream &gridfile) override
void Identify(Double_t x, Double_t y, KVIdentificationResult *) const override
virtual double DeduceAfromPID(KVIdentificationResult *idr) const
KVUniqueNameList fFits
multi-gaussian fits for mass id
KVMultiGaussIsotopeFit * GetMultiGaussFit(int z) const
void SetOnlyZId(Bool_t=kTRUE) override
KVList * GetIntervalSets()
void ReadFromAsciiFile(std::ifstream &gridfile) override
void Copy(TObject &obj) const override
int is_inside(double pid) const
void Copy(TObject &) const override
Copy this to 'obj'.
Definition: KVIDZAGrid.cpp:66
Bool_t fHasMassIDRegion
Definition: KVIDZAGrid.h:95
Int_t fICode
code de retour
Definition: KVIDZAGrid.h:85
void Initialize() override
void Identify(Double_t x, Double_t y, KVIdentificationResult *) const override
Base class for graphical cuts used in particle identification.
Definition: KVIDentifier.h:28
Full result of one attempted particle identification.
Bool_t IDOK
general quality of identification, =kTRUE if acceptable identification made
void SetComment(const Char_t *c)
Bool_t Aident
= kTRUE if A of particle established
Double_t PID
= "real" Z if Zident==kTRUE and Aident==kFALSE, "real" A if Zident==Aident==kTRUE
Int_t A
A of particle found (if Aident==kTRUE)
Int_t Z
Z of particle found (if Zident==kTRUE)
Int_t IDquality
specific quality code returned by identification procedure
Bool_t HasFlag(std::string grid_name, TString flag)
Bool_t Zident
=kTRUE if Z of particle established
Function for fitting PID mass spectra.
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
void SetValue(const Char_t *name, value_type value)
void RemoveParameter(const Char_t *name)
const Char_t * GetStringValue(const Char_t *name) const
bool Set(const KVString &)
TString GetTStringValue(const Char_t *name) const
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:85
const Char_t * AsString(Int_t maxchars=0) const
void Add(Int_t)
Add value 'n' to the list.
void Copy(TObject &obj) const override
void Add(TObject *obj) override
void AddLast(TObject *obj) override
void Clear(Option_t *option="") override
TObject * At(Int_t idx) const override
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
void Begin(TString delim) const
Definition: KVString.cpp:565
Bool_t End() const
Definition: KVString.cpp:634
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:695
void Add(TObject *obj) override
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Double_t * GetX() const
virtual Double_t Eval(Double_t x, TSpline *spline=nullptr, Option_t *option="") const
virtual void SetName(const char *name)
Ssiz_t Length() const
Int_t Atoi() const
Double_t Atof() const
const char * Data() const
Bool_t IsWhitespace() const
TString & Remove(EStripType s, char c)
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
TString GetListOfMasses()
bool is_inside(double pid)
bool is_above(double pid)
void add(int aa, double pid, double pidmin=-1., double pidmax=-1.)
double eval(KVIdentificationResult *idr)
interval_set(int zz, int type)
KVList * GetIntervals()
bool is_right_of(double pid)
double GetPID()
bool is_left_of(double pid)
double GetPIDmin()
double GetPIDmax()
bool is_inside(double pid)
Double_t y[n]
Double_t x[n]
void init()
double min(double x, double y)
Int_t Nint(T x)
ClassImp(TPyArg)