KaliVeda
Toolkit for HIC analysis
KVDataQualityAuditReporting_INDRAFAZIA.cpp
1 #include "KVDataQualityAuditReporting_INDRAFAZIA.h"
2 #include "KVDataSetManager.h"
3 #include "KVINDRA.h"
4 #include "KVFAZIA.h"
5 #include <KVFAZIADetector.h>
6 #include <KVGeoDNTrajectory.h>
7 #include <TColor.h>
8 #include <TLatex.h>
9 #include <TLegend.h>
10 #include <TMarker.h>
11 #include "TLegend.h"
12 
14 
15 
16 
20 {
21  // Make an A4-size canvas
22  Double_t w = 297 * 6;
23  Double_t h = 210 * 6;
24  if (style == canvas_t::kPortrait) std::swap(w, h);
25  myCanvas = ::new TCanvas("c", "c", w, h);
26  myCanvas->SetWindowSize(w + (w - myCanvas->GetWw()), h + (h - myCanvas->GetWh()));
27 }
28 
29 
30 
33 
35 {
36  // sort fazia telescopes into bins of theta
37 
38  TIter itdet(gFazia->GetDetectors());
39  KVDetector* det;
40 
41  while ((det = (KVDetector*)itdet())) {
42  if (det->IsLabelled("SI1")) {
43  int bin = TMath::Nint(det->GetTheta() * (1 / theta_bin));
44  fazia_map[bin].push_back(det);
45  }
46  }
47  // Now sort detectors in vectors in order of increasing theta
48  for (auto& p : fazia_map) {
49  std::sort(std::begin(p.second), std::end(p.second), [](KVDetector * a, KVDetector * b) {
50  return a->GetTheta() < b->GetTheta();
51  });
52  }
53 }
54 
55 
56 
58 
60 {
61  if (!gMultiDetArray)
63 
64  indra_csi_idtype = gIndra->GetCsIIDType();
66 
67  make_fazia_map(0.4);
68 
69  make_canvas();
70 
71  TString pdf_file = Form("%s_DataQualityAudit_%s.pdf", dataset_name.Data(), fAudit->GetName());
72  auto first_page = pdf_file + "(";
73  auto last_page = pdf_file + ")";
74  current_page = first_page;
75  // INDRA Z identification bilan
76  int nx(2), ny(3);
77  myCanvas->Divide(nx, ny);
78  int pad = 1;
79  std::vector<TString> indra_id_types = {"SI", "SI_CSI", indra_csi_idtype};
80  for (auto& id : indra_id_types) {
81  for (int ring = 6; ring <= 17; ++ring) {
82  if ((id == "SI_CSI" || id == "SI") && ring > 9) break;
83 
84  myCanvas->cd(pad)->SetFillStyle(4000);//transparent pad
85  if(INDRA_ring_reporting_Z(ring, id))
86  {
87  ++pad;
88  if (pad > nx * ny) {
90  myCanvas->Print(current_page, Form("Title:INDRA Z %s ID Quality by Ring", id.Data()));
91  if (current_page == first_page) current_page = pdf_file;
92  myCanvas->Clear();
93  myCanvas->Divide(nx, ny);
94  pad = 1;
95  }
96  }
97  }
98  if (pad != 1) {
100  myCanvas->Print(current_page, Form("Title:INDRA Z %s ID Quality by Ring", id.Data()));
101  if (current_page == first_page) current_page = pdf_file;
102  myCanvas->Clear();
103  myCanvas->Divide(nx, ny);
104  pad = 1;
105  }
106  nx = 4;
107  ny = 3;
108  myCanvas->Clear();
109  myCanvas->Divide(nx, ny);
110  pad = 1;
111  for (int ring = 6; ring <= 17; ++ring) {
112  if ((id == "SI_CSI" || id == "SI")) {
113  if (ring > 9) break; // no si after ring 9
114  }
115 
116  INDRA_ring_mean_A_vs_Z(ring, id, pad, nx, ny);
117  }
118  if (pad != 1) {
119  myCanvas->Print(current_page, Form("Title:INDRA <A> vs Z %s by Ring", id.Data()));
120  myCanvas->Clear();
121  myCanvas->Divide(nx, ny);
122  pad = 1;
123  }
124  for (int ring = 6; ring <= 17; ++ring) {
125  INDRA_ring_Z_threshold_vs_Z(ring, id, pad, nx, ny);
126  }
127  if (pad != 1) {
128  myCanvas->Print(current_page, Form("Title:INDRA Z thresh. vs Z %s by Ring", id.Data()));
129  myCanvas->Clear();
130  myCanvas->Divide(nx, ny);
131  pad = 1;
132  }
133  }
134 
135  // FAZIA Z identification bilan
136  std::vector<TString> fazia_id_types = {"SI1", "SI1_SI2", fazia_si_csi_idtype, "CSI"};
137  nx = 4;
138  ny = 3;
139  for (auto id : fazia_id_types) {
140  myCanvas->Clear();
141  myCanvas->Divide(nx, ny);
142  pad = 1;
143  int group_num = 1;
144  for (auto& p : fazia_map) {
145  myCanvas->cd(pad);
146 
147  if(FAZIA_group_reporting_Z(group_num, p.second, id))
148  {
149  ++group_num;
150  ++pad;
151  if (pad > nx * ny) {
153  myCanvas->Print(pdf_file, Form("Title:FAZIA Z %s ID Quality by Group", id.Data()));
154  myCanvas->Clear();
155  myCanvas->Divide(nx, ny);
156  pad = 1;
157  }
158  }
159  }
160  if (pad != 1) {
162  myCanvas->Print(pdf_file, Form("Title:FAZIA Z %s ID Quality by Group", id.Data()));
163  }
164 
165  myCanvas->Clear();
166  myCanvas->Divide(nx, ny);
167  pad = 1;
168  group_num = 1;
169  for (auto& p : fazia_map) {
170  myCanvas->cd(pad);
171 
172  if(FAZIA_group_mean_A_vs_Z(group_num, p.second, id))
173  {
174  ++group_num;
175  ++pad;
176  if (pad > nx * ny) {
177  myCanvas->Print(pdf_file, Form("Title:FAZIA <A> vs Z %s by Group", id.Data()));
178  myCanvas->Clear();
179  myCanvas->Divide(nx, ny);
180  pad = 1;
181  }
182  }
183  }
184  if (pad != 1) myCanvas->Print(pdf_file, Form("Title:FAZIA <A> vs Z %s by Group", id.Data()));
185 
186  myCanvas->Clear();
187  myCanvas->Divide(nx, ny);
188  pad = 1;
189  group_num = 1;
190  for (auto& p : fazia_map) {
191  myCanvas->cd(pad);
192 
193  if(FAZIA_group_Z_threshold_vs_Z(group_num, p.second, id))
194  {
195  ++group_num;
196  ++pad;
197  if (pad > nx * ny) {
198  myCanvas->Print(pdf_file, Form("Title:FAZIA Z thresh. vs Z %s by Group", id.Data()));
199  myCanvas->Clear();
200  myCanvas->Divide(nx, ny);
201  pad = 1;
202  }
203  }
204  }
205  if (pad != 1) myCanvas->Print(pdf_file, Form("Title:FAZIA Z thresh. vs Z %s by Group", id.Data()));
206 
207  myCanvas->Clear();
208  myCanvas->Divide(nx, ny);
209  pad = 1;
210  group_num = 1;
211  for (auto& p : fazia_map) {
212  myCanvas->cd(pad);
213 
214  if(FAZIA_group_A_threshold_vs_Z(group_num, p.second, id))
215  {
216  ++group_num;
217  ++pad;
218  if (pad > nx * ny) {
219  myCanvas->Print(pdf_file, Form("Title:FAZIA A thresh. vs Z %s by Group", id.Data()));
220  myCanvas->Clear();
221  myCanvas->Divide(nx, ny);
222  pad = 1;
223  }
224  }
225  }
226  if (pad != 1) myCanvas->Print(pdf_file, Form("Title:FAZIA A thresh. vs Z %s by Group", id.Data()));
227  }
228  myCanvas->Clear();
229  myCanvas->Print(last_page, "Title:Last page");
230 }
231 
232 
233 
235 
236 std::optional<KVUnownedList> KVDataQualityAuditReporting_INDRAFAZIA::get_indra_telescopes(int ring, const TString& idtype)
237 {
238  bool ok = false;
239  KVUnownedList tels;
240  for (int mod = 1; mod <= 24; ++mod) {
241  TString name = Form("%s_%02d%02d", idtype.Data(), ring, mod);
242  auto tel = gIndra->GetIDTelescope(name);
243  if (tel) {
244  tels.Add(tel);
245  if(fAudit->HasTelescope(name))
246  ok = true;
247  }
248  }
249  if(ok)
250  return tels;
251  return {};
252 }
253 
254 
255 
257 
259 {
260  auto tels = get_indra_telescopes(ring, idtype);
261 
262  if(!tels)
263  return false;
264 
265  auto get_module_number = [](const KVIDTelescope * idt) {
266  return (int)dynamic_cast<KVINDRADetector*>(idt->GetDetector(1))->GetModuleNumber();
267  };
268  auto zmean = fReport.get_mean_Z_for_telescopes(&(*tels), get_module_number);
269  auto zmax = fReport.get_max_Z_for_telescopes(&(*tels), get_module_number, 24);
270  auto zmax_iso = fReport.get_max_Z_with_isotopes_for_telescopes(&(*tels), get_module_number, 30);
271  zmax_iso->SetMarkerColor(kOrange - 2);
272  auto zmin = fReport.get_min_Z_for_telescopes(&(*tels), get_module_number, 25);
273  TMultiGraph* mg = ::new TMultiGraph;
274  mg->SetTitle(Form("Ring %d %s Min/Mean/Max Z vs. Module", ring, idtype.Data()));
275  mg->Add(zmin);
276  mg->Add(zmean);
277  mg->Add(zmax);
278  mg->Add(zmax_iso);
279  mg->Draw("ap");
280 
281  return true;
282 }
283 
284 
286 
287 void KVDataQualityAuditReporting_INDRAFAZIA::INDRA_ring_mean_A_vs_Z(int ring, const TString& idtype, int& pad, int nx, int ny)
288 {
289  auto tels = get_indra_telescopes(ring,idtype);
290  if(!tels)
291  return;
292 
293  int nmods = tels->GetEntries();
294  // 24 => 4 sets of 6, 16 => 2 sets of 8, 8 => 1 set of 8
295  int mod_set = (nmods > 16 ? 6 : 8);
296  int color_step = TColor::GetPalette().GetSize() / (mod_set + 1);
297  KVIDTelescope* idt;
298  nmods = 0;
299  while (tels->GetEntries()) {
300  TMultiGraph* mg = ::new TMultiGraph;
301  int i = 1;
302  while (i <= mod_set) {
303  idt = (KVIDTelescope*)tels->Remove(tels->First()); // "pop" the first one in the list
304  if (fAudit->HasTelescope(idt->GetName())) {
305  auto gr = fReport[idt->GetName()].get_mean_isotopic_mass_by_Z();
306  gr->SetMarkerColor(TColor::GetPalette()[color_step * i]);
308  gr->SetLineWidth(0);
309  mg->Add(gr);
310  }
311  ++i;
312  ++nmods;
313  }
314  myCanvas->cd(pad);
315  ++pad;
316  // change title of graph => title of pad
317  mg->SetTitle(Form("INDRA <A> vs. Z %s Ring %d [%d-%d]", idtype.Data(), ring, nmods - mod_set + 1, nmods));
318  mg->Draw("ap");
319  TLegend* leg;
320  if (i > 12) {
321  leg = gPad->BuildLegend(.11, .89, .61, .69);
322  leg->SetNColumns(3);
323  }
324  else if (i > 6) {
325  leg = gPad->BuildLegend(.11, .89, .61, .69);
326  leg->SetNColumns(2);
327  }
328  else {
329  leg = gPad->BuildLegend(.11, .89, .61, .69);
330  }
331  leg->SetBorderSize(0);
332  leg->SetFillColorAlpha(kWhite, 1.00);
333  if (pad > nx * ny) {
334  myCanvas->Print(current_page, Form("Title:INDRA <A> vs Z %s by Ring", idtype.Data()));
335  myCanvas->Clear();
336  myCanvas->Divide(nx, ny);
337  pad = 1;
338  }
339  }
340 }
341 
342 
343 
345 
346 void KVDataQualityAuditReporting_INDRAFAZIA::INDRA_ring_Z_threshold_vs_Z(int ring, const TString& idtype, int& pad, int nx, int ny)
347 {
348  auto tels = get_indra_telescopes(ring,idtype);
349 
350  if(!tels)
351  return;
352 
353  int nmods = tels->GetEntries();
354  // 24 => 4 sets of 6, 16 => 2 sets of 8, 8 => 1 set of 8
355  int mod_set = (nmods > 16 ? 6 : 8);
356  int color_step = TColor::GetPalette().GetSize() / (mod_set + 1);
357  KVIDTelescope* idt;
358  nmods = 0;
359  while (tels->GetEntries()) {
360  TMultiGraph* mg = ::new TMultiGraph;
361  int i = 1;
362  while (i <= mod_set) {
363  idt = (KVIDTelescope*)tels->Remove(tels->First()); // "pop" the first one in the list
364  if (fAudit->HasTelescope(idt->GetName())) {
365  auto gr = fReport[idt->GetName()].get_element_thresholds_by_Z_mev_per_nuc(TColor::GetPalette()[color_step * i]);
367  gr->SetLineWidth(0);
368  mg->Add(gr);
369  if (gr->GetMean(2) < 0) {
370  // threshold = -1 => uncalibrated particles
371  std::cout << idt->GetName() << " : particles are UNCALIBRATED";
372  // check detector calibrations
373  KVGeoDNTrajectory* traj;
374  if (idt->GetSize() > 1) traj = (KVGeoDNTrajectory*)idt->GetDetector(2)->GetNode()->GetForwardTrajectories()->First();
375  else traj = (KVGeoDNTrajectory*)(idt->GetDetector(1)->GetNode()->GetForwardTrajectories() ? idt->GetDetector(1)->GetNode()->GetForwardTrajectories()->First() : nullptr);
376  if(traj){
377  int nuncal = 0;
378  traj->IterateFrom();
379  KVGeoDetectorNode* dn;
380  while ((dn = traj->GetNextNode())) nuncal += (!dn->GetDetector()->IsCalibrated());
381  if (nuncal == traj->GetN()) std::cout << " ... just like ALL DETECTORS";
382  else {
383  std::cout << " ... just like ";
384  traj->IterateFrom();
385  KVGeoDetectorNode* dn;
386  while ((dn = traj->GetNextNode())) {
387  if (!dn->GetDetector()->IsCalibrated()) {
388  std::cout << dn->GetName() << " ";
389  }
390  }
391  }
392  }
393  std::cout << std::endl;
394  }
395  }
396  ++i;
397  ++nmods;
398  }
399  myCanvas->cd(pad);
400  ++pad;
401  // change title of graph => title of pad
402  mg->SetTitle(Form("INDRA Z thresh. [MeV/u] vs. Z %s Ring %d [%d-%d]", idtype.Data(), ring, nmods - mod_set + 1, nmods));
403  mg->Draw("ap");
404  double x1, x2;
405  // if(idtype == "CSI")
406  // {
407  // return;// no legend - doesn't fit
408  // //x1 = .39; x2 = .89;
409  // }
410  // else
411  // {
412  x1 = .11;
413  x2 = .61;
414  // }
415  TLegend* leg;
416  if (mod_set > 12) {
417  leg = gPad->BuildLegend(x1, .89, x2, .69);
418  leg->SetNColumns(3);
419  }
420  else if (mod_set > 6) {
421  leg = gPad->BuildLegend(x1, .89, x2, .69);
422  leg->SetNColumns(2);
423  }
424  else {
425  leg = gPad->BuildLegend(x1, .89, x2, .69);
426  }
427  leg->SetBorderSize(0);
428  leg->SetFillColorAlpha(kWhite, 1.00);
429  if (pad > nx * ny) {
430  myCanvas->Print(current_page, Form("Title:INDRA Z thresh. vs Z %s by Ring", idtype.Data()));
431  myCanvas->Clear();
432  myCanvas->Divide(nx, ny);
433  pad = 1;
434  }
435  }
436 }
437 
438 
439 
442 
443 bool KVDataQualityAuditReporting_INDRAFAZIA::fill_telescopes_of_group(TList& tels, std::vector<KVDetector*>& dets, const TString& idtype, double& theta_min, double& theta_max)
444 {
445  // \returns true if telescopes added to list, false if not a single telescope appears in the audit
446  bool ok = false;
447  for (auto d : dets) {
448  int index = dynamic_cast<KVFAZIADetector*>(d)->GetIndex();
449  TString name = Form("ID_%s_%d", idtype.Data(), index);
450  auto tel = gFazia->GetIDTelescope(name);
451  if (tel) {
452  if (fAudit->HasTelescope(tel->GetName())) {
453  tels.Add(tel);
454  ok = true;
455  }
456  else {
457  if (tel->IsReadyForID()) std::cout << tel->GetName() << " is absent from audit - BUT IS READY TO IDENTIFY!" << std::endl;
458  }
459  }
460  theta_min = std::min(theta_min, d->GetTheta());
461  theta_max = std::max(theta_min, d->GetTheta());
462  }
463  return ok;
464 }
465 
466 
467 
470 
471 bool KVDataQualityAuditReporting_INDRAFAZIA::FAZIA_group_reporting_Z(int group_num, std::vector<KVDetector*>& dets, const TString& idtype)
472 {
473  // \returns false if nothing was displayed
474 
475  TList tels;
476  double theta_min{360}, theta_max{0};
477  Info("FAZIA_group_reporting_Z", "id=%s Group %d", idtype.Data(), group_num);
478  if(!fill_telescopes_of_group(tels, dets, idtype, theta_min, theta_max))
479  {
480  Info("FAZIA_group_reporting_Z", "id=%s Group %d: NO TELESCOPES IN AUDIT", idtype.Data(), group_num);
481  return false;
482  }
483 
484  int index = 0;
485  auto get_index = [&](const KVIDTelescope*) {
486  return index++;
487  };
488  auto zmean = fReport.get_mean_Z_for_telescopes(&tels, get_index);
489  index = 0;
490  auto zmax = fReport.get_max_Z_for_telescopes(&tels, get_index, 24);
491  index = 0;
492  auto zmax_iso = fReport.get_max_Z_with_isotopes_for_telescopes(&tels, get_index, 30);
493  zmax_iso->SetMarkerColor(kOrange - 2);
494  index = 0;
495  auto zmin = fReport.get_min_Z_for_telescopes(&tels, get_index, 25);
496  TMultiGraph* mg = ::new TMultiGraph;
497  mg->SetTitle(Form("Group %d [%.2f#leq#theta#leq%.2f] %s Min/Mean/Max Z", group_num, theta_min, theta_max, idtype.Data()));
498  mg->Add(zmin);
499  mg->Add(zmean);
500  mg->Add(zmax);
501  mg->Add(zmax_iso);
502  relabel_FAZIA_telescope_axis(mg, &tels);
503  mg->Draw("ap");
504  return true;
505 }
506 
507 
508 
510 
511 bool KVDataQualityAuditReporting_INDRAFAZIA::FAZIA_group_mean_A_vs_Z(int group_num, std::vector<KVDetector*>& dets, const TString& idtype)
512 {
513  TList tels;
514  double theta_min{360}, theta_max{0};
515  if(!fill_telescopes_of_group(tels, dets, idtype, theta_min, theta_max))
516  return false;
517 
518  int color_step = TColor::GetPalette().GetSize() / (dets.size() + 1);
519  int i = 1;
520  TIter next(&tels);
521  KVIDTelescope* idt;
522  TMultiGraph* mg = ::new TMultiGraph;
523  while ((idt = (KVIDTelescope*)next())) {
524  auto gr = fReport[idt->GetName()].get_mean_isotopic_mass_by_Z();
525  gr->SetMarkerColor(TColor::GetPalette()[color_step * i]);
527  gr->SetLineWidth(0);
528  mg->Add(gr);
529  ++i;
530  }
531  // change title of graph => title of pad
532  mg->SetTitle(Form("FAZIA <A> vs. Z %s Group %d [%.2f#leq#theta#leq%.2f]", idtype.Data(), group_num, theta_min, theta_max));
533  mg->Draw("ap");
534  TLegend* leg;
535  if (i > 12) {
536  leg = gPad->BuildLegend(.11, .89, .61, .69);
537  leg->SetNColumns(3);
538  }
539  else if (i > 6) {
540  leg = gPad->BuildLegend(.11, .89, .61, .69);
541  leg->SetNColumns(2);
542  }
543  else {
544  leg = gPad->BuildLegend(.11, .89, .61, .69);
545  }
546  leg->SetBorderSize(0);
547  leg->SetFillColorAlpha(kWhite, 1.00);
548  return true;
549 }
550 
551 
552 
554 
555 bool KVDataQualityAuditReporting_INDRAFAZIA::FAZIA_group_Z_threshold_vs_Z(int group_num, std::vector<KVDetector*>& dets, const TString& idtype)
556 {
557  TList tels;
558  double theta_min{360}, theta_max{0};
559  if(!fill_telescopes_of_group(tels, dets, idtype, theta_min, theta_max))
560  return false;
561 
562  int color_step = TColor::GetPalette().GetSize() / (dets.size() + 1);
563  int i = 1;
564  TIter next(&tels);
565  KVIDTelescope* idt;
566  TMultiGraph* mg = ::new TMultiGraph;
567  while ((idt = (KVIDTelescope*)next())) {
568  auto gr = fReport[idt->GetName()].get_element_thresholds_by_Z_mev_per_nuc(TColor::GetPalette()[color_step * i]);
570  gr->SetLineWidth(0);
571  mg->Add(gr);
572  if (gr->GetMean(2) < 0) {
573  // threshold = -1 => uncalibrated particles
574  std::cout << idt->GetName() << " : particles are UNCALIBRATED";
575  // check detector calibrations
576  KVGeoDNTrajectory* traj = nullptr;
577  if (idt->GetSize() > 1) traj = (KVGeoDNTrajectory*)idt->GetDetector(2)->GetNode()->GetForwardTrajectories()->First();
578  else {
579  if (idt->GetDetector(1)->GetNode()->GetForwardTrajectories())
581  else {
582  if (idt->GetDetector(1)->IsCalibrated()) std::cout << " ... but " << idt->GetDetector(1)->GetName() << " is calibrated...";
583  else std::cout << " ... just like " << idt->GetDetector(1)->GetName();
584  }
585  }
586  if (traj) {
587  int nuncal = 0;
588  traj->IterateFrom();
589  KVGeoDetectorNode* dn;
590  while ((dn = traj->GetNextNode())) nuncal += (!dn->GetDetector()->IsCalibrated());
591  if (nuncal == traj->GetN()) std::cout << " ... just like ALL DETECTORS";
592  else {
593  std::cout << " ... just like ";
594  traj->IterateFrom();
595  KVGeoDetectorNode* dn;
596  while ((dn = traj->GetNextNode())) {
597  if (!dn->GetDetector()->IsCalibrated()) {
598  std::cout << dn->GetName() << " ";
599  }
600  }
601  }
602  }
603  std::cout << std::endl;
604  }
605  ++i;
606  }
607  // change title of graph => title of pad
608  mg->SetTitle(Form("FAZIA Z thresh. [MeV/u] vs. Z %s Group %d [%.2f#leq#theta#leq%.2f]", idtype.Data(), group_num, theta_min, theta_max));
609  mg->Draw("ap");
610  double x1, x2;
611  if (idtype == "CSI") {
612  return true;// no legend - doesn't fit
613  //x1 = .39; x2 = .89;
614  }
615  else {
616  x1 = .11;
617  x2 = .61;
618  }
619  TLegend* leg;
620  if (i > 12) {
621  leg = gPad->BuildLegend(x1, .89, x2, .69);
622  leg->SetNColumns(3);
623  }
624  else if (i > 6) {
625  leg = gPad->BuildLegend(x1, .89, x2, .69);
626  leg->SetNColumns(2);
627  }
628  else {
629  leg = gPad->BuildLegend(x1, .89, x2, .69);
630  }
631  leg->SetBorderSize(0);
632  leg->SetFillColorAlpha(kWhite, 1.00);
633  return true;
634 }
635 
636 
637 
639 
640 bool KVDataQualityAuditReporting_INDRAFAZIA::FAZIA_group_A_threshold_vs_Z(int group_num, std::vector<KVDetector*>& dets, const TString& idtype)
641 {
642  TList tels;
643  double theta_min{360}, theta_max{0};
644  if(!fill_telescopes_of_group(tels, dets, idtype, theta_min, theta_max))
645  return false;
646 
647  int color_step = TColor::GetPalette().GetSize() / (dets.size() + 1);
648  int i = 1;
649  TIter next(&tels);
650  KVIDTelescope* idt;
651  TMultiGraph* mg = ::new TMultiGraph;
652  while ((idt = (KVIDTelescope*)next())) {
653  auto gr = fReport[idt->GetName()].get_isotope_thresholds_by_Z_mev_per_nuc(TColor::GetPalette()[color_step * i]);
655  gr->SetLineWidth(0);
656  mg->Add(gr);
657  ++i;
658  }
659  // change title of graph => title of pad
660  mg->SetTitle(Form("FAZIA A thresh. [MeV/u] vs. Z %s Group %d [%.2f#leq#theta#leq%.2f]", idtype.Data(), group_num, theta_min, theta_max));
661  mg->Draw("ap");
662  double x1, x2;
663  if (idtype == "CSI") {
664  return true;// no legend - doesn't fit
665  //x1 = .39; x2 = .89;
666  }
667  else {
668  x1 = .11;
669  x2 = .61;
670  }
671  TLegend* leg;
672  if (i > 12) {
673  leg = gPad->BuildLegend(x1, .89, x2, .69);
674  leg->SetNColumns(3);
675  }
676  else if (i > 6) {
677  leg = gPad->BuildLegend(x1, .89, x2, .69);
678  leg->SetNColumns(2);
679  }
680  else {
681  leg = gPad->BuildLegend(x1, .89, x2, .69);
682  }
683  leg->SetBorderSize(0);
684  leg->SetFillColorAlpha(kWhite, 1.00);
685  return true;
686 }
687 
688 
689 
691 
693 {
694  auto N = dynamic_cast<TGraph*>(graf->GetListOfGraphs()->First())->GetN();
695  for (int i = 0; i < N; ++i) {
696  auto bin = graf->GetXaxis()->FindBin(i);
697  graf->GetXaxis()->SetBinLabel(bin, tels->At(i)->GetName());
698  }
699  graf->GetXaxis()->LabelsOption("u");// automatic rotation of labels to fit
700 }
701 
702 
703 
707 
709 {
710  //=========Macro generated from canvas: c1/c1
711  //========= (Tue Jul 27 11:21:28 2021) by ROOT version 6.24/02
712 
713  myCanvas->cd();
714  // make sure all pads are transparent
716  TObject* obj;
717  while ((obj = it())) {
718  if (obj->InheritsFrom("TPad")) {
719  dynamic_cast<TPad*>(obj)->SetFillStyle(4000);
720  }
721  }
722  TLatex* tex = ::new TLatex(0.02336825, 0.86604, "Max. Z");
723  tex->SetTextSize(0.028);
724  tex->SetTextAngle(90);
725  tex->SetLineWidth(2);
726  tex->Draw();
727  TMarker* marker = ::new TMarker(0.01611604, 0.8519389, 4);
728 
729  Int_t ci; // for color index setting
730  TColor* color; // for color definition with alpha
731  ci = TColor::GetColor("#333399");
732  marker->SetMarkerColor(ci);
733  marker->SetMarkerStyle(4);
734  marker->SetMarkerSize(1.7);
735  marker->Draw();
736  marker = ::new TMarker(0.01772764, 0.4453584, 30);
737 
738  ci = TColor::GetColor("#ffcc33");
739  marker->SetMarkerColor(ci);
740  marker->SetMarkerStyle(30);
741  marker->SetMarkerSize(1.7);
742  marker->Draw();
743  tex = ::new TLatex(0.02497985, 0.4629847, "Max. Z with identified A");
744  tex->SetTextSize(0.028);
745  tex->SetTextAngle(90);
746  tex->SetLineWidth(2);
747  tex->Draw();
748  tex = ::new TLatex(0.02659146, 0.2632197, "Mean Z");
749  tex->SetTextSize(0.028);
750  tex->SetTextAngle(90);
751  tex->SetLineWidth(2);
752  tex->Draw();
753  tex = ::new TLatex(0.02820306, 0.1045828, "Min. Z");
754  tex->SetTextSize(0.028);
755  tex->SetTextAngle(90);
756  tex->SetLineWidth(2);
757  tex->Draw();
758  marker = ::new TMarker(0.01853344, 0.2479436, 20);
759 
760  ci = TColor::GetColor("#333399");
761  marker->SetMarkerColor(ci);
762  marker->SetMarkerStyle(20);
763  marker->SetMarkerSize(1.7);
764  marker->Draw();
765  marker = ::new TMarker(0.02095085, 0.08578143, 25);
766 
767  ci = TColor::GetColor("#333399");
768  marker->SetMarkerColor(ci);
769  marker->SetMarkerStyle(25);
770  marker->SetMarkerSize(1.7);
771  marker->Draw();
772 }
773 
774 
775 
779 
781 {
782  //=========Macro generated from canvas: c1/c1
783  //========= (Tue Jul 27 12:21:27 2021) by ROOT version 6.24/02
784  myCanvas->cd();
785  // make sure all pads are transparent
787  TObject* obj;
788  while ((obj = it())) {
789  if (obj->InheritsFrom("TPad")) {
790  dynamic_cast<TPad*>(obj)->SetFillStyle(4000);
791  }
792  }
793  TLatex* tex = new TLatex(0.9813084, 0.9470954, "Max. Z");
794  tex->SetTextSize(0.028);
795  tex->SetTextAngle(270);
796  tex->SetLineWidth(2);
797  tex->Draw();
798  TMarker* marker = new TMarker(0.9853972, 0.9688797, 4);
799 
800  Int_t ci; // for color index setting
801  TColor* color; // for color definition with alpha
802  ci = TColor::GetColor("#333399");
803  marker->SetMarkerColor(ci);
804  marker->SetMarkerStyle(4);
805  marker->SetMarkerSize(1.7);
806  marker->Draw();
807  marker = new TMarker(0.984229, 0.8309129, 30);
808 
809  ci = TColor::GetColor("#ffcc33");
810  marker->SetMarkerColor(ci);
811  marker->SetMarkerStyle(30);
812  marker->SetMarkerSize(1.7);
813  marker->Draw();
814  tex = new TLatex(0.9807243, 0.8143154, "Max. Z with identified A");
815  tex->SetTextSize(0.028);
816  tex->SetTextAngle(270);
817  tex->SetLineWidth(2);
818  tex->Draw();
819  tex = new TLatex(0.9813084, 0.4419087, "Mean Z");
820  tex->SetTextSize(0.028);
821  tex->SetTextAngle(270);
822  tex->SetLineWidth(2);
823  tex->Draw();
824  tex = new TLatex(0.978972, 0.2479253, "Min. Z");
825  tex->SetTextSize(0.028);
826  tex->SetTextAngle(270);
827  tex->SetLineWidth(2);
828  tex->Draw();
829  marker = new TMarker(0.9853972, 0.4564315, 20);
830 
831  ci = TColor::GetColor("#333399");
832  marker->SetMarkerColor(ci);
833  marker->SetMarkerStyle(20);
834  marker->SetMarkerSize(1.7);
835  marker->Draw();
836  marker = new TMarker(0.9830607, 0.2645228, 25);
837 
838  ci = TColor::GetColor("#333399");
839  marker->SetMarkerColor(ci);
840  marker->SetMarkerStyle(25);
841  marker->SetMarkerSize(1.7);
842  marker->Draw();
843 
844 }
845 
846 
int Int_t
#define d(i)
double Double_t
kOrange
kWhite
#define N
winID w
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t SetFillStyle
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t b
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 style
char name[80]
char * Form(const char *fmt,...)
#define gPad
Bool_t IsLabelled(const Char_t *l) const
Definition: KVBase.h:206
TGraph * get_max_Z_for_telescopes(TSeqCollection *idtels, TelescopeIndexFunction F, Marker_t marker_style=20) const
TGraph * get_mean_Z_for_telescopes(TSeqCollection *idtels, TelescopeIndexFunction F, Marker_t marker_style=20) const
TGraph * get_max_Z_with_isotopes_for_telescopes(TSeqCollection *idtels, TelescopeIndexFunction F, Marker_t marker_style=20) const
TGraph * get_min_Z_for_telescopes(TSeqCollection *idtels, TelescopeIndexFunction F, Marker_t marker_style=20) const
Prepare PDF report on data quality audits for INDRA-FAZIA experiments.
void relabel_FAZIA_telescope_axis(TMultiGraph *graf, const TList *tels) const
bool FAZIA_group_reporting_Z(int group_num, std::vector< KVDetector * > &, const TString &idtype)
std::map< double, std::vector< KVDetector * > > fazia_map
bool FAZIA_group_Z_threshold_vs_Z(int group_num, std::vector< KVDetector * > &dets, const TString &idtype)
bool INDRA_ring_reporting_Z(int ring, const TString &idtype)
bool fill_telescopes_of_group(TList &tels, std::vector< KVDetector * > &dets, const TString &idtype, double &theta_min, double &theta_max)
void make_fazia_map(double theta_bin)
sort fazia telescopes into bins of theta
void INDRA_ring_Z_threshold_vs_Z(int ring, const TString &idtype, int &pad, int nx, int ny)
bool FAZIA_group_mean_A_vs_Z(int group_num, std::vector< KVDetector * > &, const TString &idtype)
void INDRA_ring_mean_A_vs_Z(int ring, const TString &idtype, int &pad, int nx, int ny)
bool FAZIA_group_A_threshold_vs_Z(int group_num, std::vector< KVDetector * > &dets, const TString &idtype)
void make_canvas(canvas_t style=canvas_t::kLandscape)
Make an A4-size canvas.
std::optional< KVUnownedList > get_indra_telescopes(int ring, const TString &idtype)
Bool_t HasTelescope(const TString &tel_name) const
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:160
Bool_t IsCalibrated(const KVNameValueList &params={}) const
Definition: KVDetector.cpp:545
KVGeoDetectorNode * GetNode()
Definition: KVDetector.h:345
Double_t GetTheta() const override
Definition: KVDetector.h:727
Base class for FAZIA detectors.
TString GetSiCsIIDType() const
Definition: KVFAZIA.h:314
Path taken by particles through multidetector geometry.
KVGeoDetectorNode * GetNextNode() const
void IterateFrom(const KVGeoDetectorNode *node0=nullptr) const
Information on relative positions of detectors & particle trajectories.
KVDetector * GetDetector() const
const KVSeqCollection * GetForwardTrajectories() const
const Char_t * GetName() const override
Name of node is same as name of associated detector.
const KVSeqCollection * GetDetectors() const
Base class for all detectors or associations of detectors in array which can identify charged particl...
Definition: KVIDTelescope.h:85
KVDetector * GetDetector(UInt_t n) const
UInt_t GetSize() const
Base class for detectors of INDRA array.
TString GetCsIIDType() const
Definition: KVINDRA.cpp:775
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray", KVExpDB *db=nullptr)
KVIDTelescope * GetIDTelescope(const Char_t *name) const
Return pointer to DeltaE-E ID Telescope with "name".
TObject * First() const override
void Add(TObject *obj) override
Extended TList class which does not own its objects by default.
Definition: KVUnownedList.h:20
Int_t GetSize() const
virtual void SetFillStyle(Style_t fstyle)
virtual void SetLineWidth(Width_t lwidth)
virtual void SetMarkerColor(Color_t mcolor=1)
virtual void SetMarkerStyle(Style_t mstyle=1)
virtual void SetMarkerSize(Size_t msize=1)
virtual void SetTextAngle(Float_t tangle=0)
virtual void SetTextSize(Float_t tsize=1)
virtual Int_t FindBin(const char *label)
virtual void LabelsOption(Option_t *option="h")
virtual void SetBinLabel(Int_t bin, const char *label)
void Clear(Option_t *option="") override
TVirtualPad * cd(Int_t subpadnumber=0) override
static const TArrayI & GetPalette()
static Int_t GetColor(const char *hexcolor)
virtual Double_t GetMean(Int_t axis=1) const
void Add(TObject *obj) override
TObject * First() const override
TObject * At(Int_t idx) const override
void Draw(Option_t *option="") override
TList * GetListOfGraphs() const
TAxis * GetXaxis()
const char * GetName() const override
virtual const char * GetName() const
virtual Bool_t InheritsFrom(const char *classname) const
virtual void Draw(Option_t *option="")
virtual void Info(const char *method, const char *msgfmt,...) const
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
TList * GetListOfPrimitives() const override
void Print(const char *filename, Option_t *option) override
const char * Data() const
void swap(RVec< T > &lhs, RVec< T > &rhs)
TGraphErrors * gr
leg
TH1 * h
Int_t Nint(T x)
TArc a
ClassImp(TPyArg)