KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVIDGCsI.cpp
Go to the documentation of this file.
1 #include "KVIDGCsI.h"
2 #include "KVIDCutLine.h"
3 #include "KVIDCsIRLLine.h"
5 
7 
8 
9 
13 {
14  //Default constructor
15  IMFLine = 0;
16  GammaLine = 0;
17 }
18 
19 
20 
23 
25 {
26  //Copy constructor
27  IMFLine = 0;
28  GammaLine = 0;
29 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
30  grid.Copy(*this);
31 #else
32  ((KVIDGCsI&) grid).Copy(*this);
33 #endif
34 }
35 
36 
37 
38 
41 
42 KVIDGCsI::~KVIDGCsI()
43 {
44  //Dtor.
45 }
46 
47 
48 
49 
55 
57 {
58  //Performs rejection of gammas - returns kTRUE if point above gamma line
59  //(note - position w.r.t. IMF line is not considered)
60  //Returns kFALSE for points below gamma line i.e. for gammas ;-)
61  //If no gamma line defined, returns kTRUE for all points
62 
63  if (!GammaLine) return kTRUE;
64  else if ((x < GammaLine->GetX()[0]) && (x < GammaLine->GetX()[GammaLine->GetN() - 1])) return kTRUE;
65  else if (GammaLine->WhereAmI(x, y, "above")) return kTRUE;
66 
67  if (reject) *reject = GammaLine->GetName();
68  return kFALSE;
69 }
70 
71 
72 
73 // small utility class which adds the IMF line to the identifier list when constructed,
74 // and removes it when it goes out of scope i.e. whenever the program leaves the Identify method
75 
77 
82  {
84  }
86  {
88  }
89 };
90 
91 
92 
104 
106 {
107  // Set Z and A of nucleus based on position in R-L grid
108  // The identification of gammas (kICODE10) and charged particles is performed
109  //
110  // ### Note
111  // for isotopically identified particles, the integer A (KVNucleus::GetA) is the mass assigned to the closest line
112  // [unless the closest line is the IMF line, in which case we use the closest identifier line],
113  // whereas the floating-point A (KVReconstructedNucleus::GetRealA) is calculated by interpolation.
114  // the integer A is not necessarily = nint(floating-point A): for example, if no 5He line is drawn in the grid
115  // (which is usually the case), there will be no isotopically-identified particle with GetA()=5, although
116  // there may be particles with GetRealA() between 4.5 and 5.5
117 
119 
120  if (!IsIdentifiable(x, y)) {
121  //point below gamma line
122  const_cast < KVIDGCsI* >(this)->fICode = kICODE10;
123  idr->IDquality = fICode;
124  idr->Z = 0;
125  idr->A = 0;
126  idr->IDOK = kTRUE;
127  idr->SetComment("gamma");
128  return;
129  }
130 
131  SetInfos(x, y, idr);
132 
133  if (!const_cast<KVIDGCsI*>(this)->FindFourEmbracingLines(x, y, "above")) {
134  //no lines corresponding to point were found
135  const_cast < KVIDGCsI* >(this)->fICode = kICODE8; // Z indetermine ou (x,y) hors limites
136  idr->IDquality = fICode;
137  idr->SetComment("no identification: (x,y) out of range covered by grid");
138  return;
139  }
140  Int_t Z;
141  Double_t A;
142  const_cast < KVIDGCsI* >(this)->IdentZA(x, y, Z, A);
143  idr->Z = Z;
144  idr->A = Aint;
145  idr->PID = A;
146  idr->IDquality = fICode;
147  switch (fICode) {
148 
149  case kICODE0:
150  idr->SetComment("ok");
151  break;
152  case kICODE1:
153  idr->SetComment("Z ok, mass may be greater than A");
154  break;
155  case kICODE2:
156  idr->SetComment("Z ok, mass may be smaller than A");
157  break;
158  case kICODE3:
159  idr->SetComment("Z ok, mass may be greater or smaller than A");
160  break;
161  case kICODE4:
162  idr->SetComment("Z ok, mass out of range, >=A");
163  break;
164  case kICODE5:
165  idr->SetComment("Z ok, mass out of range, <=A");
166  break;
167  case kICODE6:
168  idr->SetComment("point above IMF line, Z is a minimum value");
169  break;
170  case kICODE7:
171  idr->SetComment("point is left of IMF line, Z is the most probable lower limit");
172  break;
173  case kICODE8:
174  idr->SetComment("no identification: (x,y) out of range covered by grid");
175  break;
176  case kICODE9:
177  idr->SetComment("no identification for this module");
178  break;
179  default:
180  idr->SetComment("no identification: (x,y) out of range covered by grid");
181  }
182 
183  if (fICode < kICODE4) {
184  idr->IDOK = kTRUE;
185  idr->Zident = kTRUE;
186  idr->Aident = kTRUE;
187  }
188 }
189 
190 
191 
192 
198 
200 {
201  //Returns first ID line in sorted list for which GetZ() returns 'z'.
202  //index=index of line found in fIDLines list (-1 if not found).
203  //This is done by looping over all the lines in the list, not by dichotomy as in base class KVIDZAGrid,
204  //because of the 8Be line being in between 6He and 6Li.
205 
206  index = 0;
207  Int_t nlines = GetNumberOfIdentifiers();
208  while ((dynamic_cast < KVIDZALine* >(GetIdentifiers()->At(index))->GetZ() != z)
209  && (index < (nlines - 1)))
210  index++;
211  KVIDZALine* line = dynamic_cast < KVIDZALine* >(GetIdentifiers()->At(index));
212  if (line->GetZ() != z) {
213  index = -1;
214  return 0;
215  }
216  return line;
217 }
218 
219 
220 
221 
226 
228 {
229  //Returns ID line corresponding to nucleus (Z,A) and its index in fIDLines.
230  //This is done by looping over all the lines in the list, not as in base class KVIDZAGrid,
231  //because of the 8Be line being in between 6He and 6Li.
232 
233  index = 0;
234  Int_t nlines = GetNumberOfIdentifiers();
235  KVIDZALine* line = dynamic_cast < KVIDZALine* >(GetIdentifiers()->At(index));
236  while ((line->GetZ() != z || line->GetA() != a)
237  && (index < (nlines - 1))) {
238  line = dynamic_cast < KVIDZALine* >(GetIdentifiers()->At(++index));
239  }
240  if (line->GetZ() != z || line->GetA() != a) {
241  index = -1;
242  return 0;
243  }
244  return line;
245 }
246 
247 
248 
249 
255 
257 {
258  //Finds Z, A and 'real A' for point (x,y) once closest lines to point have been found.
259  // Double_t A = mass calculated by interpolation
260  //This is a line-for-line copy of the latter part of IdnCsOr, even the same
261  //variable names and comments have been used (as much as possible).
262 
263  fICode = kICODE0;
264  A = -1.;
265  Aint = 0;
266 
267 // if(fIdxClosest==ksups) cout << "*** ";
268 // cout << "ksups = " << ksups << " Zsups = " << Zsups << " Asups = " << Asups << " wsups = " << wsups << " dsups = " << dsups << endl;
269 // if(fIdxClosest==ksup) cout << "*** ";
270 // cout << "ksup = " << ksup << " Zsup = " << Zsup << " Asup = " << Asup << " wsup = " << wsup << " dsup = " << dsup << endl;
271 // if(fIdxClosest==kinf) cout << "*** ";
272 // cout << "kinf = " << kinf << " Zinf = " << Zinf << " Ainf = " << Ainf << " winf = " << winf << " dinf = " << dinf << endl;
273 // if(fIdxClosest==kinfi) cout << "*** ";
274 // cout << "kinfi = " << kinfi << " Zinfi = " << Zinfi << " Ainfi = " << Ainfi << " winfi = " << winfi << " dinfi = " << dinfi << endl;
275  Int_t ibif = 0;
276  Int_t k = -1;
277  Double_t yy, y1, y2;
278  Int_t ix1, ix2;
279  yy = y1 = y2 = 0;
280  ix1 = ix2 = 0;
281  if (ksup > -1) {
282  if (kinf > -1) {
283  //cout << " /******************* 2 lignes encadrant le point ont ete trouvees ************************/" << endl;
284  Double_t dt = dinf + dsup; //distance between the 2 lines
285  if (Zinf == Zsup) {
286  // cout << " /****************meme Z**************/" << endl;
287  Z = Zinf;
288  Int_t dA = Asup - Ainf;
289  Double_t dist = dt / dA; //distance between the 2 lines normalised to difference in A of lines
290  /*** A = Asup ***/
291  if (dinf > dsup) { //point is closest to upper line, 'sup' line
292  ibif = 1;
293  k = ksup;
294  yy = -dsup;
295  A = Asup;
296  Aint = Asup;
297  if (ksups > -1) { // there is a 'sups' line above the 2 which encadrent le point
298  y2 = dsups - dsup;
299  if (Zsups == Zsup) {
300  ibif = 0;
301  y2 /= 2.;
302  ix2 = Asups - Asup;
303  }
304  else {
305  if (Zsups > 0)
306  y2 /= 2.; // 'sups' line is not IMF line
307  Double_t x2 = wsup;
308  x2 = 0.5 * TMath::Max(x2, dist);
309  y2 = TMath::Min(y2, x2);
310  ix2 = 1;
311  }
312  }
313  else { // ksups == -1 i.e. no 'sups' line
314  y2 = wsup;
315  y2 = 0.5 * TMath::Max(y2, dist);
316  ix2 = 1;
317  }
318  y1 = -dt / 2.;
319  ix1 = -dA;
320  }
321  /*** A = Ainf ***/
322  else { //point is closest to lower line, 'inf' line
323  ibif = 2;
324  k = kinf;
325  yy = dinf;
326  A = Ainf;
327  Aint = Ainf;
328  if (kinfi > -1) { // there is a 'infi' line below the 2 which encadrent le point
329  y1 = 0.5 * (dinfi - dinf);
330  if (Zinfi == Zinf) {
331  ibif = 0;
332  ix1 = Ainfi - Ainf;
333  y1 = -y1;
334  }
335  else {
336  Double_t x1 = winf;
337  x1 = 0.5 * TMath::Max(x1, dist);
338  y1 = -TMath::Min(y1, x1);
339  ix1 = -1;
340  }
341  }
342  else { // kinfi = -1 i.e. no 'infi' line
343  y1 = winf;
344  y1 = -0.5 * TMath::Max(y1, dist);
345  ix1 = -1;
346  }
347  y2 = dt / 2.;
348  ix2 = dA;
349  }
350  }
351  else {
352  //cout << " /****************Z differents**************/ " << endl;
353  if (Zsup == -1) { //'sup' is the IMF line
354  dt *= 2.;
355  dsup = dt - dinf;
356  }
357  /*** Z = Zsup ***/
358  ibif = 3;
359  if (dinf > dsup) { // closest to upper 'sup' line
360  k = ksup;
361  yy = -dsup;
362  Z = Zsup;
363  A = Asup;
364  Aint = Asup;
365  y1 = 0.5 * wsup;
366  if (ksups > -1) { // there is a 'sups' line above the 2 which encadrent the point
367  y2 = dsups - dsup;
368  if (Zsups == Zsup) {
369  ibif = 2;
370  ix2 = Asups - Asup;
371  Double_t x1 = y2 / ix2 / 2.;
372  y1 = TMath::Max(y1, x1);
373  y1 = -TMath::Min(y1, dt / 2.);
374  ix1 = -1;
375  y2 /= 2.;
376  }
377  else {
378  if (Zsups > 0)
379  y2 /= 2.; // 'sups" is not IMF line
380  y2 = TMath::Min(y1, y2);
381  ix2 = 1;
382  y1 = -TMath::Min(y1, dt / 2.);
383  ix1 = -1;
384  }
385  }
386  else { // ksups == -1, i.e. no 'sups' line
387  fICode = kICODE7; //a gauche de la ligne fragment, Z est alors un Zmin et le plus probable
388  y2 = y1;
389  ix2 = 1;
390  y1 = -TMath::Min(y1, dt / 2.);
391  ix1 = -1;
392  }
393  }
394  /*** Z = Zinf ***/
395  else { // closest to lower 'inf' line
396  k = kinf;
397  yy = dinf;
398  Z = Zinf;
399  A = Ainf;
400  Aint = Ainf;
401  y2 = 0.5 * winf;
402  if (kinfi > -1) { // there is a 'infi' line below the 2 which encadrent the point
403  y1 = dinfi - dinf;
404  if (Zinfi == Zinf) {
405  ibif = 1;
406  ix1 = Ainfi - Ainf;
407  Double_t x2 = -y1 / ix1 / 2.;
408  y2 = TMath::Max(y2, x2);
409  y2 = TMath::Min(y2, dt / 2.);
410  ix2 = 1;
411  y1 /= -2.;
412  }
413  else {
414  y1 = -TMath::Min(y2, y1 / 2.);
415  ix1 = -1;
416  y2 = TMath::Min(y2, dt / 2.);
417  ix2 = 1;
418  }
419  }
420  else { // no kinfi line, kinfi==-1
421  y1 = -y2;
422  ix1 = -1;
423  y2 = TMath::Min(y2, dt / 2.);
424  ix1 = 1;
425  }
426  }
427  }
428  }//if(kinf>-1)...
429  else if (Zsup > 0) { // 'sup' is not IMF line
430  //cout<<" /****************** Seule une ligne superieure a ete trouvee *********************/" << endl;
431  ibif = 3;
432  k = ksup;
433  yy = -dsup;
434  Z = Zsup;
435  A = Asup;
436  Aint = Asup;
437  y1 = 0.5 * wsup;
438  if (ksups > -1) { // there is a 'sups' line above the closest line to the point
439  y2 = dsups - dsup;
440  if (Zsups == Zsup) {
441  ibif = 2;
442  ix2 = Asups - Asup;
443  Double_t x1 = y2 / ix2 / 2.;
444  y1 = -TMath::Max(y1, x1);
445  ix1 = -1;
446  y2 /= 2.;
447  }
448  else {
449  if (Zsups > 0)
450  y2 /= 2.; // 'sups' is not IMF line
451  y2 = TMath::Min(y1, y2);
452  ix2 = 1;
453  y1 = -y1;
454  ix1 = -1;
455  }
456  }
457  else { // no 'sups' line above closest line
458  fICode = kICODE7; //a gauche de la ligne fragment, Z est alors un Zmin et le plus probable
459  y2 = y1;
460  ix2 = 1;
461  y1 = -y1;
462  ix1 = -1;
463  }
464  }
465  else {
466  fICode = kICODE8; // Z indetermine ou (x,y) hors limites
467  }
468  }
469  else if (kinf > -1) {
470  //cout <<"/****************** Seule une ligne inferieure a ete trouvee ***********************/" << endl;
471  /*** Sep. fragment ***/
472  if (Zinf == -1) { // 'inf' is IMF line
473  //point is above IMF line. Z = Z of last line in grid, A = -1
474  k = -1;
475  Z = GetZmax();
476  A = -1;
477  Aint = 0;
478  fICode = kICODE6; // au-dessus de la ligne fragment, Z est alors un Zmin
479  }
480  /*** Ligne de crete (Z,A line)***/
481  else {
482  ibif = 3;
483  k = kinf;
484  Z = Zinf;
485  A = Ainf;
486  Aint = Ainf;
487  yy = dinf;
488  y2 = 0.5 * winf;
489  if (kinfi > -1) {
490  y1 = dinfi - dinf;
491  if (Zinfi == Zinf) {
492  ibif = 1;
493  ix1 = Ainfi - Ainf;
494  Double_t x2 = -y1 / ix1 / 2.;
495  y2 = TMath::Max(y2, x2);
496  ix2 = 1;
497  y1 /= -2.;
498  }
499  else {
500  y1 = -TMath::Min(y2, y1 / 2.);
501  ix1 = -1;
502  ix2 = 1;
503  }
504  }
505  else {
506  y1 = -y2;
507  ix1 = -1;
508  ix2 = 1;
509  }
510  fICode = kICODE7; // a gauche de la ligne fragment, Z est alors un Zmin et le plus probable
511  }
512  }
513  /*****************Aucune ligne n'a ete trouvee*********************************/
514  else {
515  fICode = kICODE8; // Z indetermine ou (x,y) hors limites
516  }
517  /****************Test des bornes********************************************/
518  if (k > -1 && fICode == kICODE0) {
519  if (yy > y2)
520  fICode = kICODE4; // Z ok, masse hors limite superieure ou egale a A
521  }
522  if (k > -1 && (fICode == kICODE0 || fICode == kICODE7)) {
523  if (yy < y1)
524  fICode = kICODE5; // Z ok, masse hors limite inferieure ou egale a A
525  }
526  if (fICode == kICODE4 || fICode == kICODE5) {
527  A = -1;
528  Aint = 0;
529  }
530 
531  /****************Interpolation de la masse: da = f*log(1+b*dy)********************/
532  if (fICode == kICODE0 || (fICode == kICODE7 && yy <= y2)) {
533  Double_t deltaA = 0.;
534  Bool_t i = kFALSE;
535  Double_t dt, dist = y1 * y2;
536  dt = 0.;
537  if (ix2 == -ix1) { //dA1 = dA2
538  if (dist != 0) {
539  dt = -(y1 + y2) / dist;
540  i = kTRUE;
541  }
542  /*else
543  Warning("IdentZA","%s : cannot calculate interpolated mass (dist=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
544  GetName(), dist, Z, Aint, fICode);*/
545  }
546  else if (ix2 == -ix1 * 2) { // dA2 = 2*dA1
547  Double_t tmp = y1 * y1 - 4. * dist;
548  if (tmp > 0. && dist != 0) {
549  dt = -(y1 + 2. * y2 -
550  TMath::Sqrt(tmp)) / dist / 2.;
551  i = kTRUE;
552  }
553  /*else
554  Warning("IdentZA","%s : cannot calculate interpolated mass (y1*y1-4*dist=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
555  GetName(), tmp, Z, Aint, fICode);*/
556  }
557  else if (ix1 == -ix2 * 2) { // dA1 = 2*dA2
558  Double_t tmp = y2 * y2 - 4. * dist;
559  if (tmp > 0. && dist != 0) {
560  dt = -(y2 + 2. * y1 +
561  TMath::Sqrt(tmp)) / dist / 2.;
562  i = kTRUE;
563  }
564  /*else
565  Warning("IdentZA","%s : cannot calculate interpolated mass (y2*y2-4*dist=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
566  GetName(), tmp, Z, Aint, fICode);*/
567  }
568  if (i) {
569  dist = dt * y2;
570  if (TMath::Abs(dist) < 0.001) {
571  if (y2 != 0)
572  deltaA = yy * ix2 / y2 / 2.;
573  /*else
574  Warning("IdentZA","%s : cannot calculate interpolated mass (y2=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
575  GetName(), y2, Z, Aint, fICode);*/
576  }
577  else {
578  if (dist > -1. && dt * yy > -1.)
579  deltaA = ix2 / 2. / TMath::Log(1. + dist) * TMath::Log(1. + dt * yy);
580  /*else
581  Warning("IdentZA","%s : cannot calculate interpolated mass (dist=%f dt*yy=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
582  GetName(), dist, dt*yy, Z, Aint, fICode);*/
583  }
584  A += deltaA;
585  }
586  }
587  /***************D'autres masses sont-elles possibles ?*************************/
588  if (fICode == kICODE0) {
589  //cout << "icode = 0, ibif = " << ibif << endl;
590  /***Masse superieure***/
591  if (ibif == 1 || ibif == 3) {
592  //We look at next line in the complete list of lines, after the closest line.
593  //If it has the same Z as the closest line, but was excluded from research for closest line
594  //because the point lies outside the endpoints, there remains a doubt about the mass:
595  //on rajoute 1 a fICode, effectivement on le met = kICODE1
596  Int_t idx = fIdxClosest;
597  if (idx > -1 && ++idx < GetNumberOfIdentifiers()) {
598  KVIDCsIRLLine* nextline =
600  if (nextline->GetZ() == Z
601  && !nextline->IsBetweenEndPoints(x, y, "x")) {
602  fICode++; // Z ok, mais les masses superieures a A sont possibles
603  //cout <<"//on rajoute 1 a fICode, effectivement on le met = kICODE1" << endl;
604  }
605  }
606  }
607  /***Masse inferieure***/
608  if (ibif == 2 || ibif == 3) {
609  //We look at line below the closest line in the complete list of lines.
610  //If it has the same Z as the closest line, but was excluded from research for closest line
611  //because the point lies outside the endpoints, there remains a doubt about the mass:
612  //on rajoute 2 a fICode, so it can be = kICODE2 or kICODE3
613  Int_t idx = fIdxClosest;
614  if (idx > -1 && --idx >= 0) {
615  KVIDCsIRLLine* nextline =
617  if (nextline->GetZ() == Z
618  && !nextline->IsBetweenEndPoints(x, y, "x")) {
619  fICode += 2;
620  //cout << "//on rajoute 2 a fICode, so it can be = kICODE2 or kICODE3" << endl;
621  }
622  }
623  }
624  }
625 
626  // cout << "Z = " << Z << " A = " << A << " icode = " << fICode << endl;
627 }
628 
629 
630 
631 
640 
642 {
643  // General initialisation method for CsI R-L identification grid.
644  // This method MUST be called once before using the grid for identifications.
645  // The ID lines are sorted.
646  // The natural line widths of all ID lines are calculated.
647  // The line with the largest Z (Zmax line) is found.
648  // IMF & Gamma line pointers are initialised
649 
650  // if grid has already been used for identification, IMF_line will be in identifiers list.
651  TObject* imfline = fIdentifiers->FindObject("IMF_line");
652  if (imfline) fIdentifiers->Remove(imfline); // remove to avoid problems with CalculateLineWidths
654  GammaLine = (KVIDLine*)GetCut("gamma_line");
655  if (!GammaLine) {
656  GetName();
657  Error("Initialize", "%s: Cut 'gamma_line' not found in grid. Not a valid CsI R-L grid. Listing existing cuts:",
658  GetName());
659  GetCuts()->ls();
660  }
661  IMFLine = (KVIDLine*)GetCut("IMF_line");
662  if (!IMFLine) {
663  Error("Initialize", "%s: Cut 'IMF_line' not found in grid. Not a valid CsI R-L grid. Listing existing cuts:",
664  GetName());
665  GetCuts()->ls();
666  }
667 }
668 
669 
670 
671 
683 
685 {
686  // Called after reading a grid from an ascii file.
687  // Tries to convert information written by an old version of the class:
688  //
689  //<PARAMETER> Ring min=... ----> <PARAMETER> IDTelescopes=...
690  //<PARAMETER> Ring max=...
691  //<PARAMETER> Mod min=...
692  //<PARAMETER> Mod max=...
693  //
694  // This will fail. The fix is no longer supported. Such files should
695  // no longer exist.
696 
698  if (fPar->HasParameter("IDTelescopes")) return;
699 
700  Fatal("BackwardsCompatibilityFix",
701  "This fix no longer works. There will be problems.");
702  GammaLine = (KVIDLine*)GetCut("gamma_line");
703  IMFLine = (KVIDLine*)GetCut("IMF_line");
704  if (GammaLine)((KVIDCutLine*)GammaLine)->SetAcceptedDirection("above");
705  if (IMFLine)((KVIDCutLine*)IMFLine)->SetAcceptedDirection("below");
706  SetVarY("CSI-R");
707  SetVarX("CSI-L");
708 }
709 
710 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
Base class for CsI R-L identification lines (A and Z identification).
Definition: KVIDCsIRLLine.h:26
Line in ID grid used to delimit regions where no identification is possible.
Definition: KVIDCutLine.h:22
Identification grids for CsI R-L (fast-slow) matrices.
Definition: KVIDGCsI.h:41
KVIDLine * GammaLine
Definition: KVIDGCsI.h:44
KVIDZALine * GetZALine(Int_t z, Int_t a, Int_t &) const
Definition: KVIDGCsI.cpp:227
void IdentZA(Double_t x, Double_t y, Int_t &Z, Double_t &A)
Definition: KVIDGCsI.cpp:256
KVIDLine * IMFLine
Definition: KVIDGCsI.h:43
KVIDZALine * GetZLine(Int_t z, Int_t &) const
Definition: KVIDGCsI.cpp:199
virtual void Initialize()
Definition: KVIDGCsI.cpp:641
virtual void BackwardsCompatibilityFix()
Definition: KVIDGCsI.cpp:684
virtual Bool_t IsIdentifiable(Double_t x, Double_t y, TString *=nullptr) const
Definition: KVIDGCsI.cpp:56
virtual void Identify(Double_t x, Double_t y, KVIdentificationResult *) const
Definition: KVIDGCsI.cpp:105
KVIDGCsI()
Default constructor.
Definition: KVIDGCsI.cpp:12
virtual void SetVarX(const char *v)
Definition: KVIDGraph.h:506
Int_t GetNumberOfIdentifiers() const
Definition: KVIDGraph.h:300
KVIDentifier * GetCut(const Char_t *name) const
Definition: KVIDGraph.h:272
KVNameValueList * fPar
parameters associated to grid
Definition: KVIDGraph.h:44
KVList * fIdentifiers
list of identification objects
Definition: KVIDGraph.h:39
KVIDentifier * GetIdentifierAt(Int_t index) const
Definition: KVIDGraph.h:263
virtual void SetVarY(const char *v)
Definition: KVIDGraph.h:510
const Char_t * GetName() const
Definition: KVIDGraph.cpp:1320
virtual void BackwardsCompatibilityFix()
Definition: KVIDGraph.cpp:1349
virtual void SetInfos(Double_t, Double_t, KVIdentificationResult *) const
loop over KVIDGraph::fInfoZones to set flags in KVIdentificationResult
Definition: KVIDGraph.cpp:1285
KVList * GetCuts() const
Definition: KVIDGraph.h:290
KVList * GetIdentifiers() const
Definition: KVIDGraph.h:285
Base class for lines/cuts used for particle identification in 2D data maps.
Definition: KVIDLine.h:143
Bool_t WhereAmI(Double_t px, Double_t py, Option_t *opt)
Definition: KVIDLine.h:348
Bool_t IsBetweenEndPoints(Double_t x, Double_t y, const Char_t *axis="") const
Definition: KVIDLine.h:442
Identification grid with lines corresponding to different nuclear isotopes (KVIDZALine)
Definition: KVIDZAGrid.h:65
Double_t dsups
Definition: KVIDZAGrid.h:87
Int_t Zsups
Definition: KVIDZAGrid.h:89
Double_t wsup
Definition: KVIDZAGrid.h:88
Int_t Zinf
Definition: KVIDZAGrid.h:89
Int_t Ainfi
Definition: KVIDZAGrid.h:90
Int_t Ainf
Definition: KVIDZAGrid.h:90
virtual void Copy(TObject &) const
Copy this to 'obj'.
Definition: KVIDZAGrid.cpp:71
Int_t Zsup
Definition: KVIDZAGrid.h:89
Int_t Asups
Definition: KVIDZAGrid.h:90
Int_t kinfi
Definition: KVIDZAGrid.h:86
Int_t Asup
Definition: KVIDZAGrid.h:90
Int_t ksup
Definition: KVIDZAGrid.h:86
Int_t ksups
used by IdentZA and IdentZ
Definition: KVIDZAGrid.h:86
Int_t GetZmax() const
Definition: KVIDZAGrid.h:141
Double_t winf
Definition: KVIDZAGrid.h:88
Double_t dinf
Definition: KVIDZAGrid.h:87
Int_t kinf
Definition: KVIDZAGrid.h:86
virtual void Initialize()
Double_t dsup
Definition: KVIDZAGrid.h:87
Int_t fIdxClosest
index of closest line in main list fIdentifiers
Definition: KVIDZAGrid.h:83
Int_t Zinfi
Definition: KVIDZAGrid.h:89
Int_t fICode
code de retour
Definition: KVIDZAGrid.h:84
Double_t dinfi
Definition: KVIDZAGrid.h:87
virtual Bool_t FindFourEmbracingLines(Double_t x, Double_t y, const Char_t *position)
Definition: KVIDZAGrid.cpp:525
Int_t Aint
mass of line used to identify particle
Definition: KVIDZAGrid.h:92
Base class for identification ridge lines corresponding to different nuclear species.
Definition: KVIDZALine.h:32
virtual Int_t GetZ() const
Definition: KVIDentifier.h:78
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 Zident
=kTRUE if Z of particle established
Extended TList class which owns its objects by default.
Definition: KVList.h:27
Bool_t HasParameter(const Char_t *name) const
virtual void AddLast(TObject *obj)
virtual TObject * At(Int_t idx) const
virtual TObject * Remove(TObject *obj)
Remove object from list.
virtual TObject * FindObject(const char *name) const
virtual void ls(Option_t *option="") const
Int_t GetN() const
Double_t * GetX() const
virtual const char * GetName() const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Fatal(const char *method, const char *msgfmt,...) const
TLine * line
Double_t y[n]
Double_t x[n]
double dist(AxisAngle const &r1, AxisAngle const &r2)
Double_t Min(Double_t a, Double_t b)
Double_t Log(Double_t x)
Double_t Sqrt(Double_t x)
Double_t Abs(Double_t d)
Double_t Max(Double_t a, Double_t b)
KVIDLine * imfline
Definition: KVIDGCsI.cpp:80
add_remove_imf_line(KVList *l, KVIDLine *i)
Definition: KVIDGCsI.cpp:81
auto * l
auto * a