KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVZALineFinder.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Tue Dec 4 09:49:15 2012
2 //Author: dgruyer
3 
4 #include "KVZALineFinder.h"
6 #include "KVSpiderLine.h"
7 #include "TSystem.h"
8 #include "KVIDZALine.h"
9 #include "KVIDGridEditor.h"
10 
12 
13 
14 
18 {
19  // constructor
20  fGrid = gg;
21  fHisto = hh;
22 
23  fGeneratedGrid = 0;
24  fLinearHisto = 0;
25  // fCanvas = 0;
26 
27  fPoints = new TGraph;
28  fPoints->SetMarkerStyle(21);
29  fNPoints = 0;
30  fLines = new TList;
31  SetAList("1,4,7,9,11,12,14,16,19,21,23,25,27,29,31,34,37,40,41");
32  fBinsByZ = 60;
33 }
34 
35 
36 
37 
44 
46 {
47  // Copy constructor
48  // This ctor is used to make a copy of an existing object (for example
49  // when a method returns an object), and it is always a good idea to
50  // implement it.
51  // If your class allocates memory in its constructor(s) then it is ESSENTIAL :-)
52 
53  obj.Copy(*this);
54 }
55 
56 
57 
60 
62 {
63  // Destructor
64 }
65 
66 
67 
75 
76 void KVZALineFinder::Copy(TObject& obj) const
77 {
78  // This method copies the current state of 'this' object into 'obj'
79  // You should add here any member variables, for example:
80  // (supposing a member variable KVZALineFinder::fToto)
81  // CastedObj.fToto = fToto;
82  // or
83  // CastedObj.SetToto( GetToto() );
84 
85  KVBase::Copy(obj);
86  //KVZALineFinder& CastedObj = (KVZALineFinder&)obj;
87 }
88 
89 
90 
92 
93 void KVZALineFinder::SetAList(const char* Alist)
94 {
95  Int_t zmax = (Int_t)(((KVIDentifier*)fGrid->GetIdentifiers()->Last())->GetPID() + 1);
96  KVNucleus nuc;
97 
98  Int_t Z = 1;
99  fAList.clear();
100  KVNumberList Al(Alist);
101  Al.Begin();
102  while (!Al.End()) {
103  Int_t A = Al.Next();
104  if (A <= 0) A = nuc.GetAFromZ(Z, fGrid->GetMassFormula());
105  fAList.push_back(A);
106  Z++;
107  }
108 
109  for (int i = Z + 1; i <= zmax; i++) fAList.push_back(nuc.GetAFromZ(i, fGrid->GetMassFormula()));
110 }
111 
112 
113 
115 
117 {
118  Double_t zmin = ((KVIDentifier*)fGrid->GetIdentifiers()->First())->GetPID() - 1.0;
119  Double_t zmax = ((KVIDentifier*)fGrid->GetIdentifiers()->Last())->GetPID() + 1.0;
120  Int_t zbins = (Int_t)(zmax - zmin) * nZbin;
121 
122  Double_t emin = fHisto->GetXaxis()->GetXmin();
123  Double_t emax = fHisto->GetXaxis()->GetXmax();
124  Int_t ebins = fHisto->GetXaxis()->GetNbins();
125 
126  if (fLinearHisto) delete fLinearHisto;
127  fLinearHisto = new TH2F("fLinearHisto", "fLinearHisto", ebins, emin, emax, zbins, zmin, zmax);
128 
129  // if(fConvertHisto) delete fConvertHisto ;
130  fConvertHisto = new TH2F("fConvertHisto", "fConvertHisto", ebins, emin, emax, zbins, zmin, zmax);
131 
133 
134  Int_t events_read = 0;
135  for (int i = 1; i <= fHisto->GetNbinsX(); i++) {
136  for (int j = 1; j <= fHisto->GetNbinsY(); j++) {
137  Stat_t poids = fHisto->GetBinContent(i, j);
138  if (poids == 0) continue;
139 
140  Axis_t x0 = fHisto->GetXaxis()->GetBinCenter(i);
141  Axis_t y0 = fHisto->GetYaxis()->GetBinCenter(j);
142  Axis_t wx = fHisto->GetXaxis()->GetBinWidth(i);
143  Axis_t wy = fHisto->GetYaxis()->GetBinWidth(j);
144 
145  Double_t x, y;
146  Int_t kmax = (Int_t) TMath::Min(20., poids);
147  Double_t weight = (kmax == 20 ? poids / 20. : 1.);
148  for (int k = 0; k < kmax; k++) {
149  x = gRandom->Uniform(x0 - .5 * wx, x0 + .5 * wx);
150  y = gRandom->Uniform(y0 - .5 * wy, y0 + .5 * wy);
151  if (fGrid->IsIdentifiable(x, y)) {
152  fGrid->Identify(x, y, idr);
153  Float_t PID = idr->PID;
154  if (idr->Aident) PID = (idr->Z + 0.1 * (idr->PID - 2. * idr->Z));
155  fLinearHisto->Fill(x, PID, weight);
156  // fConvertHisto->SetPoint(fConvertHisto->GetN(),x0, nuc.GetPID(), y0);
157  }
158  events_read += (Int_t) poids;
159  IncrementLinear((Float_t) events_read);
161  }
162  if (fGrid->IsIdentifiable(x0, y0)) {
163  fGrid->Identify(x0, y0, idr);
164  Float_t PID = idr->PID;
165  if (idr->Aident) PID = (idr->Z + 0.1 * (idr->PID - 2. * idr->Z));
166  Int_t ibin = fConvertHisto->FindBin(x0, PID);
167  fConvertHisto->SetBinContent(ibin, y0);
168  }
169  }
170  }
171 
172  delete idr;
173 
174  // new KVCanvas;
175  // fConvertHisto->Draw("colz");
176  // cout << fConvertHisto->GetN() << endl;
177 
178  return fLinearHisto;
179 }
180 
181 
182 
184 
186 {
187  fLinearHisto->SetAxisRange(zz - 0.5, zz + 0.5, "Y");
188 
189  KVIDLine* line = (KVIDLine*)fGrid->GetIdentifier(zz, 2 * zz + 1); // A=2*zz+1 : dummy, A is ignored in this case
190  if (!line) {
191  int i = 1;
192  while (!(line = (KVIDLine*)fGrid->GetIdentifier(zz + i, 2 * zz + 1))) i++;
193  }
194  if (!line) return;
195 
196  Double_t lX, lY;
197  line->GetStartPoint(lX, lY);
198  Int_t xbmin = 1;
199  line->GetEndPoint(lX, lY);
200  Int_t xbmax = fLinearHisto->GetXaxis()->FindBin(lX);
201  Int_t width = (Int_t)((xbmax - xbmin) * 1.0 / 100.);
202 
203 
204  Int_t widthmax = (Int_t)((xbmax - xbmin) * 1.0 / 30.);
205 
206  KVIDZALine* TheLine = 0;
207  TheLine = (KVIDZALine*)((KVIDZAGrid*)fGeneratedGrid)->NewLine("ID");
208  TheLine->SetZ(zz);
209  TheLine->SetA(fAList.at(zz - 1));
210 
211  TH1* projey = 0;
212  Int_t i = 0;
213  Double_t lasty = zz;
214  for (int xx = xbmin; xx < xbmax; xx += width) {
215  projey = fLinearHisto->ProjectionY("ProjectionAfterLin", TMath::Max(xx - width / 2, xbmin), xx + width / 2);
216  Double_t xline = fLinearHisto->GetBinCenter(xx);
217  Double_t yline = projey->GetMean();
218  if ((yline > zz - 0.5) && (yline < zz + 0.5)) {
219  if (i == 0) lasty = yline;
220  if (TMath::Abs(yline - lasty) < 0.2) {
221  TheLine->SetPoint(i, xline, yline);
222  lasty = yline;
223  i++;
224  if (width < widthmax) width *= 1.2;
225  }
226  }
227  delete projey;
228  }
229 
230  if (TheLine->GetN() > 5) fGeneratedGrid->Add("ID", TheLine);
231 }
232 
233 
234 
236 
238 {
239  fLinearHisto->SetAxisRange(zz - 0.5, zz + 0.5, "Y");
240 
241  KVIDLine* line = (KVIDLine*)fGrid->GetIdentifier(zz, 2 * zz + 1); // A=2*zz+1 : dummy, A is ignored in this case
242  if (!line) {
243  int i = 1;
244  while (!(line = (KVIDLine*)fGrid->GetIdentifier(zz + i, 2 * zz + 1))) i++;
245  }
246  if (!line) return;
247 
248  Double_t lX, lY;
249  line->GetStartPoint(lX, lY);
250  Int_t xbmin = 0;//fLinearHisto->GetYaxis()->FindBin(lX);
251  line->GetEndPoint(lX, lY);
252  Int_t xbmax = fLinearHisto->GetXaxis()->FindBin(lX);
253 
254  // create lines
255  TList Lines;
256  KVSpiderLine* tmp = 0;
257 
258  fLinearHisto->SetAxisRange(fLinearHisto->GetXaxis()->GetBinCenter(50), lX, "X"); //fLinearHisto->GetXaxis()->GetXmax(),"X");
259  TH1* tmph = fLinearHisto->ProjectionX(Form("tmph%d", zz));
260  Int_t startBin = (Int_t)(tmph->GetMaximumBin() * 0.95);
261  delete tmph;
262 
263  TH1* projey = 0;
264  if (startBin) {
265  projey = fLinearHisto->ProjectionY("ProjectionAfterLin", startBin - width * 3, startBin + width * 3);
266  int nfound = fSpectrum.Search(projey, 0.05, "goff", 0.0001);
267  Info("FindALine", "%d peack found...", nfound);
268 #if ROOT_VERSION_CODE > ROOT_VERSION(5,99,01)
269  Double_t* xpeaks = fSpectrum.GetPositionX();
270  Double_t* ypeaks = fSpectrum.GetPositionY();
271 #else
272  Float_t* xpeaks = fSpectrum.GetPositionX();
273  Float_t* ypeaks = fSpectrum.GetPositionY();
274 #endif
275  for (int p = 0; p < nfound; p++) {
276  if (p > 8) break;
277  if (ypeaks[p] < 10) continue;
278  Double_t xline = fLinearHisto->GetXaxis()->GetBinCenter(startBin);
279  Double_t yline = xpeaks[p];
280  KVSpiderLine* tmp = 0;
281  TIter next(&Lines);
282  while ((tmp = (KVSpiderLine*)next())) {
283  Info("FindALine", "line found but I don't know why...");
284  if (TMath::Abs(tmp->GetY() - yline) < 0.05) continue;
285  }
286  tmp = new KVSpiderLine(zz, -1);
287  Lines.AddLast(tmp);
288  tmp->AddPoint(xline, yline);
289  fPoints->SetPoint(fNPoints, xline, yline);
290  fNPoints++;
291  }
292  if (projey) delete projey;
293  }
294  else Error("FindALine", "not starting bin indicated...");
295  SortLines(&Lines);
296 
297  Int_t nLines = Lines.GetSize();
298  tmp = 0;
299  for (int xx = startBin - width; xx > xbmin; xx -= width) {
300  projey = fLinearHisto->ProjectionY("ProjectionAfterLin", xx - width / 2, xx + width / 2);
301  int nfound = fSpectrum.Search(projey, 0.05, "goff", 0.02);
302 #if ROOT_VERSION_CODE > ROOT_VERSION(5,99,01)
303  Double_t* xpeaks = fSpectrum.GetPositionX();
304  Double_t* ypeaks = fSpectrum.GetPositionY();
305 #else
306  Float_t* xpeaks = fSpectrum.GetPositionX();
307  Float_t* ypeaks = fSpectrum.GetPositionY();
308 #endif
309  for (int p = 0; p < nfound; p++) {
310  if (p >= nLines + 1) continue;
311  if (ypeaks[p] < 5) continue;
312  Double_t xline = fLinearHisto->GetXaxis()->GetBinCenter(xx);
313  Double_t yline = xpeaks[p];
314  KVSpiderLine* tmp = 0;
315  TIter next(&Lines);
316  while ((tmp = (KVSpiderLine*)next())) {
317  if ((TMath::Abs(tmp->GetY() - yline) < 0.05)) break;
318  }
319  if (tmp) {
320  if ((TMath::Abs(tmp->GetX() - xline) < 10 * width)) tmp->AddPoint(xline, yline);
321  }
322  }
323  if (projey) delete projey;
324  }
325 
326  TIter nextli(&Lines);
327  while ((tmp = (KVSpiderLine*)nextli()))tmp->Sort(true);
328 
329  tmp = 0;
330  for (int xx = startBin + width; xx <= xbmax - width / 2; xx += width) {
331  projey = fLinearHisto->ProjectionY("ProjectionAfterLin", xx - width / 2, xx + width / 2);
332  int nfound = fSpectrum.Search(projey, 0.05, "goff", 0.02);
333 #if ROOT_VERSION_CODE > ROOT_VERSION(5,99,01)
334  Double_t* xpeaks = fSpectrum.GetPositionX();
335  Double_t* ypeaks = fSpectrum.GetPositionY();
336 #else
337  Float_t* xpeaks = fSpectrum.GetPositionX();
338  Float_t* ypeaks = fSpectrum.GetPositionY();
339 #endif
340  for (int p = 0; p < nfound; p++) {
341  if (p >= nLines + 1) continue;
342  if (ypeaks[p] < 5) continue;
343  Double_t xline = fLinearHisto->GetXaxis()->GetBinCenter(xx);
344  Double_t yline = xpeaks[p];
345  KVSpiderLine* tmp = 0;
346  TIter next(&Lines);
347  while ((tmp = (KVSpiderLine*)next())) {
348  if (TMath::Abs(tmp->GetY() - yline) < 0.05) break;
349  }
350  if (tmp) {
351  if ((TMath::Abs(tmp->GetX() - xline) < 10 * width)) tmp->AddPoint(xline, yline);
352  }
353  }
354  if (projey) delete projey;
355  }
356 
357  fLines->AddAll(&Lines);
358 }
359 
360 
361 
363 
365 {
366  int nn = Lines->GetSize();
367  if (!nn) return;
368 
369  Int_t zz = ((KVSpiderLine*)Lines->At(0))->GetZ();
370 
371  double* yy = new double[nn];
372  Int_t* ii = new int[nn];
373  for (int i = 0; i < nn; i++) yy[i] = ((KVSpiderLine*)Lines->At(i))->GetY(0);
374 
375  KVNucleus nuc;
376  TMath::Sort(nn, yy, ii, kFALSE);
377  Int_t iMostProb = 0;
378  for (int i = 0; i < nn; i++) {
379  if (ii[i] == 0) iMostProb = i;
380  }
381 
382  Int_t aMostProb = fAList.at(zz - 1);
383  // if(zz<=fAList.size()) aMostProb = fAList.at(zz-1);
384  // else aMostProb = 2*zz+1;
385 
386  for (int i = 0; i < nn; i++) {
387  Int_t aa = aMostProb - (iMostProb - i);
388  nuc.SetZandA(zz, aa);
389  if ((nuc.GetLifeTime() < 1.e-06)) aa -= 1;
390  ((KVSpiderLine*)Lines->At(ii[i]))->SetA(aa);
391  }
392 
393  return;
394 }
395 
396 
397 
398 
400 
402 {
404  fGeneratedGrid->SetName(Form("%s_Up", fGrid->GetName()));
408 
409  KVIDZALine* TheLine = 0;
410  KVIDZALine* TheZLine = 0;
411  KVSpiderLine* spline = 0;
412  TIter next_line(fLines);
413  while ((spline = (KVSpiderLine*)next_line())) { // generate KVLines from KVSpiderLines
414  spline->Sort();
415  if ((spline->GetN() > 5)) {
416  TheLine = (KVIDZALine*)((KVIDZAGrid*)fGeneratedGrid)->NewLine("ID");
417  TheLine->SetZ(spline->GetZ());
418  TheLine->SetA(spline->GetA());
419  for (int i = 0; i < spline->GetN(); i++) {
420  Int_t ibin = fConvertHisto->FindBin(spline->GetX(i), spline->GetY(i));
422  Int_t di = 1;
423  while (!y) {
424  y = fConvertHisto->GetBinContent(ibin + di);
425  di++;
426  }
427  if ((i > spline->GetN() - 3) && (std::abs(y - TheLine->GetY()[TheLine->GetN() - 1]) > 20)) continue;
428  TheLine->SetPoint(i, spline->GetX(i), y);//spline->GetY(i));
429  }
430 
431  TheZLine = (KVIDZALine*)fGrid->GetIdentifier(TheLine->GetZ(), 200);
432  if (!TheZLine) continue;
433  Double_t x0 = TheLine->GetX()[0];
434  Double_t x0ref = TheZLine->GetX()[0];
435  Double_t x1ref = TheZLine->GetX()[TheZLine->GetN() - 1];
436  if ((x0 > 2.*x0ref) && (x0 < 0.15 * x1ref)) {
437  Double_t dy = (TheLine->GetY()[0]) - (TheZLine->Eval(x0));
438  for (int i = 0; i <= TheZLine->GetN(); i++) {
439  Double_t x = TheZLine->GetX()[i];
440  Double_t y = dy + TheZLine->GetY()[i];
441  if (x >= x0) break;
442  TheLine->SetPoint(TheLine->GetN(), x, y);
443  }
444  TheLine->Sort(&TGraph::CompareX, kTRUE);
445  }
446 
447  Double_t x1 = TheLine->GetX()[TheLine->GetN() - 1];
448  if ((x1 < x1ref)) {
449  Double_t dy = (TheLine->GetY()[TheLine->GetN() - 1]) - (TheZLine->Eval(x1));
450  for (int i = TheZLine->GetN() - 1; i > 0; i--) {
451  Double_t x = TheZLine->GetX()[i];
452  Double_t y = dy + TheZLine->GetY()[i];
453  if (x <= x1) break;
454  TheLine->SetPoint(TheLine->GetN(), x, y);
455  }
456  TheLine->Sort(&TGraph::CompareX, kTRUE);
457  }
458 
459  fGeneratedGrid->Add("ID", TheLine);
460  fLines->Remove(spline);
461  }
462  }
463 
464  // return;
465 
466  // TheLine = 0;
467  // spline = 0;
468  // TIter next(fLines);
469  // while((spline = (KVSpiderLine*)next())) // scan rejected lines
470  // {
471  // Int_t zl = spline->GetZ();
472  // if(zl<4) continue;
473  // Int_t aMostProb = 0;
474  // if(zl<=fAList.size()) aMostProb = fAList.at(zl-1);
475  // else aMostProb = 2*zl+1;
476  // int index = 0;
477  // KVIDZALine* oldLine = 0;
478  // oldLine=(fGeneratedGrid->GetZALine(spline->GetZ(),aMostProb,index));
479  // if(!oldLine)
480  // {
481  // continue;
482  // }
483  // TheLine = (KVIDZALine*)((KVIDZAGrid*)fGeneratedGrid)->NewLine("ID");
484  // TheLine->SetZ(spline->GetZ());
485  // TheLine->SetA(spline->GetA());
486 
487  // Int_t ibin = fConvertHisto->FindBin(spline->GetX(),spline->GetY());
488  // Double_t y = fConvertHisto->GetBinContent(ibin);
489  // Int_t di = 1;
490  // while(!y)
491  // {
492  // y = fConvertHisto->GetBinContent(ibin+di);
493  // di++;
494  // }
495  // Double_t dy = y - oldLine->Eval(spline->GetX());
496  // for(int i=0; i<oldLine->GetN(); i++)
497  // {
498  // TheLine->SetPoint(i, (oldLine->GetX()[i]), (oldLine->GetY()[i])+dy);
500  // }
502  // }
503 
504 
505 }
506 
507 
508 
509 
511 
513 {
517 }
518 
519 
520 
521 
523 
525 {
526  Info("ProcessIdentification", "histo linearization using Z grid...");
528  Info("ProcessIdentification", "linearization done!");
529 
530  if (zmin < 0) zmin = ((KVIDentifier*)fGrid->GetIdentifiers()->First())->GetZ();
531  if (zmax < 0) zmax = ((KVIDentifier*)fGrid->GetIdentifiers()->Last())->GetZ();
532 
533  int ww = 10;
534  for (int z = zmin; z <= zmax; z++) {
535  if (z > 4) ww = 20;
536  if (z > 6) ww = 30;
537  if (z > 10) ww = 40;
538  if (z > 12) ww = 50;
539  FindALine(z, ww);
540  Increment(z - zmin);
542  }
543  MakeGrid();
544 
545  Double_t zmGrid = ((KVIDentifier*)fGrid->GetIdentifiers()->Last())->GetPID();
546  if (zmax >= zmGrid) return;
547 
548  for (int z = zmax + 1; z <= zmGrid; z++) {
549  FindZLine(z);
550  Increment(z - zmin);
552  }
553 
554 }
555 
556 
557 
558 
559 
560 
561 
int Int_t
KVIDGridEditor * gIDGridEditor
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
const Bool_t kFALSE
double Axis_t
double Double_t
double Stat_t
float Float_t
const Bool_t kTRUE
include TDocParser_001 C image html pict1_TDocParser_001 png width
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
Base class for KaliVeda framework.
Definition: KVBase.h:141
virtual void Copy(TObject &) const
Make a copy of this object.
Definition: KVBase.cpp:394
void Add(TString, KVIDentifier *)
Definition: KVIDGraph.cpp:838
virtual void SetVarX(const char *v)
Definition: KVIDGraph.h:524
const TList * GetIDTelescopes() const
Definition: KVIDGraph.h:244
virtual void SetName(const char *name)
Definition: KVIDGraph.h:139
Int_t GetMassFormula() const
Definition: KVIDGraph.h:437
KVIDentifier * GetIdentifier(Int_t Z, Int_t A) const
Definition: KVIDGraph.cpp:310
virtual Bool_t IsIdentifiable(Double_t, Double_t, TString *rejected_by=nullptr) const
Definition: KVIDGraph.cpp:1269
virtual void SetVarY(const char *v)
Definition: KVIDGraph.h:528
const Char_t * GetName() const
Definition: KVIDGraph.cpp:1332
void AddIDTelescopes(const TList *)
Associate this graph with all ID telescopes in list.
Definition: KVIDGraph.cpp:1516
const KVList * GetIdentifiers() const
Definition: KVIDGraph.h:297
void SetGrid(KVIDGraph *gg, Bool_t histo=true)
void StartViewer()
Close();.
void SetHisto(TH2 *hh)
Base class for lines/cuts used for particle identification in 2D data maps.
Definition: KVIDLine.h:142
Identification grid with lines corresponding to different nuclear isotopes (KVIDZALine)
Definition: KVIDZAGrid.h:65
virtual void Identify(Double_t x, Double_t y, KVIdentificationResult *) const
Base class for identification ridge lines corresponding to different nuclear species.
Definition: KVIDZALine.h:32
Base class for graphical cuts used in particle identification.
Definition: KVIDentifier.h:27
virtual Int_t GetZ() const
Definition: KVIDentifier.h:78
virtual void SetA(Int_t atnum)
Definition: KVIDentifier.h:87
virtual void SetZ(Int_t ztnum)
Definition: KVIDentifier.h:82
Full result of one attempted particle identification.
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 Z
Z of particle found (if Zident==kTRUE)
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
void SetZandA(Int_t z, Int_t a)
Set atomic number and mass number.
Definition: KVNucleus.cpp:721
static Int_t GetAFromZ(Double_t, Char_t mt)
Definition: KVNucleus.cpp:564
Double_t GetLifeTime(Int_t z=-1, Int_t a=-1) const
Definition: KVNucleus.cpp:1040
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
Bool_t End(void) const
Definition: KVNumberList.h:197
void Begin(void) const
Int_t Next(void) const
virtual TObject * Last() const
virtual TObject * First() const
Part of Spider Identification.
Definition: KVSpiderLine.h:18
void Sort(bool ascending_=true)
double GetX(int n_) const
double GetY(int n_) const
bool AddPoint(double x_, double y_, bool test_=false, int n_=-1)
(try to) find mass lines from charge lines
void FindALine(Int_t zz, Int_t width=10)
TH2 * LinearizeHisto(Int_t nZbin=40)
void FindZLine(Int_t zz)
void IncrementLinear(Float_t x)
void Copy(TObject &) const
std::vector< int > fAList
virtual ~KVZALineFinder()
Destructor.
TSpectrum fSpectrum
void SetAList(const char *Alist)
void SortLines(TList *Lines)
TGraph * fPoints
void Increment(Float_t x)
void ProcessIdentification(Int_t zmin=-1, Int_t zmax=-1)
KVZALineFinder(KVIDZAGrid *gg, TH2 *hh)
constructor
KVIDZAGrid * fGeneratedGrid
KVIDZAGrid * fGrid
virtual void SetMarkerStyle(Style_t mstyle=1)
virtual Int_t FindBin(const char *label)
virtual Double_t GetBinCenter(Int_t bin) const
Double_t GetXmax() const
Double_t GetXmin() const
Int_t GetNbins() const
virtual Double_t GetBinWidth(Int_t bin) const
virtual void AddAll(const TCollection *col)
virtual Int_t GetSize() const
const char * GetVarY() const
const char * GetVarX() const
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
static Bool_t CompareX(const TGraph *gr, Int_t left, Int_t right)
Int_t GetN() const
virtual void Sort(Bool_t(*greater)(const TGraph *, Int_t, Int_t)=&TGraph::CompareX, Bool_t ascending=kTRUE, Int_t low=0, Int_t high=-1111)
Double_t * GetX() const
virtual Double_t Eval(Double_t x, TSpline *spline=nullptr, Option_t *option="") const
Double_t * GetY() const
virtual Double_t GetBinCenter(Int_t bin) const
virtual Int_t GetNbinsY() const
TAxis * GetXaxis()
virtual Double_t GetMean(Int_t axis=1) const
TAxis * GetYaxis()
virtual Int_t GetNbinsX() const
virtual Int_t GetMaximumBin() const
virtual void SetAxisRange(Double_t xmin, Double_t xmax, Option_t *axis="X")
virtual Int_t FindBin(Double_t x, Double_t y=0, Double_t z=0)
TH1D * ProjectionY(const char *name="_py", Int_t firstxbin=0, Int_t lastxbin=-1, Option_t *option="") const
virtual Int_t Fill(const char *namex, const char *namey, Double_t w)
TH1D * ProjectionX(const char *name="_px", Int_t firstybin=0, Int_t lastybin=-1, Option_t *option="") const
virtual Double_t GetBinContent(Int_t bin) const
virtual void SetBinContent(Int_t bin, Double_t content)
TObject * Remove(const TObjLinkPtr_t &lnk)
virtual TObject * At(Int_t idx) const
virtual void AddLast(TObject *obj)
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
virtual Double_t Uniform(Double_t x1, Double_t x2)
virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
Double_t * GetPositionX() const
Double_t * GetPositionY() const
virtual Bool_t ProcessEvents()
TLine * line
Double_t y[n]
Double_t x[n]
Double_t Min(Double_t a, Double_t b)
Double_t Abs(Double_t d)
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Double_t Max(Double_t a, Double_t b)
void spline(double x[], double y[], int n, double yp1, double ypn, double *y2)