KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVIDGCsI.cpp
1#include "KVIDGCsI.h"
2#include "KVIDCutLine.h"
3#include "KVIDCsIRLLine.h"
4#include "KVIdentificationResult.h"
5
7
8
9
11
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
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
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 SetInfos(x, y, idr);
121
122 if (!const_cast<KVIDGCsI*>(this)->FindFourEmbracingLines(x, y, "above")) {
123 //no lines corresponding to point were found
124 const_cast < KVIDGCsI* >(this)->fICode = kICODE8; // Z indetermine ou (x,y) hors limites
125 idr->IDquality = fICode;
126 idr->SetComment("no identification: (x,y) out of range covered by grid");
127 return;
128 }
129 Int_t Z;
130 Double_t A;
131 const_cast < KVIDGCsI* >(this)->IdentZA(x, y, Z, A);
132 idr->Z = Z;
133 idr->A = Aint;
134 idr->PID = A;
135 idr->IDquality = fICode;
136 switch (fICode) {
137
138 case kICODE0:
139 idr->SetComment("ok");
140 break;
141 case kICODE1:
142 idr->SetComment("Z ok, mass may be greater than A");
143 break;
144 case kICODE2:
145 idr->SetComment("Z ok, mass may be smaller than A");
146 break;
147 case kICODE3:
148 idr->SetComment("Z ok, mass may be greater or smaller than A");
149 break;
150 case kICODE4:
151 idr->SetComment("Z ok, mass out of range, >=A");
152 break;
153 case kICODE5:
154 idr->SetComment("Z ok, mass out of range, <=A");
155 break;
156 case kICODE6:
157 idr->SetComment("point above IMF line, Z is a minimum value");
158 break;
159 case kICODE7:
160 idr->SetComment("point is left of IMF line, Z is the most probable lower limit");
161 break;
162 case kICODE8:
163 idr->SetComment("no identification: (x,y) out of range covered by grid");
164 break;
165 case kICODE9:
166 idr->SetComment("no identification for this module");
167 break;
168 default:
169 idr->SetComment("no identification: (x,y) out of range covered by grid");
170 }
171
172 if (fICode < kICODE4) {
173 idr->IDOK = kTRUE;
174 idr->Zident = kTRUE;
175 idr->Aident = kTRUE;
176 }
177}
178
179
180
181
187
189{
190 //Returns first ID line in sorted list for which GetZ() returns 'z'.
191 //index=index of line found in fIDLines list (-1 if not found).
192 //This is done by looping over all the lines in the list, not by dichotomy as in base class KVIDZAGrid,
193 //because of the 8Be line being in between 6He and 6Li.
194
195 index = 0;
196 Int_t nlines = GetNumberOfIdentifiers();
197 while ((dynamic_cast < KVIDZALine* >(GetIdentifiers()->At(index))->GetZ() != z)
198 && (index < (nlines - 1)))
199 index++;
200 KVIDZALine* line = dynamic_cast < KVIDZALine* >(GetIdentifiers()->At(index));
201 if (line->GetZ() != z) {
202 index = -1;
203 return 0;
204 }
205 return line;
206}
207
208
209
210
215
217{
218 //Returns ID line corresponding to nucleus (Z,A) and its index in fIDLines.
219 //This is done by looping over all the lines in the list, not as in base class KVIDZAGrid,
220 //because of the 8Be line being in between 6He and 6Li.
221
222 index = 0;
223 Int_t nlines = GetNumberOfIdentifiers();
224 KVIDZALine* line = dynamic_cast < KVIDZALine* >(GetIdentifiers()->At(index));
225 while ((line->GetZ() != z || line->GetA() != a)
226 && (index < (nlines - 1))) {
227 line = dynamic_cast < KVIDZALine* >(GetIdentifiers()->At(++index));
228 }
229 if (line->GetZ() != z || line->GetA() != a) {
230 index = -1;
231 return 0;
232 }
233 return line;
234}
235
236
237
238
244
246{
247 //Finds Z, A and 'real A' for point (x,y) once closest lines to point have been found.
248 // Double_t A = mass calculated by interpolation
249 //This is a line-for-line copy of the latter part of IdnCsOr, even the same
250 //variable names and comments have been used (as much as possible).
251
252 fICode = kICODE0;
253 A = -1.;
254 Aint = 0;
255
256// if(fIdxClosest==ksups) cout << "*** ";
257// cout << "ksups = " << ksups << " Zsups = " << Zsups << " Asups = " << Asups << " wsups = " << wsups << " dsups = " << dsups << endl;
258// if(fIdxClosest==ksup) cout << "*** ";
259// cout << "ksup = " << ksup << " Zsup = " << Zsup << " Asup = " << Asup << " wsup = " << wsup << " dsup = " << dsup << endl;
260// if(fIdxClosest==kinf) cout << "*** ";
261// cout << "kinf = " << kinf << " Zinf = " << Zinf << " Ainf = " << Ainf << " winf = " << winf << " dinf = " << dinf << endl;
262// if(fIdxClosest==kinfi) cout << "*** ";
263// cout << "kinfi = " << kinfi << " Zinfi = " << Zinfi << " Ainfi = " << Ainfi << " winfi = " << winfi << " dinfi = " << dinfi << endl;
264 Int_t ibif = 0;
265 Int_t k = -1;
266 Double_t yy, y1, y2;
267 Int_t ix1, ix2;
268 yy = y1 = y2 = 0;
269 ix1 = ix2 = 0;
270 if (ksup > -1) {
271 if (kinf > -1) {
272 //cout << " /******************* 2 lignes encadrant le point ont ete trouvees ************************/" << endl;
273 Double_t dt = dinf + dsup; //distance between the 2 lines
274 if (Zinf == Zsup) {
275 // cout << " /****************meme Z**************/" << endl;
276 Z = Zinf;
277 Int_t dA = Asup - Ainf;
278 Double_t dist = dt / dA; //distance between the 2 lines normalised to difference in A of lines
279 /*** A = Asup ***/
280 if (dinf > dsup) { //point is closest to upper line, 'sup' line
281 ibif = 1;
282 k = ksup;
283 yy = -dsup;
284 A = Asup;
285 Aint = Asup;
286 if (ksups > -1) { // there is a 'sups' line above the 2 which encadrent le point
287 y2 = dsups - dsup;
288 if (Zsups == Zsup) {
289 ibif = 0;
290 y2 /= 2.;
291 ix2 = Asups - Asup;
292 }
293 else {
294 if (Zsups > 0)
295 y2 /= 2.; // 'sups' line is not IMF line
296 Double_t x2 = wsup;
297 x2 = 0.5 * TMath::Max(x2, dist);
298 y2 = TMath::Min(y2, x2);
299 ix2 = 1;
300 }
301 }
302 else { // ksups == -1 i.e. no 'sups' line
303 y2 = wsup;
304 y2 = 0.5 * TMath::Max(y2, dist);
305 ix2 = 1;
306 }
307 y1 = -dt / 2.;
308 ix1 = -dA;
309 }
310 /*** A = Ainf ***/
311 else { //point is closest to lower line, 'inf' line
312 ibif = 2;
313 k = kinf;
314 yy = dinf;
315 A = Ainf;
316 Aint = Ainf;
317 if (kinfi > -1) { // there is a 'infi' line below the 2 which encadrent le point
318 y1 = 0.5 * (dinfi - dinf);
319 if (Zinfi == Zinf) {
320 ibif = 0;
321 ix1 = Ainfi - Ainf;
322 y1 = -y1;
323 }
324 else {
325 Double_t x1 = winf;
326 x1 = 0.5 * TMath::Max(x1, dist);
327 y1 = -TMath::Min(y1, x1);
328 ix1 = -1;
329 }
330 }
331 else { // kinfi = -1 i.e. no 'infi' line
332 y1 = winf;
333 y1 = -0.5 * TMath::Max(y1, dist);
334 ix1 = -1;
335 }
336 y2 = dt / 2.;
337 ix2 = dA;
338 }
339 }
340 else {
341 //cout << " /****************Z differents**************/ " << endl;
342 if (Zsup == -1) { //'sup' is the IMF line
343 dt *= 2.;
344 dsup = dt - dinf;
345 }
346 /*** Z = Zsup ***/
347 ibif = 3;
348 if (dinf > dsup) { // closest to upper 'sup' line
349 k = ksup;
350 yy = -dsup;
351 Z = Zsup;
352 A = Asup;
353 Aint = Asup;
354 y1 = 0.5 * wsup;
355 if (ksups > -1) { // there is a 'sups' line above the 2 which encadrent the point
356 y2 = dsups - dsup;
357 if (Zsups == Zsup) {
358 ibif = 2;
359 ix2 = Asups - Asup;
360 Double_t x1 = y2 / ix2 / 2.;
361 y1 = TMath::Max(y1, x1);
362 y1 = -TMath::Min(y1, dt / 2.);
363 ix1 = -1;
364 y2 /= 2.;
365 }
366 else {
367 if (Zsups > 0)
368 y2 /= 2.; // 'sups" is not IMF line
369 y2 = TMath::Min(y1, y2);
370 ix2 = 1;
371 y1 = -TMath::Min(y1, dt / 2.);
372 ix1 = -1;
373 }
374 }
375 else { // ksups == -1, i.e. no 'sups' line
376 fICode = kICODE7; //a gauche de la ligne fragment, Z est alors un Zmin et le plus probable
377 y2 = y1;
378 ix2 = 1;
379 y1 = -TMath::Min(y1, dt / 2.);
380 ix1 = -1;
381 }
382 }
383 /*** Z = Zinf ***/
384 else { // closest to lower 'inf' line
385 k = kinf;
386 yy = dinf;
387 Z = Zinf;
388 A = Ainf;
389 Aint = Ainf;
390 y2 = 0.5 * winf;
391 if (kinfi > -1) { // there is a 'infi' line below the 2 which encadrent the point
392 y1 = dinfi - dinf;
393 if (Zinfi == Zinf) {
394 ibif = 1;
395 ix1 = Ainfi - Ainf;
396 Double_t x2 = -y1 / ix1 / 2.;
397 y2 = TMath::Max(y2, x2);
398 y2 = TMath::Min(y2, dt / 2.);
399 ix2 = 1;
400 y1 /= -2.;
401 }
402 else {
403 y1 = -TMath::Min(y2, y1 / 2.);
404 ix1 = -1;
405 y2 = TMath::Min(y2, dt / 2.);
406 ix2 = 1;
407 }
408 }
409 else { // no kinfi line, kinfi==-1
410 y1 = -y2;
411 ix1 = -1;
412 y2 = TMath::Min(y2, dt / 2.);
413 ix1 = 1;
414 }
415 }
416 }
417 }//if(kinf>-1)...
418 else if (Zsup > 0) { // 'sup' is not IMF line
419 //cout<<" /****************** Seule une ligne superieure a ete trouvee *********************/" << endl;
420 ibif = 3;
421 k = ksup;
422 yy = -dsup;
423 Z = Zsup;
424 A = Asup;
425 Aint = Asup;
426 y1 = 0.5 * wsup;
427 if (ksups > -1) { // there is a 'sups' line above the closest line to the point
428 y2 = dsups - dsup;
429 if (Zsups == Zsup) {
430 ibif = 2;
431 ix2 = Asups - Asup;
432 Double_t x1 = y2 / ix2 / 2.;
433 y1 = -TMath::Max(y1, x1);
434 ix1 = -1;
435 y2 /= 2.;
436 }
437 else {
438 if (Zsups > 0)
439 y2 /= 2.; // 'sups' is not IMF line
440 y2 = TMath::Min(y1, y2);
441 ix2 = 1;
442 y1 = -y1;
443 ix1 = -1;
444 }
445 }
446 else { // no 'sups' line above closest line
447 fICode = kICODE7; //a gauche de la ligne fragment, Z est alors un Zmin et le plus probable
448 y2 = y1;
449 ix2 = 1;
450 y1 = -y1;
451 ix1 = -1;
452 }
453 }
454 else {
455 fICode = kICODE8; // Z indetermine ou (x,y) hors limites
456 }
457 }
458 else if (kinf > -1) {
459 //cout <<"/****************** Seule une ligne inferieure a ete trouvee ***********************/" << endl;
460 /*** Sep. fragment ***/
461 if (Zinf == -1) { // 'inf' is IMF line
462 //point is above IMF line. Z = Z of last line in grid, A = -1
463 k = -1;
464 Z = GetZmax();
465 A = -1;
466 Aint = 0;
467 fICode = kICODE6; // au-dessus de la ligne fragment, Z est alors un Zmin
468 }
469 /*** Ligne de crete (Z,A line)***/
470 else {
471 ibif = 3;
472 k = kinf;
473 Z = Zinf;
474 A = Ainf;
475 Aint = Ainf;
476 yy = dinf;
477 y2 = 0.5 * winf;
478 if (kinfi > -1) {
479 y1 = dinfi - dinf;
480 if (Zinfi == Zinf) {
481 ibif = 1;
482 ix1 = Ainfi - Ainf;
483 Double_t x2 = -y1 / ix1 / 2.;
484 y2 = TMath::Max(y2, x2);
485 ix2 = 1;
486 y1 /= -2.;
487 }
488 else {
489 y1 = -TMath::Min(y2, y1 / 2.);
490 ix1 = -1;
491 ix2 = 1;
492 }
493 }
494 else {
495 y1 = -y2;
496 ix1 = -1;
497 ix2 = 1;
498 }
499 fICode = kICODE7; // a gauche de la ligne fragment, Z est alors un Zmin et le plus probable
500 }
501 }
502 /*****************Aucune ligne n'a ete trouvee*********************************/
503 else {
504 fICode = kICODE8; // Z indetermine ou (x,y) hors limites
505 }
506 /****************Test des bornes********************************************/
507 if (k > -1 && fICode == kICODE0) {
508 if (yy > y2)
509 fICode = kICODE4; // Z ok, masse hors limite superieure ou egale a A
510 }
511 if (k > -1 && (fICode == kICODE0 || fICode == kICODE7)) {
512 if (yy < y1)
513 fICode = kICODE5; // Z ok, masse hors limite inferieure ou egale a A
514 }
515 if (fICode == kICODE4 || fICode == kICODE5) {
516 A = -1;
517 Aint = 0;
518 }
519
520 /****************Interpolation de la masse: da = f*log(1+b*dy)********************/
521 if (fICode == kICODE0 || (fICode == kICODE7 && yy <= y2)) {
522 Double_t deltaA = 0.;
523 Bool_t i = kFALSE;
524 Double_t dt, dist = y1 * y2;
525 dt = 0.;
526 if (ix2 == -ix1) { //dA1 = dA2
527 if (dist != 0) {
528 dt = -(y1 + y2) / dist;
529 i = kTRUE;
530 }
531 /*else
532 Warning("IdentZA","%s : cannot calculate interpolated mass (dist=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
533 GetName(), dist, Z, Aint, fICode);*/
534 }
535 else if (ix2 == -ix1 * 2) { // dA2 = 2*dA1
536 Double_t tmp = y1 * y1 - 4. * dist;
537 if (tmp > 0. && dist != 0) {
538 dt = -(y1 + 2. * y2 -
539 TMath::Sqrt(tmp)) / dist / 2.;
540 i = kTRUE;
541 }
542 /*else
543 Warning("IdentZA","%s : cannot calculate interpolated mass (y1*y1-4*dist=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
544 GetName(), tmp, Z, Aint, fICode);*/
545 }
546 else if (ix1 == -ix2 * 2) { // dA1 = 2*dA2
547 Double_t tmp = y2 * y2 - 4. * dist;
548 if (tmp > 0. && dist != 0) {
549 dt = -(y2 + 2. * y1 +
550 TMath::Sqrt(tmp)) / dist / 2.;
551 i = kTRUE;
552 }
553 /*else
554 Warning("IdentZA","%s : cannot calculate interpolated mass (y2*y2-4*dist=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
555 GetName(), tmp, Z, Aint, fICode);*/
556 }
557 if (i) {
558 dist = dt * y2;
559 if (TMath::Abs(dist) < 0.001) {
560 if (y2 != 0)
561 deltaA = yy * ix2 / y2 / 2.;
562 /*else
563 Warning("IdentZA","%s : cannot calculate interpolated mass (y2=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
564 GetName(), y2, Z, Aint, fICode);*/
565 }
566 else {
567 if (dist > -1. && dt * yy > -1.)
568 deltaA = ix2 / 2. / TMath::Log(1. + dist) * TMath::Log(1. + dt * yy);
569 /*else
570 Warning("IdentZA","%s : cannot calculate interpolated mass (dist=%f dt*yy=%f), Areal will equal Aint (Z=%d Aint=%d fICode=%d)",
571 GetName(), dist, dt*yy, Z, Aint, fICode);*/
572 }
573 A += deltaA;
574 }
575 }
576 /***************D'autres masses sont-elles possibles ?*************************/
577 if (fICode == kICODE0) {
578 //cout << "icode = 0, ibif = " << ibif << endl;
579 /***Masse superieure***/
580 if (ibif == 1 || ibif == 3) {
581 //We look at next line in the complete list of lines, after the closest line.
582 //If it has the same Z as the closest line, but was excluded from research for closest line
583 //because the point lies outside the endpoints, there remains a doubt about the mass:
584 //on rajoute 1 a fICode, effectivement on le met = kICODE1
585 Int_t idx = fIdxClosest;
586 if (idx > -1 && ++idx < GetNumberOfIdentifiers()) {
587 KVIDCsIRLLine* nextline =
589 if (nextline->GetZ() == Z
590 && !nextline->IsBetweenEndPoints(x, y, "x")) {
591 fICode++; // Z ok, mais les masses superieures a A sont possibles
592 //cout <<"//on rajoute 1 a fICode, effectivement on le met = kICODE1" << endl;
593 }
594 }
595 }
596 /***Masse inferieure***/
597 if (ibif == 2 || ibif == 3) {
598 //We look at line below the closest line in the complete list of lines.
599 //If it has the same Z as the closest line, but was excluded from research for closest line
600 //because the point lies outside the endpoints, there remains a doubt about the mass:
601 //on rajoute 2 a fICode, so it can be = kICODE2 or kICODE3
602 Int_t idx = fIdxClosest;
603 if (idx > -1 && --idx >= 0) {
604 KVIDCsIRLLine* nextline =
606 if (nextline->GetZ() == Z
607 && !nextline->IsBetweenEndPoints(x, y, "x")) {
608 fICode += 2;
609 //cout << "//on rajoute 2 a fICode, so it can be = kICODE2 or kICODE3" << endl;
610 }
611 }
612 }
613 }
614
615 // cout << "Z = " << Z << " A = " << A << " icode = " << fICode << endl;
616}
617
618
619
620
629
631{
632 // General initialisation method for CsI R-L identification grid.
633 // This method MUST be called once before using the grid for identifications.
634 // The ID lines are sorted.
635 // The natural line widths of all ID lines are calculated.
636 // The line with the largest Z (Zmax line) is found.
637 // IMF & Gamma line pointers are initialised
638
639 // if grid has already been used for identification, IMF_line will be in identifiers list.
640 TObject* imfline = fIdentifiers.FindObject("IMF_line");
641 if (imfline) fIdentifiers.Remove(imfline); // remove to avoid problems with CalculateLineWidths
643 GammaLine = (KVIDLine*)GetCut("gamma_line");
644 if (!GammaLine) {
645 GetName();
646 Error("Initialize", "%s: Cut 'gamma_line' not found in grid. Not a valid CsI R-L grid. Listing existing cuts:",
647 GetName());
648 GetCuts()->ls();
649 }
650 IMFLine = (KVIDLine*)GetCut("IMF_line");
651 if (!IMFLine) {
652 Error("Initialize", "%s: Cut 'IMF_line' not found in grid. Not a valid CsI R-L grid. Listing existing cuts:",
653 GetName());
654 GetCuts()->ls();
655 }
656}
657
658
659
660
672
674{
675 // Called after reading a grid from an ascii file.
676 // Tries to convert information written by an old version of the class:
677 //
678 //<PARAMETER> Ring min=... ----> <PARAMETER> IDTelescopes=...
679 //<PARAMETER> Ring max=...
680 //<PARAMETER> Mod min=...
681 //<PARAMETER> Mod max=...
682 //
683 // This will fail. The fix is no longer supported. Such files should
684 // no longer exist.
685
687 if (GetParameters()->HasParameter("IDTelescopes")) return;
688
689 Fatal("BackwardsCompatibilityFix",
690 "This fix no longer works. There will be problems.");
691 GammaLine = (KVIDLine*)GetCut("gamma_line");
692 IMFLine = (KVIDLine*)GetCut("IMF_line");
693 if (GammaLine)((KVIDCutLine*)GammaLine)->SetAcceptedDirection("above");
694 if (IMFLine)((KVIDCutLine*)IMFLine)->SetAcceptedDirection("below");
695 SetVarY("CSI-R");
696 SetVarX("CSI-L");
697}
698
699
int Int_t
bool Bool_t
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 index
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
Base class for CsI R-L identification lines (A and Z identification).
Line in ID grid used to delimit regions where no identification is possible.
Definition KVIDCutLine.h:23
Identification grids for CsI R-L (fast-slow) matrices.
Definition KVIDGCsI.h:42
KVIDLine * GammaLine
Definition KVIDGCsI.h:45
virtual ~KVIDGCsI()
Dtor.
Definition KVIDGCsI.cpp:42
KVIDZALine * GetZALine(Int_t z, Int_t a, Int_t &) const
Definition KVIDGCsI.cpp:216
void IdentZA(Double_t x, Double_t y, Int_t &Z, Double_t &A)
Definition KVIDGCsI.cpp:245
KVIDLine * IMFLine
Definition KVIDGCsI.h:44
KVIDZALine * GetZLine(Int_t z, Int_t &) const
Definition KVIDGCsI.cpp:188
virtual void Initialize()
Definition KVIDGCsI.cpp:630
virtual void BackwardsCompatibilityFix()
Definition KVIDGCsI.cpp:673
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
const KVNameValueList * GetParameters() const
Definition KVIDGraph.h:288
virtual void SetVarX(const char *v)
Definition KVIDGraph.h:525
KVList fIdentifiers
list of identification objects
Definition KVIDGraph.h:40
Int_t GetNumberOfIdentifiers() const
Definition KVIDGraph.h:328
const KVList * GetIdentifiers() const
Definition KVIDGraph.h:298
virtual void SetVarY(const char *v)
Definition KVIDGraph.h:529
const Char_t * GetName() const
virtual void BackwardsCompatibilityFix()
virtual void SetInfos(Double_t, Double_t, KVIdentificationResult *) const
loop over KVIDGraph::fInfoZones to set flags in KVIdentificationResult
const KVList * GetCuts() const
Definition KVIDGraph.h:308
KVIDentifier * GetCut(const Char_t *name) const
Definition KVIDGraph.h:280
KVIDentifier * GetIdentifierAt(Int_t index) const
Definition KVIDGraph.h:271
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)
Bool_t IsBetweenEndPoints(Double_t x, Double_t y, const Char_t *axis="") const
Identification grid with lines corresponding to different nuclear isotopes (KVIDZALine)
Definition KVIDZAGrid.h:66
Double_t dsups
Definition KVIDZAGrid.h:88
Int_t Zsups
Definition KVIDZAGrid.h:90
Double_t wsup
Definition KVIDZAGrid.h:89
Int_t Zinf
Definition KVIDZAGrid.h:90
Int_t Ainfi
Definition KVIDZAGrid.h:91
Int_t Ainf
Definition KVIDZAGrid.h:91
virtual void Copy(TObject &) const
Copy this to 'obj'.
Int_t Zsup
Definition KVIDZAGrid.h:90
Int_t Asups
Definition KVIDZAGrid.h:91
Int_t kinfi
Definition KVIDZAGrid.h:87
Int_t Asup
Definition KVIDZAGrid.h:91
Int_t ksup
Definition KVIDZAGrid.h:87
Int_t ksups
used by IdentZA and IdentZ
Definition KVIDZAGrid.h:87
Int_t GetZmax() const
Definition KVIDZAGrid.h:133
Double_t winf
Definition KVIDZAGrid.h:89
Double_t dinf
Definition KVIDZAGrid.h:88
Int_t kinf
Definition KVIDZAGrid.h:87
virtual void Initialize()
Double_t dsup
Definition KVIDZAGrid.h:88
Int_t fIdxClosest
index of closest line in main list fIdentifiers
Definition KVIDZAGrid.h:84
Int_t Zinfi
Definition KVIDZAGrid.h:90
Int_t fICode
code de retour
Definition KVIDZAGrid.h:85
Double_t dinfi
Definition KVIDZAGrid.h:88
virtual Bool_t FindFourEmbracingLines(Double_t x, Double_t y, const Char_t *position)
Int_t Aint
mass of line used to identify particle
Definition KVIDZAGrid.h:93
Base class for identification ridge lines corresponding to different nuclear species.
Definition KVIDZALine.h:33
virtual Int_t GetZ() const
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:28
virtual TObject * FindObject(const char *name) const
virtual void AddLast(TObject *obj)
virtual TObject * At(Int_t idx) const
virtual TObject * Remove(TObject *obj)
Remove object from list.
void ls(Option_t *option="") const override
Int_t GetN() const
Double_t * GetX() const
const char * GetName() const override
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(const KVList &l, KVIDLine *i)
Definition KVIDGCsI.cpp:81
TLine l
TArc a
ClassImp(TPyArg)