KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVItvFinderDialog.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Mon Jan 23 10:03:13 2017
2 //Author: Diego Gruyer
3 
4 #include "KVItvFinderDialog.h"
5 #include "TRootEmbeddedCanvas.h"
6 #include "TStyle.h"
7 #include "TSystem.h"
8 #include "TROOT.h"
9 #include "TGMsgBox.h"
10 #include "TGFileDialog.h"
11 #include "KVTestIDGridDialog.h"
12 #include "KVIdentificationResult.h"
13 #include "KVNameValueListGUI.h"
14 #include <KVMultiGaussIsotopeFit.h>
15 #include "KeySymbols.h"
16 #include <thread>
17 
19 //ClassImp(interval_painter)
20 
21 // Set default values for mass-fit parameters
22 
23 
26  {"Limit range of fit", false},
27  {"PID min for fit", 0.},
28  {"PID max for fit", 10.},
29  {"Minimum probability [%]", 50.},
30  {"Minimum #sigma", 1.e-2},
31  {"Maximum #sigma", 5.e-2}
32 };
33 
34 
35 
39 
41 {
42  // remove painter from list and modify the 'left_painter' and 'right_painter' references
43  // in any adjacent painters/intervals, then delete painter
44 
45  std::unique_ptr<KVPIDIntervalPainter> _p(p);
46  auto pleft = p->get_left_interval();
47  auto pright = p->get_right_interval();
48  if (pleft) pleft->set_right_interval(pright);
49  if (pright) pright->set_left_interval(pleft);
50  fItvPaint.Remove(p);
51 }
52 
53 
54 
56 
58 {
59  fGrid = gg;
60  fHisto = hh;
61 
62  fPoints = new TGraph;
64  fNPoints = 0;
65 
66  gStyle->SetOptStat(0);
67  gStyle->SetOptTitle(0);
68 
69  gROOT->ProcessLine(Form("KVItvFinderDialog* _dummy_itv=(KVItvFinderDialog*)%p", this));
70 
71  fMain = new TGTransientFrame(gClient->GetDefaultRoot(), gClient->GetDefaultRoot(), 10, 10);
72  // Here is the recipe for cleanly closing a window
73  // (see https://root-forum.cern.ch/t/error-in-rootx11errorhandler-baddrawable-quot/8095/6)
74  fMain->Connect("CloseWindow()", "KVItvFinderDialog", this, "DoClose()");
77  // afterwards, in the destructor do fMain->CloseWindo(), and in DoClose(), do 'delete this'
78 
79  // Default constructor
80  TGHorizontalFrame* fCanvasFrame = new TGHorizontalFrame(fMain, 627, 7000, kHorizontalFrame);
81  // fCanvasFrame->SetBackgroundColor(fColor);
82 
83 
84  TRootEmbeddedCanvas* fRootEmbeddedCanvas615 = new TRootEmbeddedCanvas(0, fCanvasFrame, 800, 440);
85  // to replace the TCanvas in a TRootEmbeddedCanvas, just delete the original like so:
86  // (see https://root-forum.cern.ch/t/error-in-rootx11errorhandler-baddrawable-quot/8095/6)
87  auto WID = fRootEmbeddedCanvas615->GetCanvasWindowId();
88  delete fRootEmbeddedCanvas615->GetCanvas();
89  fCanvas = new TCanvas("c123", 10, 10, WID);
90  fCanvas->AddExec("toto", "if(_dummy_itv)_dummy_itv->HandleKey();");
91  fRootEmbeddedCanvas615->AdoptCanvas(fCanvas);
92  fPad = fCanvas->cd();
93 // fPad->AddExec("test","printf('coucou\n')");
94  fCanvas->SetRightMargin(0.02);
95  fCanvas->SetTopMargin(0.02);
96  fCanvas->SetLeftMargin(0.08);
97  fCanvas->SetBottomMargin(0.07);
98 
99  fCanvasFrame->AddFrame(fRootEmbeddedCanvas615, new TGLayoutHints(kLHintsExpandX | kLHintsExpandY, 2, 2, 2, 2));
100  fMain->AddFrame(fCanvasFrame, new TGLayoutHints(kLHintsExpandX | kLHintsExpandY, 0, 0, 0, 0));
101 
102  TGVerticalFrame* fControlOscillo = new TGVerticalFrame(fCanvasFrame, 2000, 7000, kVerticalFrame);
103  // fControlOscillo->SetBackgroundColor(fColor);
104 
105 
106  {
107  const char* xpms[] = {
108  "filesaveas.xpm",
109  "profile_t.xpm",
110  "bld_copy.png",
111  "ed_new.png",
112  "refresh2.xpm",
113  "sm_delete.xpm",
114  "h1_t.xpm",
115  "query_new.xpm",
116  "tb_back.xpm",
117  "bld_colorselect.png",
118  "latex.xpm",
119  "move_cursor.png",
120  0
121  };
122  // toolbar tool tip text
123  const char* tips[] = {
124  "Save intervals in current grid",
125  "Find intervals",
126  "Create a new interval set",
127  "Create a new interval",
128  "Update interval lists",
129  "Remove selected intervals",
130  "Multigauss fit to isotopes in interval set",
131  "Set parameters for fit",
132  "Remove fit from selected interval set",
133  "Test identification",
134  "Set log scale on y axis",
135  "Unzoom the histogram",
136  0
137  };
138  int spacing[] = {
139  5,
140  20,
141  0,
142  0,
143  0,
144  20,
145  50,
146  0,
147  0,
148  50,
149  80,
150  0,
151  0
152  };
153  const char* method[] = {
154  "SaveGrid()",
155  "Identify()",
156  "NewIntervalSet()",
157  "NewInterval()",
158  "UpdateLists()",
159  "RemoveInterval()",
160  "FitIsotopes()",
161  "SetFitParameters()",
162  "RemoveFit()",
163  "TestIdent()",
164  "SetLogy()",
165  "UnzoomHisto()",
166  0
167  };
168  fNbButtons = 0;
169  ToolBarData_t t[50];
170  fToolBar = new TGToolBar(fControlOscillo, 450, 80);
171  int i = 0;
172  while (xpms[i]) {
173  t[i].fPixmap = xpms[i];
174  t[i].fTipText = tips[i];
175  t[i].fStayDown = kFALSE;
176  t[i].fId = i + 1;
177  t[i].fButton = NULL;
178  TGButton* bb = fToolBar->AddButton(fControlOscillo, &t[i], spacing[i]);
179  bb->Connect("Clicked()", "KVItvFinderDialog", this, method[i]);
180  fNbButtons++;
181  i++;
182  }
183  fControlOscillo->AddFrame(fToolBar, new TGLayoutHints(kLHintsTop | kLHintsExpandX));
184  }
185 
186  fIntervalSetListView = new KVListView(interval_set::Class(), fControlOscillo, 450, 200);
188  fIntervalSetListView->SetDataColumn(0, "Z", "GetZ", kTextLeft);
189  fIntervalSetListView->SetDataColumn(1, "PIDs", "GetNPID", kTextCenterX);
190  fIntervalSetListView->SetDataColumn(2, "Masses", "GetListOfMasses", kTextLeft);
191  fIntervalSetListView->Connect("SelectionChanged()", "KVItvFinderDialog", this, "DisplayPIDint()");
192  fIntervalSetListView->SetDoubleClickAction("KVItvFinderDialog", this, "ZoomOnCanvas()");
194  fControlOscillo->AddFrame(fIntervalSetListView, new TGLayoutHints(kLHintsTop | kLHintsExpandX | kLHintsExpandY, 2, 2, 2, 2));
195 
196  fIntervalListView = new KVListView(interval::Class(), fControlOscillo, 450, 180);
198  fIntervalListView->SetDataColumn(0, "Z", "GetZ", kTextLeft);
199  fIntervalListView->SetDataColumn(1, "A", "GetA", kTextCenterX);
200  fIntervalListView->SetDataColumn(2, "min", "GetPIDmin", kTextCenterX);
201  fIntervalListView->SetDataColumn(3, "pid", "GetPID", kTextCenterX);
202  fIntervalListView->SetDataColumn(4, "max", "GetPIDmax", kTextCenterX);
203 
204 
205 
206  {
207  const char* xpms[] = {
208  "arrow_down.xpm",
209  "arrow_up.xpm",
210  0
211  };
212  const char* tips[] = {
213  "Decrease A by one",
214  "Increase A by one",
215  0
216  };
217  int spacing[] = {
218  120,
219  0,
220  0
221  };
222  const char* method[] = {
223  "MassesDown()",
224  "MassesUp()",
225  0
226  };
227  fNbButtons = 0;
228  ToolBarData_t t[50];
229  fToolBar2 = new TGToolBar(fControlOscillo, 450, 80);
230  int i = 0;
231  while (xpms[i]) {
232  t[i].fPixmap = xpms[i];
233  t[i].fTipText = tips[i];
234  t[i].fStayDown = kFALSE;
235  t[i].fId = i + 1;
236  t[i].fButton = NULL;
237  TGButton* bb = fToolBar2->AddButton(fControlOscillo, &t[i], spacing[i]);
238  bb->Connect("Clicked()", "KVItvFinderDialog", this, method[i]);
239  fNbButtons++;
240  i++;
241  }
242  fControlOscillo->AddFrame(fToolBar2, new TGLayoutHints(kLHintsTop | kLHintsExpandX));
243  }
244 
245  // fCurrentView->ActivateSortButtons();
246  fIntervalListView->Connect("SelectionChanged()", "KVItvFinderDialog", this, "SelectionITVChanged()");
247  // fCurrentView->SetDoubleClickAction("FZCustomFrameManager",this,"ChangeParValue()");
249 
250  fControlOscillo->AddFrame(fIntervalListView, new TGLayoutHints(kLHintsTop | kLHintsExpandX, 2, 2, 2, 2));
251 
252  fCanvasFrame->AddFrame(fControlOscillo, new TGLayoutHints(kLHintsExpandY, 0, 0, 0, 0));
253 
254 
255  //layout and display window
256  fMain->MapSubwindows();
258 
259  // position relative to the parent's window
261 
262  fMain->SetWindowName("Masses Identification");
263  fMain->RequestFocus();
264  fMain->MapWindow();
265 
266  fIntervalSetListView->Display(((KVIDZAFromZGrid*)fGrid)->GetIntervalSets());
267  current_interval_set = nullptr;
268  fPad->cd();
269 
270  LinearizeHisto(100);
273  fLinearHisto->Draw("hist");
274 
275  int tmp[30] = {3, 3, 3, 4, 4, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6, 6};
276  for (int ii = 0; ii < 30; ii++) fNpeaks[ii] = tmp[ii];
277 
278  fSig = 0.1;
279  fRat = 0.0001;
280 
281  DrawIntervals();
282 
283 }
284 
285 
286 
287 
289 
291 {
292  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
293  Int_t nSelected = list->GetSize();
294  if (nSelected == 1) {
295  current_interval_set = (interval_set*)list->At(0);
297  }
298  else {
299  current_interval_set = nullptr;
301  }
303 }
304 
305 
306 
308 
310 {
311  fPad->cd();
312  fItvPaint.Execute("HighLight", "0");
313 
315  int zz = ((interval*)fIntervalListView->GetLastSelectedObject())->GetZ();
316  int aa = ((interval*)fIntervalListView->GetLastSelectedObject())->GetA();
317  KVPIDIntervalPainter* painter = (KVPIDIntervalPainter*)fItvPaint.FindObject(Form("%d_%d", zz, aa));
318  if (!painter) Info("SelectionITVChanged", "%d %d not found...", zz, aa);
319  painter->HighLight();
320  }
321  fCanvas->Modified();
322  fCanvas->Update();
323 }
324 
325 
326 
328 
330 {
331  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
332  Int_t nSelected = list->GetSize();
333  if (nSelected == 1) {
334  current_interval_set = (interval_set*)list->At(0);
336  }
337  else {
338  current_interval_set = nullptr;
340  }
341 }
342 
343 
344 
347 
349 {
350  // Display the interval set for a given Z when the user double clicks on it
351 
353  current_interval_set = nullptr;
354  return;
355  }
357  auto zz = current_interval_set->GetZ();
358 
359  fLinearHisto->GetXaxis()->SetRangeUser(zz - 0.5, zz + 0.5);
360 
361  fItvPaint.Execute("SetDisplayLabel", "0");
362 
363  std::unique_ptr<KVSeqCollection> tmp(fItvPaint.GetSubListWithMethod(Form("%d", zz), "GetZ"));
364 
365  tmp->Execute("SetDisplayLabel", "1");
366 
367  // Check - if no mass fit is displayed, we look to see if the grid has a saved
368  // multigauss fit from a previous session, and if so we display it
370  if (fGrid->GetParameters()->HasParameter(Form("MASSFIT_%d", zz))) {
371  KVString massfit = fGrid->GetParameters()->GetTStringValue(Form("MASSFIT_%d", zz));
372  massfit.ReplaceAll(":", "=");
373  KVNameValueList fitparams;
374  fitparams.Set(massfit);
375  KVMultiGaussIsotopeFit fitfunc(zz, fitparams);
376  fitfunc.DrawFitWithGaussians("same");
377  // deactivate intervals in painter
378  tmp->Execute("DeactivateIntervals", "");
379  }
380  }
381 
382  fCanvas->Modified();
383  fCanvas->Update();
384 }
385 
386 
387 
389 
391 {
392  interval_set* itvs = 0;
393  last_drawn_interval = nullptr;
394  TIter it(fGrid->GetIntervalSets());
395  while ((itvs = (interval_set*)it())) {
396  DrawInterval(itvs);
397  }
398 }
399 
400 
401 
403 
405 {
406  fPad->cd();
407  interval* itv = nullptr;
408  // give a different colour to each interval
409  auto nitv = itvs->GetIntervals()->GetEntries();
410  auto cstep = TColor::GetPalette().GetSize() / (nitv + 1);
411  int i = 1;
412  auto deactivate_intervals = fGrid->GetParameters()->HasParameter(Form("MASSFIT_%d", itvs->GetZ()));
413 
414  TIter itt(itvs->GetIntervals());
415  while ((itv = (interval*)itt())) {
417  ++i;
418  if (label) dummy->SetDisplayLabel();
419  if (deactivate_intervals) dummy->DeactivateIntervals();
420  dummy->Draw();
421  dummy->SetCanvas(fCanvas);
422  dummy->Connect("IntMod()", "KVItvFinderDialog", this, "UpdatePIDList()");
423  fItvPaint.Add(dummy);
424  last_drawn_interval = dummy;
425  }
426 }
427 
428 
429 
434 
436 {
437  // empty an interval set, effectively removing it from the interval sets which will be saved with the grid.
438  //
439  // we also remove any previous fits from the grid's parameters
440 
441  std::vector<int> alist;
442  for (int ii = 0; ii < itvs->GetNPID(); ii++) {
443  interval* itv = (interval*)itvs->GetIntervals()->At(ii);
444  KVPIDIntervalPainter* pid = (KVPIDIntervalPainter*)fItvPaint.FindObject(Form("%d_%d", itv->GetZ(), itv->GetA()));
446  alist.push_back(itv->GetA());
447  }
448  itvs->GetIntervals()->Clear();
449  UpdateLists();
450  // remove any fit displayed in pad
451  KVMultiGaussIsotopeFit fitfunc(itvs->GetZ(), alist);
452  fitfunc.UnDraw(fPad);
453  // mop up any stray gaussians (from intervals which have been removed)
455 
456  // remove from grid parameters
457  if (fGrid->GetParameters()->HasParameter(Form("MASSFIT_%d", itvs->GetZ()))) {
458  fGrid->GetParameters()->RemoveParameter(Form("MASSFIT_%d", itvs->GetZ()));
459  KVNumberList zlist(fGrid->GetParameters()->GetStringValue("MASSFITS"));
460  zlist.Remove(itvs->GetZ());
461  fGrid->GetParameters()->SetValue("MASSFITS", zlist.AsString());
462  }
463  fCanvas->Modified();
464  fCanvas->Update();
465 }
466 
467 
468 
470 
472 {
473  Double_t zmin = ((KVIDentifier*)fGrid->GetIdentifiers()->First())->GetPID() - 1.0;
474  Double_t zmax = 0;
475 
476  for (int iz = 1; iz < fGrid->GetIdentifiers()->GetSize() + 1; iz++) {
478  if (tmp && tmp->GetPID() > zmax) zmax = tmp->GetPID();
479  }
480 
481  Int_t zbins = (Int_t)(zmax - zmin) * nbins;
482 
483  fLinearHisto = new TH1F("fLinearHisto", "fLinearHisto", zbins, zmin, zmax);
484 
485  fGrid->SetOnlyZId();
486 
487  // use multi-threading capacities
488  auto available_cpu = WITH_MULTICORE_CPU;
489  Int_t xbins_per_cpu = fHisto->GetNbinsX() / available_cpu;
490 
491  std::vector<std::thread> jobs; // threads to do the work
492 
493  // to clean up copies of grid
494  KVList grid_copies;
495 
496  // do not add copies of grid to ID grid manager
497  auto save_auto_add = KVIDGraph::GetAutoAdd();
498  KVIDGraph::SetAutoAdd(false);
499 
500  std::cout << "Will run " << available_cpu << " threads, each for " << xbins_per_cpu << " bins in X" << std::endl;
501 
502  int nthreads = available_cpu;
503  // join threads
504  std::cout << "Histo linearization using " << nthreads << " threads..." << std::endl;
505 
506  for (int job = 0; job < available_cpu; ++job) {
507  auto imin = 1 + job * xbins_per_cpu;
508  auto imax = (job + 1) * xbins_per_cpu;
509  if (job == available_cpu - 1) imax = fHisto->GetNbinsX();
510 
511  // make new copy of grid
512  auto grid_copy = new KVIDZAFromZGrid(*fGrid);
513  grid_copy->Initialize();
514  grid_copies.Add(grid_copy);
515 
516  // start new thread
517  jobs.push_back(std::thread([ =, &nthreads]() {
518 
520 
521  bool no_mass_id_zone_defined = (grid_copy->GetInfos()->FindObject("MassID") == nullptr);
522 
523  for (int i = imin; i <= imax; ++i) {
524  for (int j = 1; j <= fHisto->GetNbinsY(); j++) {
525  Stat_t poids = fHisto->GetBinContent(i, j);
526  if (poids == 0) continue;
527 
528  Axis_t x0 = fHisto->GetXaxis()->GetBinCenter(i);
529  Axis_t y0 = fHisto->GetYaxis()->GetBinCenter(j);
530  Axis_t wx = fHisto->GetXaxis()->GetBinWidth(i);
531  Axis_t wy = fHisto->GetYaxis()->GetBinWidth(j);
532 
533  if (x0 < 4) continue;
534 
535  Double_t x, y;
536  Int_t kmax = (Int_t) TMath::Min(20., poids);
537  Double_t weight = (kmax == 20 ? poids / 20. : 1.);
538  for (int k = 0; k < kmax; k++) {
539  x = gRandom->Uniform(x0 - .5 * wx, x0 + .5 * wx);
540  y = gRandom->Uniform(y0 - .5 * wy, y0 + .5 * wy);
541  if (grid_copy->IsIdentifiable(x, y)) {
542  idr.Clear();
543  grid_copy->KVIDZAGrid::Identify(x, y, &idr);
544  if (no_mass_id_zone_defined || idr.HasFlag(grid_copy->GetName(), "MassID")) {
545  Float_t PID = idr.PID;
546  fLinearHisto->Fill(PID, weight);
547  }
548  }
549  }
550  }
551  }
552  --nthreads;
553  std::cout << "...remaining threads: " << nthreads << std::endl;
554  }));
555  }
556  for (auto& j : jobs) {
557  if (j.joinable()) j.join();
558  }
559 
560  // reset automatic grid adding to previous state
561  KVIDGraph::SetAutoAdd(save_auto_add);
562 }
563 
564 
565 
568 
570 {
571  // KVBase::OpenContextMenu("Identify(double,double)",this);
572  Identify(0.1, 0.001);
573 }
574 
575 
576 
578 
579 void KVItvFinderDialog::Identify(double sigma, double ratio)
580 {
581  fSig = sigma;
582  fRat = ratio;
583 
584  fPad->cd();
585  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
586  if (!list->GetSize()) {
588  for (int ii = 0; ii < fGrid->GetIntervalSets()->GetSize(); ii++) DrawInterval((interval_set*)fGrid->GetIntervalSets()->At(ii), 0);
589  }
590  else {
591  for (int ii = 0; ii < list->GetSize(); ii++) {
592  interval_set* itvs = (interval_set*) list->At(ii);
593  ProcessIdentification(itvs->GetZ(), itvs->GetZ());
594  DrawInterval(itvs, 0);
595  }
596  }
597 
598  fCanvas->Modified();
599  fCanvas->Update();
600 
601 }
602 
603 
604 
606 
608 {
609  ExportToGrid();
610 
611  TString currentdir(gSystem->ExpandPathName("."));
612 
613  TString fn = fHisto->GetName();
614  fn += ".dat";
615 
616  Int_t ret_code;
617  new TGMsgBox(
618  gClient->GetRoot(),
619  gClient->GetDefaultRoot(),
620  "KVIDGridEditor::SaveCurrentGrid", Form("Do you wat to save the grid here : %s", fn.Data()),
621  kMBIconExclamation, kMBYes | kMBNo, &ret_code
622  );
623 
624  if (ret_code == kMBYes) {
625  fGrid->WriteAsciiFile(Form("%s", fn.Data()));
626  return;
627  }
628 
629  static TString dir(".");
630  const char* filetypes[] = {
631  "ID Grid files", "*.dat",
632  "All files", "*",
633  0, 0
634  };
635  TGFileInfo fi;
636  fi.fFileTypes = filetypes;
637  fi.fIniDir = StrDup(dir);
638  new TGFileDialog(gClient->GetDefaultRoot(), gClient->GetDefaultRoot(), kFDSave, &fi);
639  if (fi.fFilename) {
640  TString filenam(fi.fFilename);
641  if (filenam.Contains("toto")) filenam.ReplaceAll("toto", fGrid->GetName());
642  if (!filenam.Contains('.')) filenam += ".dat";
643  fGrid->WriteAsciiFile(filenam.Data());
644  }
645  dir = fi.fIniDir;
647 
648  fIntervalSetListView->Display(((KVIDZAFromZGrid*)fGrid)->GetIntervalSets());
649  current_interval_set = nullptr;
651 
652  fItvPaint.Clear();
653  DrawIntervals();
654 }
655 
656 
657 
660 
662 {
663  // Write all PID intervals in grid parameters "PIDRANGE", "PIDRANGE%d", etc.
664 
666  KVNumberList pids;
667  interval_set* itvs = 0;
668  TIter npid(fGrid->GetIntervalSets());
669  while ((itvs = (interval_set*)npid())) {
670  if (!itvs->GetNPID()) continue;
671  pids.Add(itvs->GetZ());
672  }
673  fGrid->GetParameters()->SetValue("PIDRANGE", pids.AsString());
674 
675  itvs = 0;
676  TIter next(fGrid->GetIntervalSets());
677  while ((itvs = (interval_set*)next())) {
678  if (!itvs->GetNPID()) continue;
679  KVString par = Form("PIDRANGE%d", itvs->GetZ());
680  KVString val = "";
681  interval* itv = 0;
682  TIter ni(itvs->GetIntervals());
683  while ((itv = (interval*)ni())) {
684  val += Form("%d:%lf,%lf,%lf|", itv->GetA(), itv->GetPIDmin(), itv->GetPID(), itv->GetPIDmax());
685  }
686  val.Remove(val.Length() - 1);
687  fGrid->GetParameters()->SetValue(par.Data(), val.Data());
688  }
689 }
690 
691 
692 
694 
696 {
697  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
698  if (!list->GetSize() || list->GetSize() > 1) {
699  current_interval_set = nullptr;
700  return;
701  }
702 
703  current_interval_set = (interval_set*)list->At(0);
704 
705  fPad->WaitPrimitive("TMarker");
706  auto mm = dynamic_cast<TMarker*>(fPad->GetListOfPrimitives()->Last());
707  assert(mm);
708  double pid = mm->GetX();
709  delete mm;
710 
711  int aa = 0;
712  int iint = 0;
713 
714  // try to guess mass from position (PID)
715  // typical isotope separation in PID is 0.12 (e.g. for carbon isotopes)
716  // assume that PID=z means A=2*Z
717  auto aa_guessstimate = TMath::Nint(current_interval_set->GetZ() * 2 + (pid - current_interval_set->GetZ()) / 0.12);
718 //Info("NewINt","pid = %f guess a=%d",pid,aa_guessstimate);
719  if (!current_interval_set->GetNPID()) {
720  aa = aa_guessstimate;
721  iint = 0;
722  }
723  else if (pid < ((interval*)current_interval_set->GetIntervals()->First())->GetPID()) { // to left of all others
724  aa = ((interval*)current_interval_set->GetIntervals()->First())->GetA() - 1;
725  // NO MORE : use guesstimate as long as it is smaller than A of previously defined interval with largest PID
726  // NOW : use 'aa_guessstimate' only for first interval
727  //Info("NewInt","before first interval, with a=%d",aa+1);
728 // if (aa_guessstimate < aa + 1) aa = aa_guessstimate;
729  iint = 0;
730  }
731  else if (pid > ((interval*)current_interval_set->GetIntervals()->Last())->GetPID()) { // to right of all others
732  aa = ((interval*)current_interval_set->GetIntervals()->Last())->GetA() + 1;
733  // NO MORE : use guesstimate as long as it is larger than A of previously defined interval with smallest PID
734  // NOW : use 'aa_guessstimate' only for first interval
735  //Info("NewInt","after last interval, with a=%d",aa-1);
736 // if (aa_guessstimate > aa - 1) aa = aa_guessstimate;
737  iint = current_interval_set->GetNPID();
738  }
739  else {
740  // look for intervals between which the new one is places
741  for (int ii = 1; ii < current_interval_set->GetNPID(); ii++) {
742  if (pid > ((interval*)current_interval_set->GetIntervals()->At(ii - 1))->GetPID()
743  && pid < ((interval*)current_interval_set->GetIntervals()->At(ii))->GetPID()) {
744  aa = ((interval*)current_interval_set->GetIntervals()->At(ii - 1))->GetA() + 1;
745  // use guesstimate if it is in between masses of adjacent intervals (in terms of PID)
746  if ((aa_guessstimate > (aa - 1)) &&
747  (aa_guessstimate < ((interval*)current_interval_set->GetIntervals()->At(ii))->GetA()))
748  aa = aa_guessstimate;
749  iint = ii;
750  break;
751  }
752  }
753  }
754 
755  interval* itv = new interval(current_interval_set->GetZ(), aa, mm->GetX(), mm->GetX() - 0.05, mm->GetX() + 0.05);
756  current_interval_set->GetIntervals()->AddAt(itv, iint);
757 
758  // find intervals which are now left (smaller mass) and right (higher mass) than this one
759  interval* left_interval = iint > 0 ? (interval*)current_interval_set->GetIntervals()->At(iint - 1) : nullptr;
760  interval* right_interval = iint < current_interval_set->GetIntervals()->GetEntries() - 1 ? (interval*)current_interval_set->GetIntervals()->At(iint + 1) : nullptr;
761  // find the corresponding painters
762  KVPIDIntervalPainter* left_painter{nullptr}, *right_painter{nullptr};
763  TIter next_painter(&fItvPaint);
764  KVPIDIntervalPainter* pidpnt;
765  while ((pidpnt = (KVPIDIntervalPainter*)next_painter())) {
766  if (pidpnt->GetInterval() == left_interval) left_painter = pidpnt;
767  else if (pidpnt->GetInterval() == right_interval) right_painter = pidpnt;
768  }
769 
770  // give a different colour to each interval
771  auto nitv = current_interval_set->GetIntervals()->GetEntries();
772  auto cstep = TColor::GetPalette().GetSize() / (nitv + 1);
773 
774  KVPIDIntervalPainter* dummy = new KVPIDIntervalPainter(itv, fLinearHisto, TColor::GetPalette()[cstep * (iint + 1)],
775  left_painter);
776  // set up links between painters
777  if (right_painter) {
778  right_painter->set_left_interval(dummy);
779  dummy->set_right_interval(right_painter);
780  }
781  fPad->cd();
782  dummy->Draw();
783  dummy->Connect("IntMod()", "KVItvFinderDialog", this, "UpdatePIDList()");
784  dummy->SetDisplayLabel(1);
785  dummy->SetCanvas(fCanvas);
786  fItvPaint.Add(dummy);
787 
789 
790  fItvPaint.Execute("Update", "");
791 
792  fCanvas->Modified();
793  fCanvas->Update();
794 }
795 
796 
797 
799 
801 {
802  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
803  if (!list->GetSize() || list->GetSize() > 1) {
804  current_interval_set = nullptr;
805  return;
806  }
807 
808  current_interval_set = (interval_set*)list->At(0);
809 
810  int aa = 0;
811  int iint = 0;
812 
813  // try to guess mass from position (PID)
814  // typical isotope separation in PID is 0.12 (e.g. for carbon isotopes)
815  // assume that PID=z means A=2*Z
816  auto aa_guessstimate = TMath::Nint(current_interval_set->GetZ() * 2 + (pid - current_interval_set->GetZ()) / 0.12);
817 //Info("NewINt","pid = %f guess a=%d",pid,aa_guessstimate);
818  if (!current_interval_set->GetNPID()) {
819  aa = aa_guessstimate;
820  iint = 0;
821  }
822  else if (pid < ((interval*)current_interval_set->GetIntervals()->First())->GetPID()) { // to left of all others
823  aa = ((interval*)current_interval_set->GetIntervals()->First())->GetA() - 1;
824  // NO MORE : use guesstimate as long as it is smaller than A of previously defined interval with largest PID
825  // NOW : use 'aa_guessstimate' only for first interval
826  //Info("NewInt","before first interval, with a=%d",aa+1);
827 // if (aa_guessstimate < aa + 1) aa = aa_guessstimate;
828  iint = 0;
829  }
830  else if (pid > ((interval*)current_interval_set->GetIntervals()->Last())->GetPID()) { // to right of all others
831  aa = ((interval*)current_interval_set->GetIntervals()->Last())->GetA() + 1;
832  // NO MORE : use guesstimate as long as it is larger than A of previously defined interval with smallest PID
833  // NOW : use 'aa_guessstimate' only for first interval
834  //Info("NewInt","after last interval, with a=%d",aa-1);
835 // if (aa_guessstimate > aa - 1) aa = aa_guessstimate;
836  iint = current_interval_set->GetNPID();
837  }
838  else {
839  // look for intervals between which the new one is places
840  for (int ii = 1; ii < current_interval_set->GetNPID(); ii++) {
841  if (pid > ((interval*)current_interval_set->GetIntervals()->At(ii - 1))->GetPID()
842  && pid < ((interval*)current_interval_set->GetIntervals()->At(ii))->GetPID()) {
843  aa = ((interval*)current_interval_set->GetIntervals()->At(ii - 1))->GetA() + 1;
844  // use guesstimate if it is in between masses of adjacent intervals (in terms of PID)
845  if ((aa_guessstimate > (aa - 1)) &&
846  (aa_guessstimate < ((interval*)current_interval_set->GetIntervals()->At(ii))->GetA()))
847  aa = aa_guessstimate;
848  iint = ii;
849  break;
850  }
851  }
852  }
853 
854  interval* itv = new interval(current_interval_set->GetZ(), aa, pid, pid - 0.05, pid + 0.05);
855  current_interval_set->GetIntervals()->AddAt(itv, iint);
856 
857  // find intervals which are now left (smaller mass) and right (higher mass) than this one
858  interval* left_interval = iint > 0 ? (interval*)current_interval_set->GetIntervals()->At(iint - 1) : nullptr;
859  interval* right_interval = iint < current_interval_set->GetIntervals()->GetEntries() - 1 ? (interval*)current_interval_set->GetIntervals()->At(iint + 1) : nullptr;
860  // find the corresponding painters
861  KVPIDIntervalPainter* left_painter{nullptr}, *right_painter{nullptr};
862  TIter next_painter(&fItvPaint);
863  KVPIDIntervalPainter* pidpnt;
864  while ((pidpnt = (KVPIDIntervalPainter*)next_painter())) {
865  if (pidpnt->GetInterval() == left_interval) left_painter = pidpnt;
866  else if (pidpnt->GetInterval() == right_interval) right_painter = pidpnt;
867  }
868 
869  // give a different colour to each interval
870  auto nitv = current_interval_set->GetIntervals()->GetEntries();
871  auto cstep = TColor::GetPalette().GetSize() / (nitv + 1);
872 
873  KVPIDIntervalPainter* dummy = new KVPIDIntervalPainter(itv, fLinearHisto, TColor::GetPalette()[cstep * (iint + 1)],
874  left_painter);
875  // set up links between painters
876  if (right_painter) {
877  right_painter->set_left_interval(dummy);
878  dummy->set_right_interval(right_painter);
879  }
880  fPad->cd();
881  dummy->Draw();
882  dummy->Connect("IntMod()", "KVItvFinderDialog", this, "UpdatePIDList()");
883  dummy->SetDisplayLabel(1);
884  dummy->SetCanvas(fCanvas);
885  fItvPaint.Add(dummy);
886 
888 
889  fItvPaint.Execute("Update", "");
890 
891  fCanvas->Modified();
892  fCanvas->Update();
893 
894 }
895 
896 
897 
899 
901 {
902  if (fGrid->GetIntervalSets()->GetSize() == 0) fNextIntervalZ = 1;
903  else fNextIntervalZ = ((interval_set*)fGrid->GetIntervalSets()->Last())->GetZ() + 1;
904  // KVBase::OpenContextMenu("SetNextIntervalZ()",this);
906  UpdateLists();
907 }
908 
909 
910 
912 
913 void KVItvFinderDialog::remove_interval_from_interval_set(interval_set* itvs, interval* itv, bool remove_fit)
914 {
915  KVPIDIntervalPainter* pid = (KVPIDIntervalPainter*)fItvPaint.FindObject(Form("%d_%d", itv->GetZ(), itv->GetA()));
916  itvs->GetIntervals()->Remove(itv);
918  if (remove_fit) {
919  // remove any fits from pad corresponding to intervals
921  }
922 }
923 
924 
925 
927 
929 {
930  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
931  Int_t nSelected = list->GetSize();
932  current_interval_set = nullptr;
933 
934  if (nSelected == 1) {
935  current_interval_set = (interval_set*)list->At(0);
936  list.reset(fIntervalListView->GetSelectedObjects());
937  nSelected = list->GetSize();
938  if (nSelected >= 1) {
939  for (int ii = 0; ii < nSelected; ii++) {
940  interval* itv = (interval*) list->At(ii);
941  remove_interval_from_interval_set(current_interval_set, itv);
942  }
944  fCanvas->Modified();
945  fCanvas->Update();
946  }
948  }
949  else if (nSelected > 1) {
950  for (int ii = 0; ii < nSelected; ii++) {
951  auto itvs = (interval_set*)list->At(ii);
952  ClearInterval(itvs);
953  }
954  }
955 }
956 
957 
958 
960 
962 {
963  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
964  Int_t nSelected = list->GetSize();
965 
966  current_interval_set = nullptr;
967  if (nSelected == 1) {
968  current_interval_set = (interval_set*)list->At(0);
969 
970  list.reset(fIntervalListView->GetSelectedObjects());
971  nSelected = list->GetSize();
972 
973  if (nSelected == 1) {
974  interval* itv = (interval*) list->At(0);
975  itv->SetA(itv->GetA() + 1);
976  fItvPaint.Execute("Update", "");
977  // change the name of any gaussian in the pad associated with this isotope
979  if (gfit) gfit->SetName(KVMultiGaussIsotopeFit::get_name_of_isotope_gaussian(itv->GetZ(), itv->GetA()));
980  fCanvas->Modified();
981  fCanvas->Update();
982  }
983  else {
984  TIter next_itv(list.get());
985  interval* itv;
986  while ((itv = (interval*)next_itv())) {
987  itv->SetA(itv->GetA() + 1);
988  // change the name of any gaussian in the pad associated with this isotope
990  if (gfit) gfit->SetName(KVMultiGaussIsotopeFit::get_name_of_isotope_gaussian(itv->GetZ(), itv->GetA()));
991  }
992  fItvPaint.Execute("Update", "");
993  fCanvas->Modified();
994  fCanvas->Update();
995  }
996  }
997 }
998 
999 
1000 
1002 
1004 {
1005  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
1006  Int_t nSelected = list->GetSize();
1007 
1008  current_interval_set = nullptr;
1009  if (nSelected == 1) {
1010  current_interval_set = (interval_set*)list->At(0);
1011  list.reset(fIntervalListView->GetSelectedObjects());
1012  nSelected = list->GetSize();
1013 
1014  if (nSelected == 1) {
1015  interval* itv = (interval*) list->At(0);
1016  itv->SetA(itv->GetA() - 1);
1017  fItvPaint.Execute("Update", "");
1018  // change the name of any gaussian in the pad associated with this isotope
1020  if (gfit) gfit->SetName(KVMultiGaussIsotopeFit::get_name_of_isotope_gaussian(itv->GetZ(), itv->GetA()));
1021  fCanvas->Modified();
1022  fCanvas->Update();
1023  }
1024  else {
1025  TIter next_itv(list.get());
1026  interval* itv;
1027  while ((itv = (interval*)next_itv())) {
1028  itv->SetA(itv->GetA() - 1);
1029  // change the name of any gaussian in the pad associated with this isotope
1031  if (gfit) gfit->SetName(KVMultiGaussIsotopeFit::get_name_of_isotope_gaussian(itv->GetZ(), itv->GetA()));
1032  }
1033  fItvPaint.Execute("Update", "");
1034  fCanvas->Modified();
1035  fCanvas->Update();
1036  }
1037  }
1038 }
1039 
1040 
1041 
1043 
1045 {
1046  fIntervalSetListView->Display(((KVIDZAFromZGrid*)fGrid)->GetIntervalSets());
1047  std::unique_ptr<TList> list(fIntervalSetListView->GetSelectedObjects());
1048  Int_t nSelected = list->GetSize();
1049  interval_set* itvs = 0;
1050  if (nSelected == 1) {
1051  current_interval_set = (interval_set*)list->At(0);
1053  }
1054  else {
1055  current_interval_set = nullptr;
1057  }
1058 }
1059 
1060 
1061 
1064 
1066 {
1067  //fGrid->SetOnlyZId(0);
1068  fGrid->Initialize();
1069  ExportToGrid();
1070 
1072 
1073  fIntervalSetListView->Display(((KVIDZAFromZGrid*)fGrid)->GetIntervalSets());
1074  current_interval_set = nullptr;
1076 
1077  fItvPaint.Clear();
1078  DrawIntervals();
1079 
1080  new KVTestIDGridDialog(gClient->GetDefaultRoot(), gClient->GetDefaultRoot(), 10, 10, fGrid, fHisto);
1081 }
1082 
1083 
1084 
1086 
1088 {
1090  fCanvas->Modified();
1091  fCanvas->Update();
1092 }
1093 
1094 
1095 
1097 
1099 {
1100  fItvPaint.Execute("SetDisplayLabel", "0");
1101  current_interval_set = nullptr;
1103  fCanvas->Modified();
1104  fCanvas->Update();
1105 }
1106 
1107 
1108 
1119 
1121 {
1122  // fit the PID spectrum for the currently selected interval set (Z).
1123  //
1124  // for an interval with N isotopes, we use N gaussians plus an exponential (decreasing)
1125  // background. each gaussian has the same width. the centroids of the gaussians are first
1126  // fixed to the positions of the PID markers, the intensity and width of the peaks
1127  // (plus the background) are fitted.
1128  // then another fit is performed without constraining the centroids.
1129  //
1130  // finally the PID markers (PID of each interval) are modified according to the fitted centroid positions.
1131 
1132  if (!current_interval_set) return;
1133 
1134  KVNumberList alist;
1135  std::vector<double> pidlist;
1137  interval* intvl = nullptr;
1138  while ((intvl = (interval*)nxt_int())) {
1139  alist.Add(intvl->GetA());
1140  pidlist.push_back(intvl->GetPID());
1141  }
1144  alist, pidlist);
1145 
1146  // check user fit parameters
1147  if (mass_fit_parameters.GetBoolValue("Limit range of fit")) {
1148  // if the user's range is not valid for the current interval set, we ignore it
1149  if (mass_fit_parameters.GetDoubleValue("PID min for fit") >= current_interval_set->GetZ() - 0.5
1150  && mass_fit_parameters.GetDoubleValue("PID max for fit") <= current_interval_set->GetZ() + 0.5)
1151  fitfunc.SetFitRange(mass_fit_parameters.GetDoubleValue("PID min for fit"),
1152  mass_fit_parameters.GetDoubleValue("PID max for fit"));
1153 
1154  }
1155  fitfunc.SetSigmaLimits(mass_fit_parameters.GetDoubleValue("Minimum #sigma"),
1156  mass_fit_parameters.GetDoubleValue("Maximum #sigma"));
1157 
1158  fLinearHisto->Fit(&fitfunc, "NR");
1159 
1160  // now release the centroids
1161  fitfunc.ReleaseCentroids();
1162 
1163  fLinearHisto->Fit(&fitfunc, "NRME");
1164 
1165  // remove any previous fit from pad
1166  fitfunc.UnDraw();
1167  // mop up any stray gaussians (from intervals which have been removed)
1169 
1170  // draw fit with individual gaussians
1171  fitfunc.DrawFitWithGaussians("same");
1172 
1173  // deactivate intervals for fitted masses
1174  std::unique_ptr<KVSeqCollection> tmp(fItvPaint.GetSubListWithMethod(Form("%d", current_interval_set->GetZ()), "GetZ"));
1175  tmp->Execute("DeactivateIntervals", "");
1176 
1177  // set interval limits according to regions of most probable mass
1178  int most_prob_A = 0;
1179  nxt_int.Reset();
1180  // minimum probability for which isotopes are taken into account
1181  double min_proba = mass_fit_parameters.GetDoubleValue("Minimum probability [%]") / 100.;
1182  double delta_pid = 0.001;
1183  TList accepted_intervals;//any intervals not in this list at the end of the procedure will be removed
1184  for (double pid = fitfunc.GetPIDmin() ; pid <= fitfunc.GetPIDmax(); pid += delta_pid) {
1185  double proba;
1186  auto Amax = fitfunc.GetMostProbableA(pid, proba);
1187  if (proba > min_proba) {
1188  if (most_prob_A) {
1189  if (Amax > most_prob_A) {
1190  //std::cout << pid << " " << Amax << " " << proba << std::endl;
1191  //std::cout << "got PIDmax for A=" << intvl->GetA() << std::endl;
1192  intvl->SetPIDmax(pid - delta_pid);
1193  most_prob_A = Amax;
1194  intvl = (interval*)nxt_int();
1195  while (intvl->GetA() < most_prob_A) {
1196  intvl = (interval*)nxt_int();
1197  }
1198  accepted_intervals.Add(intvl);
1199  intvl->SetPIDmin(pid);
1200  //std::cout << "got PIDmin for A=" << intvl->GetA() << std::endl;
1201  }
1202  }
1203  else {
1204  most_prob_A = Amax;
1205  intvl = (interval*)nxt_int();
1206  while (intvl->GetA() < most_prob_A) {
1207  intvl = (interval*)nxt_int();
1208  }
1209  accepted_intervals.Add(intvl);
1210  intvl->SetPIDmin(pid);
1211  //std::cout << pid << " " << Amax << " " << proba << std::endl;
1212  //std::cout << "got PIDmin for A=" << intvl->GetA() << std::endl;
1213  }
1214  }
1215  else if (most_prob_A) {
1216  intvl->SetPIDmax(pid - delta_pid);
1217  //std::cout << pid << " " << Amax << " " << proba << std::endl;
1218  //std::cout << "got PIDmax for A=" << intvl->GetA() << std::endl;
1219  most_prob_A = 0;
1220  }
1221  }
1222  nxt_int.Reset();
1223 // TList intervals_to_remove;
1224 // while ((intvl = (interval*)nxt_int())) {
1225 // if (!accepted_intervals.FindObject(intvl)) intervals_to_remove.Add(intvl);
1226 // }
1227 // if (intervals_to_remove.GetEntries()) {
1228 // // remove intervals below minimum probability (leave gaussians on display)
1229 // TIter it_rem(&intervals_to_remove);
1230 // while ((intvl = (interval*)it_rem())) remove_interval_from_interval_set(current_interval_set, intvl, false);
1231 // }
1232 // intervals_to_remove.Clear();
1233  int ig(1);
1234  // update PID positions from fitted centroids
1235  nxt_int.Reset();
1236  KVNumberList remaining_gaussians, remaining_alist;
1237  auto vec_alist = alist.GetArray();
1238  while ((intvl = (interval*)nxt_int())) {
1239  // in case we removed some peaks
1240  while (vec_alist[ig - 1] < intvl->GetA()) {
1241  ++ig;
1242  }
1243  intvl->SetPID(fitfunc.GetCentroid(ig));
1244  remaining_gaussians.Add(ig);
1245  remaining_alist.Add(intvl->GetA());
1246  ++ig;
1247  }
1248 // if (intervals_to_remove.GetEntries()) {
1249 // // remove intervals with centroids outside PID limits
1250 // TIter it_rem(&intervals_to_remove);
1251 // while ((intvl = (interval*)it_rem())) remove_interval_from_interval_set(current_interval_set, intvl, false);
1252 // }
1253  UpdatePIDList();
1254  fItvPaint.Execute("Update", "");
1255 
1256  fPad->Modified();
1257  fPad->Update();
1258 
1259  // save results in grid parameters
1260  KVNumberList zlist;
1261  if (fGrid->GetParameters()->HasStringParameter("MASSFITS"))
1262  zlist.Set(fGrid->GetParameters()->GetStringValue("MASSFITS"));
1263  zlist.Add(current_interval_set->GetZ());
1264  fGrid->GetParameters()->SetValue("MASSFITS", zlist.AsString());
1265  TString massfit = Form("MASSFIT_%d", current_interval_set->GetZ());
1266  KVNameValueList fitparams;
1267  fitparams.SetValue("Ng", alist.GetNValues());
1268  fitparams.SetValue("Alist", alist.AsQuotedString());
1269  fitparams.SetValue("PIDmin", fitfunc.GetPIDmin());
1270  fitparams.SetValue("PIDmax", fitfunc.GetPIDmax());
1271  fitparams.SetValue("Bkg_cst", fitfunc.GetBackgroundConstant());
1272  fitparams.SetValue("Bkg_slp", fitfunc.GetBackgroundSlope());
1273  fitparams.SetValue("GausWid", fitfunc.GetGaussianWidth(0));
1274  fitparams.SetValue("PIDvsA_a0", fitfunc.GetPIDvsAfit_a0());
1275  fitparams.SetValue("PIDvsA_a1", fitfunc.GetPIDvsAfit_a1());
1276  fitparams.SetValue("PIDvsA_a2", fitfunc.GetPIDvsAfit_a2());
1277  for (ig = 1; ig <= alist.GetNValues(); ++ig) {
1278  fitparams.SetValue(Form("Norm_%d", ig), fitfunc.GetGaussianNorm(ig));
1279  }
1280  auto sanitized = fitparams.Get().ReplaceAll("=", ":");
1281  fGrid->GetParameters()->SetValue(massfit, sanitized);
1282 }
1283 
1284 
1285 
1288 
1290 {
1291  // Open dialog to modify parameters for multigauss mass fit
1292 
1293  bool cancel = false;
1294  auto dialog = new KVNameValueListGUI(fMain, &mass_fit_parameters, &cancel);
1295  dialog->EnableDependingOnBool("PID min for fit", "Limit range of fit");
1296  dialog->EnableDependingOnBool("PID max for fit", "Limit range of fit");
1297  dialog->DisplayDialog();
1298 }
1299 
1300 
1301 
1304 
1306 {
1307  // Remove fit of currently selected interval set from pad
1308 
1309  if (!current_interval_set) return;
1310 
1311  std::vector<int> alist;
1313  interval* intvl = nullptr;
1314  while ((intvl = (interval*)nxt_int())) {
1315  alist.push_back(intvl->GetA());
1316  }
1318 
1319  fitfunc.UnDraw(fPad);
1320  // mop up any stray gaussians (from intervals which have been removed)
1322 
1323  fPad->Modified();
1324  fPad->Update();
1325 }
1326 
1327 
1328 
1330 
1332 {
1333 
1334  if (!fCanvas) return;
1335  if (fCanvas->GetEvent() == kMouseMotion) return;
1336 
1337  switch (fCanvas->GetEvent()) {
1338  case kButton2Up:
1339  FitIsotopes();
1340  break;
1341  case kButton1Double:
1343  break;
1344  case kButton1Shift:
1345 // std::cout << fCanvas->GetEventX() << " " << fCanvas->AbsPixeltoX(fCanvas->GetEventX()) << std::endl;
1347  break;
1348  }
1349 
1350 }
1351 
1352 
1353 
1355 
1357 {
1358  interval_set* itvs = fGrid->GetIntervalSet(zz);
1359  if (!itvs) {
1360  itvs = new interval_set(zz, KVIDZAFromZGrid::kIntType);
1361  fGrid->GetIntervalSets()->Add(itvs);
1362  }
1363  else ClearInterval(itvs);
1364 
1365 
1366  if (zz == 1) fLinearHisto->SetAxisRange(0.9, zz + 0.5, "X");
1367  else fLinearHisto->SetAxisRange(zz - 0.5, zz + 0.5, "X");
1368 
1369  int nfound = fSpectrum.Search(fLinearHisto, fSig, "goff", ((zz == 2) ? 0.1 * fRat : fRat));
1370 
1371 #if ROOT_VERSION_CODE > ROOT_VERSION(5,99,01)
1372  Double_t* xpeaks = fSpectrum.GetPositionX();
1373 // Double_t* ypeaks = fSpectrum.GetPositionY();
1374 #else
1375  Float_t* xpeaks = fSpectrum.GetPositionX();
1376 // Float_t* ypeaks = fSpectrum.GetPositionY();
1377 #endif
1378 
1379  nfound = TMath::Min(fNpeaks[zz - 1], nfound);
1380 
1381  int idx[nfound];
1382  TMath::Sort(nfound, xpeaks, idx, 0);
1383 
1384  int zrefs[] = {1, 4, 7, 9, 11, 12, 15, 16, 19, 21, 23, 25, 27, 29, 31, 34, 35, 38, 40, 42, 44, 47, 49, 51, 53};
1385  int zref = zrefs[zz - 1];
1386 
1387  int idref = -1;
1388  for (int p = 0; p < nfound; p++) {
1389  if (abs(xpeaks[idx[p]] - xpeaks[0]) < .0001) idref = p;
1390  }
1391  Info("FindPIDIntervals", "Z=%d : idref = %d ", zz, idref);
1392 
1393  for (int p = 0; p < nfound; p++) {
1394 // Info("FindPIDIntervals","Z=%d : (%.2lf %.2lf) (%d %d) ",zz,xpeaks[p],ypeaks[p],idx[p],p);
1395  double pid = xpeaks[idx[p]];//ff->GetParameter(3 * ii + 1);
1396  itvs->add(zref + (p - idref), pid, pid - 0.05, pid + 0.05);
1397  }
1398 
1399 
1400 }
1401 
1402 
1403 
1405 
1407 {
1408  Double_t result = 0;
1409  int np = par[0];
1410  for (Int_t p = 0; p < np; p++) {
1411  Double_t norm = par[3 * p + 3];
1412  Double_t mean = par[3 * p + 1];
1413  Double_t sigma = par[3 * p + 2];
1414  result += norm * TMath::Gaus(x[0], mean, sigma);
1415  }
1416  return result;
1417 }
1418 
1419 
1421 
1423 {
1424 
1425  if (zmin < 0) zmin = ((KVIDentifier*)fGrid->GetIdentifiers()->First())->GetZ();
1426  if (zmax < 0) zmax = ((KVIDentifier*)fGrid->GetIdentifiers()->Last())->GetZ();
1427 
1428  for (int z = zmin; z <= zmax; z++) FindPIDIntervals(z);
1429 
1430 }
1431 
1432 
1433 
1434 
1435 
1436 
1437 
1438 
1439 
1440 
kMouseMotion
kButton1Double
kButton1Shift
kButton2Up
int Int_t
kVerticalFrame
kHorizontalFrame
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
const Bool_t kFALSE
double Axis_t
double Double_t
double Stat_t
float Float_t
kGray
kBlack
#define gClient
kFDSave
kDeepCleanup
kLHintsExpandY
kLHintsTop
kLHintsExpandX
kMBNo
kMBYes
kMBIconExclamation
kTextCenterX
kTextLeft
#define gROOT
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
char * StrDup(const char *str)
R__EXTERN TStyle * gStyle
R__EXTERN TSystem * gSystem
static void SetAutoAdd(Bool_t yes=kTRUE)
Definition: KVIDGraph.h:99
static Bool_t GetAutoAdd()
Definition: KVIDGraph.h:105
const KVNameValueList * GetParameters() const
Definition: KVIDGraph.h:287
const Char_t * GetName() const
Definition: KVIDGraph.cpp:1332
void WriteAsciiFile(const Char_t *filename)
Open, write and close ascii file containing this grid.
Definition: KVIDGraph.cpp:404
const KVList * GetIdentifiers() const
Definition: KVIDGraph.h:297
Hybrid charge & mass identification grid.
interval_set * GetIntervalSet(int zint) const
void SetOnlyZId(Bool_t=kTRUE)
KVList * GetIntervalSets()
Base class for graphical cuts used in particle identification.
Definition: KVIDentifier.h:27
virtual Double_t GetPID() const
Full result of one attempted particle identification.
void Clear(Option_t *opt="")
Reset to initial values.
Double_t PID
= "real" Z if Zident==kTRUE and Aident==kFALSE, "real" A if Zident==Aident==kTRUE
Bool_t HasFlag(std::string grid_name, TString flag)
GUI for finding/fixing mass identification intervals.
KVListView * fIntervalListView
void DrawInterval(interval_set *itvs, bool label=0)
KVIDZAFromZGrid * fGrid
TGTransientFrame * fMain
void ClearInterval(interval_set *itvs)
void ProcessIdentification(Int_t zmin=-1, Int_t zmax=-1)
void AddInterval(double pid)
void delete_painter_from_painter_list(KVPIDIntervalPainter *)
TVirtualPad * fPad
void Identify()
KVBase::OpenContextMenu("Identify(double,double)",this);.
void SetFitParameters()
Open dialog to modify parameters for multigauss mass fit.
interval_set * current_interval_set
void FindPIDIntervals(Int_t zz)
KVItvFinderDialog(KVIDZAFromZGrid *gg, TH2 *hh)
KVListView * fIntervalSetListView
static KVNameValueList mass_fit_parameters
for user control of multi-gaussian fit
void RemoveFit()
Remove fit of currently selected interval set from pad.
void ExportToGrid()
Write all PID intervals in grid parameters "PIDRANGE", "PIDRANGE%d", etc.
void TestIdent()
fGrid->SetOnlyZId(0);
void LinearizeHisto(int nbins)
Double_t fpeaks(Double_t *x, Double_t *par)
KVPIDIntervalPainter * last_drawn_interval
void ZoomOnCanvas()
Display the interval set for a given Z when the user double clicks on it.
Enhanced version of ROOT TGListView widget.
Definition: KVListView.h:145
virtual void SetDataColumns(Int_t ncolumns)
Definition: KVListView.cpp:91
void SetDoubleClickAction(const char *receiver_class, void *receiver, const char *slot)
Definition: KVListView.cpp:210
virtual void Display(const TCollection *l)
Definition: KVListView.h:172
TObject * GetLastSelectedObject() const
Definition: KVListView.h:230
TList * GetSelectedObjects() const
Definition: KVListView.h:244
virtual void RemoveAll()
Definition: KVListView.h:189
void AllowContextMenu(Bool_t on=kTRUE)
Definition: KVListView.h:282
virtual void SetDataColumn(Int_t index, const Char_t *name, const Char_t *method="", Int_t mode=kTextCenterX)
Definition: KVListView.cpp:106
Extended TList class which owns its objects by default.
Definition: KVList.h:27
Function for fitting PID mass spectra.
static TString get_name_of_isotope_gaussian(int z, int a)
static TString get_name_of_multifit(int z)
double GetGaussianWidth(int) const
double GetCentroid(int i) const
void SetSigmaLimits(double smin, double smax)
void SetFitRange(double min, double max)
Change range of fit.
int GetMostProbableA(double PID, double &P) const
double GetGaussianNorm(int i) const
static void UnDrawAnyGaussian(int z, TVirtualPad *pad=gPad)
void UnDraw(TVirtualPad *pad=gPad) const
Remove the graphical representation of this fit from the given pad.
static void UnDrawGaussian(int z, int a, TVirtualPad *pad=gPad)
void DrawFitWithGaussians(Option_t *opt="") const
Draw the overall fit plus the individual gaussians for each isotope.
GUI for setting KVNameValueList parameters.
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Double_t GetDoubleValue(const Char_t *name) const
void SetValue(const Char_t *name, value_type value)
void RemoveParameter(const Char_t *name)
Bool_t HasStringParameter(const Char_t *name) const
Bool_t GetBoolValue(const Char_t *name) const
KVString Get() const
const Char_t * GetStringValue(const Char_t *name) const
bool Set(const KVString &)
Bool_t HasParameter(const Char_t *name) const
TString GetTStringValue(const Char_t *name) const
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
const Char_t * AsQuotedString() const
const Char_t * AsString(Int_t maxchars=0) const
void Remove(Int_t)
Remove value 'n' from the list.
Int_t GetNValues() const
void Add(Int_t)
Add value 'n' to the list.
void Set(const TString &l)
Definition: KVNumberList.h:133
IntArray GetArray() const
Graphical representation of a PID interval in the KVIDZAFromZGrid mass assignation GUI.
void SetCanvas(TCanvas *cc)
KVPIDIntervalPainter * get_right_interval() const
void Draw(Option_t *option="")
interval * GetInterval() const
void SetDisplayLabel(bool dis=true)
void HighLight(bool hi=true)
KVPIDIntervalPainter * get_left_interval() const
void set_right_interval(KVPIDIntervalPainter *i)
KVSeqCollection * GetSubListWithMethod(const Char_t *retvalue, const Char_t *method) const
virtual TObject * Last() const
virtual void Clear(Option_t *option="")
virtual Int_t GetSize() const
virtual TObject * At(Int_t idx) const
virtual TObject * First() const
virtual void Execute(const char *method, const char *params, Int_t *error=0)
virtual void Add(TObject *obj)
virtual void AddAt(TObject *obj, Int_t idx)
virtual TObject * Remove(TObject *obj)
Remove object from list.
virtual TObject * FindObject(const char *name) const
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
GUI for testing identification grids.
Int_t GetSize() const
virtual void SetFillColor(Color_t fcolor)
virtual void SetLineColor(Color_t lcolor)
virtual void SetMarkerStyle(Style_t mstyle=1)
virtual void SetBottomMargin(Float_t bottommargin)
virtual void SetLeftMargin(Float_t leftmargin)
virtual void SetRightMargin(Float_t rightmargin)
virtual void SetTopMargin(Float_t topmargin)
virtual void UnZoom()
virtual void SetRangeUser(Double_t ufirst, Double_t ulast)
Int_t GetEventX() const override
TVirtualPad * cd(Int_t subpadnumber=0) override
void Update() override
Int_t GetEvent() const override
virtual Int_t GetEntries() const
static const TArrayI & GetPalette()
virtual void AddFrame(TGFrame *f, TGLayoutHints *l=0)
virtual void SetCleanup(Int_t mode=kLocalCleanup)
virtual TGDimension GetDefaultSize() const
virtual void MapSubwindows()
char * fFilename
const char ** fFileTypes
char * fIniDir
virtual void Resize(TGDimension size)
virtual void MapWindow()
void DontCallClose()
void SetWindowName(const char *name=0)
virtual TGButton * AddButton(const TGWindow *w, TGPictureButton *button, Int_t spacing=0)
virtual void CenterOnParent(Bool_t croot=kTRUE, EPlacement pos=kCenter)
virtual void RequestFocus()
virtual Int_t GetNbinsY() const
TAxis * GetXaxis()
TAxis * GetYaxis()
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
virtual Int_t GetNbinsX() const
virtual Int_t Fill(const char *name, Double_t w)
virtual void Draw(Option_t *option="")
virtual void SetAxisRange(Double_t xmin, Double_t xmax, Option_t *axis="X")
virtual Double_t GetBinContent(Int_t bin) const
void Reset()
virtual void Add(TObject *obj)
virtual TObject * Last() const
virtual const char * GetName() const
void AddExec(const char *name, const char *command) override
Double_t AbsPixeltoX(Int_t px) override
void Modified(Bool_t flag=1) override
void SetLogy(Int_t value=1) override
Int_t GetLogy() const override
Bool_t Connect(const char *signal, const char *receiver_class, void *receiver, const char *slot)
virtual Double_t Uniform(Double_t x1, Double_t x2)
void AdoptCanvas(TCanvas *c)
TCanvas * GetCanvas() const
Int_t GetCanvasWindowId() const
virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
Double_t * GetPositionX() const
Ssiz_t Length() const
const char * Data() 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)
void SetOptTitle(Int_t tit=1)
void SetOptStat(Int_t stat=1)
virtual char * ExpandPathName(const char *path)
virtual TList * GetListOfPrimitives() const=0
virtual void Modified(Bool_t flag=1)=0
virtual TObject * GetPrimitive(const char *name) const=0
virtual void Update()=0
virtual TObject * WaitPrimitive(const char *pname="", const char *emode="")=0
virtual TVirtualPad * cd(Int_t subpadnumber=0)=0
void add(int aa, double pid, double pidmin=-1., double pidmax=-1.)
KVList * GetIntervals()
bool SetPID(double pid)
double GetPID()
double GetPIDmin()
void SetPIDmin(double pidmin)
void SetA(int aa)
double GetPIDmax()
void SetPIDmax(double pidmax)
const Double_t sigma
Double_t y[n]
Double_t x[n]
const long double mm
Definition: KVUnits.h:69
void Info(const char *location, const char *va_(fmt),...)
Double_t Min(Double_t a, Double_t b)
Double_t Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
Int_t Nint(T x)
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
const char * fTipText
TGButton * fButton
Bool_t fStayDown
const char * fPixmap