KaliVeda
Toolkit for HIC analysis
KVIDTelescope.cpp
1 /***************************************************************************
2 $Id: KVIDTelescope.cpp,v 1.52 2009/05/05 15:54:04 franklan Exp $
3 Author : $Author: franklan $
4  KVIDTelescope.cpp - description
5  -------------------
6  begin : Wed Jun 18 2003
7  copyright : (C) 2003 by J.D Frankland
8  email : frankland@ganil.fr
9  ***************************************************************************/
10 
11 /***************************************************************************
12  * *
13  * This program is free software; you can redistribute it and/or modify *
14  * it under the terms of the GNU General Public License as published by *
15  * the Free Software Foundation; either version 2 of the License, or *
16  * (at your option) any later version. *
17  * *
18  ***************************************************************************/
19 #include "TROOT.h"
20 #include "KVIDTelescope.h"
21 #include "KVGroup.h"
22 #include "KVNucleus.h"
23 #include "KVReconstructedNucleus.h"
24 #include "KVIDGraph.h"
25 #include "KVIDGrid.h"
26 #include "Riostream.h"
27 #include "TPluginManager.h"
28 #include "KVMultiDetArray.h"
29 #include "KVDataSet.h"
30 #include "KVIDGridManager.h"
31 #include "KVIDZALine.h"
32 #include "KVIDCutLine.h"
33 #include "KVIdentificationResult.h"
34 #include "TMath.h"
35 #include "TClass.h"
36 #include "TH2.h"
37 #include "KVParticleCondition.h"
38 
39 #include <KVDetectorSignalExpression.h>
40 #include <KVIDZAGrid.h>
41 
42 using namespace std;
43 
46 
47 
49 
51  : fGroup(nullptr)
52 {
53  init();
54 }
55 
56 
57 
60 
62 {
63  //default init
66 }
67 
68 
69 
101 
103 {
104  // Default initialisation for ID telescopes.
105  // If telescope has at least 1 grid then it is ready to identify
106  // particles after initialising the grid(s) (kReadyForID=true);
107  // otherwise kReadyForID is set to kFALSE, unless the current dataset (if defined)
108  // has been declared to have no associated identification/calibration parameters,
109  // in which case kReadyForID is by default set to kTRUE (for filtering simulations).
110  //
111  // In order to enable mass identification for certain telescopes without a dedicated
112  // implementation (e.g. for simulating array response), put the following lines
113  // in your .kvrootrc:
114  //
115  // [dataset].[telescope label].MassID: yes
116  //
117  // If you want to limit mass identification to certain values of Z and/or A,
118  // add the following line:
119  //
120  // [dataset].[telescope label].MassID.Validity: [expression]
121  //
122  // where [expression] is some valid C++ boolean expression involving Z and/or A,
123  // for example
124  //
125  // [dataset].[telescope label].MassID.Validity: (Z>3)&&(A<20)
126  //
127  //For identifications using more than one grid, the default behaviour is to try identification
128  //with each grid in turn until a successful identification is obtained. The order in which
129  //the grids should be tried should be specified by a variable with the following format:
130  //
131  //~~~~~~~~~~~~~~~~
132  //[Dataset].[telescope label].GridOrder: [Grid1],[Grid2],...
133  //~~~~~~~~~~~~~~~~
134 
136 
137  // for datasets with no calib/ident infos, all id telescopes work
138  if (gDataSet && !gDataSet->HasCalibIdentInfos()) {
140  }
141  else { // for datasets with calib/ident infos, we need a grid & all detectors working
142  // looping over detectors to check they are working
143  // if one of them is not -> set kReadyForID to false
144  TIter it(GetDetectors());
145  KVDetector* det = 0;
146  while ((det = (KVDetector*)it())) if (!det->IsOK()) {
148  return;
149  }
150 
151  if (GetIDGrid()) {
152  KVIDGraph* gr;
153  TIter it(GetListOfIDGrids());
154  bool ok = kTRUE;
155  KVUniqueNameList tmp_list;// for re-ordering grids
156  bool mass_id = false;
157  while ((gr = (KVIDGraph*)it())) {
158  tmp_list.Add(gr);
159  if (gr->HasMassIDCapability()) mass_id = true;
160  gr->Initialize();
161  // make sure both x & y axes' signals are well set up
162  if (!fGraphCoords[gr].fVarX || !fGraphCoords[gr].fVarY) {
163  ok = kFALSE;
164  Warning("Initialize",
165  "ID tel. %s: grid %s has undefined VarX(%s:%p) or VarY(%s:%p) - WILL NOT USE",
166  GetName(), gr->GetName(), gr->GetVarX(), fGraphCoords[gr].fVarX, gr->GetVarY(), fGraphCoords[gr].fVarY);
167  }
168  }
169  // set to true if at least one grid can provide mass identification
170  SetHasMassID(mass_id);
171  // if more than one grid, need to re-order them according to [Dataset].[telescope label].GridOrder
172  if (GetListOfIDGrids()->GetEntries() > 1 && gDataSet) {
173  KVString grid_list = gDataSet->GetDataSetEnv(Form("%s.GridOrder", GetLabel()));
174  ok = kFALSE;
175  if (grid_list == "")
176  Warning("Initialize", "ID telescope %s has %d grids but no %s variable defined",
177  GetName(), GetListOfIDGrids()->GetEntries(), Form("%s.GridOrder", GetLabel()));
178  else if (grid_list.GetNValues(",") != GetListOfIDGrids()->GetEntries())
179  Warning("Initialize", "ID telescope %s has %d grids but %d grids appear in variable %s",
180  GetName(), GetListOfIDGrids()->GetEntries(), grid_list.GetNValues(","), Form("%s.GridOrder", GetLabel()));
181  else {
182  fIDGrids.Clear();
183  grid_list.Begin(",");
184  while (!grid_list.End()) {
185  auto gr_name = grid_list.Next();
186  auto gr_ob = tmp_list.FindObject(gr_name);
187  if (!gr_ob) {
188  Info("Initialize", "IDtel=%s grid %s missing", GetName(), gr_name.Data());
189  }
190  else {
191  fIDGrids.Add(gr_ob);
192  }
193  }
194  ok = kTRUE;
195  }
196  }
197  if (ok) SetBit(kReadyForID);
198  }
199  }
200 
201  if (gDataSet) {
202  SetHasMassID(gDataSet->GetDataSetEnv(Form("%s.MassID", GetLabel()), kFALSE));
203  KVString valid;
204  if ((valid = gDataSet->GetDataSetEnv(Form("%s.MassID.Validity", GetLabel()), "")) != "") {
205  valid.ReplaceAll("Z", "_NUC_->GetZ()");
206  valid.ReplaceAll("A", "_NUC_->GetA()");
207  fMassIDValidity.reset(new KVParticleCondition(valid));
208  }
209  }
210 }
211 
212 
213 
222 
223 void KVIDTelescope::AddDetector(KVDetector* d)
224 {
225  // Add a detector to the telescope.
226  //
227  // Detectors must be added in the order they will be hit by impinging particles,
228  // with the last detector being the one particles stopped in the telescope will stop in.
229  // i.e. dE1, dE2, ..., Eres
230  //
231  // Update name of telescope to "ID_[name of 1st detector]_[name of 2nd detector]_ ... _[name of last detector]"
232 
233  if (d) {
234  fDetectors.Add(d);
235  if (GetSize() > 1) {
236  TString old_name = GetName();
237  old_name += Form("_%s", GetDetectors()->Last()->GetName());
238  SetName(old_name);
239  }
240  else SetName(Form("ID_%s", GetDetector(1)->GetName()));
241  //d->AddIDTelescope(this); <= caused multiple copies to exist in detector's list
242  }
243  else {
244  Warning("AddDetector", "Called with null pointer");
245  }
246 }
247 
248 
249 
253 
255 {
256  // print out telescope structure
257  //if opt="fired" only fired detectors are printed
258 
259  TIter next(GetDetectors());
260  KVDetector* obj;
261 
262  if (!strcmp(opt, "fired")) {
263  while ((obj = (KVDetector*) next())) {
264 
265  if (obj->Fired() || obj->GetEnergy())
266  obj->Print("data");
267  }
268  }
269  else {
270  cout << "\n" << opt << "Structure of KVIDTelescope object: " <<
271  GetName() << " " << GetType() << endl;
272  cout << opt <<
273  "--------------------------------------------------------" <<
274  endl;
275  while ((obj = (KVDetector*) next())) {
276  cout << opt << "Detector: " << obj->GetName() << endl;
277  }
278  }
279 }
280 
281 
282 
283 
286 
287 KVDetector* KVIDTelescope::GetDetector(const Char_t* name) const
288 {
289  // Return a pointer to the detector in the telescope with the name "name".
290 
291  KVDetector* tmp = (KVDetector*) GetDetectors()->FindObject(name);
292  if (!tmp)
293  Warning("GetDetector(const Char_t *name)",
294  "Detector %s not found in telescope %s", name, GetName());
295  return tmp;
296 }
297 
298 
299 
300 
302 
304 {
305  return fGroup;
306 }
307 
308 
309 
310 
312 
314 {
315  fGroup = kvg;
316 }
317 
318 
319 
321 
323 {
324  return (GetGroup() ? GetGroup()->GetNumber() : 0);
325 }
326 
327 
328 
351 
353 {
354  //Default identification method.
355  //
356  //Works for ID telescopes for which one or more identification grids are defined, each
357  //with VARX/VARY parameters corresponding to a KVDetectorSignal or KVDetectorSignalExpression
358  //associated with one or other of the detectors constituting the telescope.
359  //
360  //For identifications using more than one grid, the default behaviour is to try identification
361  //with each grid in turn until a successful identification is obtained. The order in which
362  //the grids should be tried should be specified by a variable with the following format:
363  //
364  //~~~~~~~~~~~~~~~~
365  //[Dataset].[Array Name].[ID type].GridOrder: [Grid1],[Grid2],...
366  //~~~~~~~~~~~~~~~~
367  //
368  //where the name of each grid is given as "VARY_VARX". If no variable defining the order is found,
369  //the grids will be tried in the order they were found in the file containing the grids for this
370  //telescope.
371  //
372  // The KVIdentificationResult is first Clear()ed; then it is filled with IDtype = GetType()
373  // of this identification telescope, IDattempted = true, and the results of the identification
374  // procedure.
375 
376  idr->Clear();
377  idr->IDattempted = true;
378  idr->SetIDType(GetType());
379 
380  KVIDGraph* grid;
381  TIter it(GetListOfIDGrids());
382  while ((grid = (KVIDGraph*)it())) { //loop over grids in order given by [Dataset].[Array Name].[ID type].GridOrder:
383  Double_t de, e;
384  GetIDGridCoords(e, de, grid, x, y);
385  idr->SetGridName(grid->GetName());
386  if (grid->IsIdentifiable(e, de, &idr->Rejecting_Cut)) {
387  grid->Identify(e, de, idr);
388  if (idr->IDOK) break; // stop on first successful identification
389  }
390  else {
391  // particle rejected by cut in grid. idr->Rejecting_Cut contains its name.
392  idr->IDOK = kFALSE;
394  }
395  }
396  idr->IDcode = GetIDCode();
397 
398  return kTRUE;
399 }
400 
401 
402 
403 
439 
441 {
442  // Add an identification grid to the list of grids used by this telescope.
443  //
444  // If the grid's VARX and VARY parameters are set and contain the names of valid
445  // detector signals (see formatting rules below) they will be used by
446  // GetIDGridXCoord() and GetIDGridYCoord() to return the coordinates
447  // needed to perform particle identification using the grid.
448  //
449  // The name of the grid is set to "VARY_VARX" (just the signal names, not the detector
450  // label part - see below). This value will be stored in the
451  // KVIdentificationResult corresponding to an attempted identification of a
452  // KVReconstructedNucleus by this grid.
453  //
454  // VARX/VARY Formatting
455  //
456  // To be valid, grid VARX/Y parameters should be set as follows:
457  //
458  //~~~~~~~~~~~~~~~~~~
459  // [signal name]
460  // or [det_label]::[signal name]
461  //~~~~~~~~~~~~~~~~~~
462  //
463  // where
464  //
465  //~~~~~~~~~~~~~~~~~~
466  // [det_label] (optional) = detector label i.e. string returned by KVDetector::GetLabel()
467  // method for detector. By default, VARX is assumed to be the Eres detector
468  // or last detector and VARY the DE detector or first detector
469  // [signal_name] = name of a signal defined for the detector, possibly depending
470  // on availability of calibration
471  //
472  // To see all available signals for a detector, use
473  //
474  // KVDetector::GetListOfDetectorSignals()
475  //~~~~~~~~~~~~~~~~~~
476 
477  if (grid) {
478  fIDGrids.Add(grid);
479  KVString det_labels_x, det_labels_y;
480  KVDetectorSignal* xx = GetSignalFromGridVar(grid->GetVarX(), "X", det_labels_x);
481  KVDetectorSignal* yy = GetSignalFromGridVar(grid->GetVarY(), "Y", det_labels_y);
482  GraphCoords gc;
483  gc.fVarX = xx;
484  gc.fDetLabelsX = det_labels_x;
485  gc.fVarY = yy;
486  gc.fDetLabelsY = det_labels_y;
487  fGraphCoords[grid] = gc;
488  TString grid_name;
489  if (xx && yy) {
490  grid_name.Form("%s_%s", yy->GetName(), xx->GetName());
491  grid->SetName(grid_name);
492  }
493  }
494 }
495 
496 
497 
543 
545 {
546  // Deduce & return pointer to detector signal from grid VARX/VARY parameter
547  //
548  // To be valid, grid VARX/Y parameters should be set using mathematical expressions
549  // which use the following references to detector signals for the telescope:
550  //
551  //~~~~~~~~~~~~~~~~~~
552  // [signal name]
553  // OR [det_label]::[signal name]
554  //~~~~~~~~~~~~~~~~~~
555  //
556  // where
557  //
558  //~~~~~~~~~~~~~~~~~~
559  // [signal_name] = name of a signal defined for 1 of the detectors of the telescope
560  // [det_label] = optional detector label i.e. string returned by
561  // KVDetector::GetLabel() method for detector
562  //~~~~~~~~~~~~~~~~~~
563  //
564  // If `[det_label]` is not given, we assume for `VARX` the last (E) detector,
565  // while for `VARY` we assume the first (dE) detector. If this telescope has only
566  // one detector, we use it for both variables.
567  //
568  // To see all available signals for a detector, use
569  //
570  //~~~~~~~~~~~~~~~~~~
571  // KVDetector::GetListOfDetectorSignals()
572  //~~~~~~~~~~~~~~~~~~
573  //
574  // #### Example ####
575  // Imagine a telescope which combines 2 detectors, with labels SI and CSI (in that order, i.e. SI is the dE detector,
576  // CSI is the residual energy detector). The following cases are valid:
577  //
578  //~~~~~
579  // VARX = Energy : use 'Energy' signal of CSI detector (default for VARX)
580  // VARY = ADC : use 'ADC' signal of SI detector (default for VARY)
581  // VARY = CSI::Energy : use 'Energy' signal of CSI detector (overrides default for VARY, SI)
582  // VARX = (Q3-Q3Fast)/(0.8*Q3) : uses 'Q3' and 'Q3Fast' signals of CSI detector (default for VARX)
583  // VARY = 1.5*SI::ADC - CSI::Q3Fast/CSI::Q3 : combination of signals from both detectors
584  //~~~~~
585  //
586  // The 'det_labels' string will be filled with a comma-separated list of the labels of each detector used
587  // in the expressions.
588  //
589  // \note if any signals are not defined, they will be evaluated as zero
590 
591  if (var == "") {
592  Warning("GetSignalFromGridVar",
593  "No VAR%s defined for grid for telescope %s. KVIDTelescope-derived class handling identification must override GetIDMapX/GetIDMapY",
594  axe.Data(), GetName());
595  return nullptr;
596  }
597 
598  KVString dum = var;
599  KVDetector* det = nullptr;
600  KVDetectorSignal* ds(nullptr);
601  KVString sig_type;
602 
603  // check if VARX/Y is a mathematical expression
604  Bool_t is_expression = var.GetNValues("+-*/()") > 1;
605  // examine each term in expression (this will split the elements in a mathematical expression,
606  // or just take the whole expression if no math operators are present)
607  var.Begin("+-*/()");
608  bool explicit_det_ref = false;
609  bool multidetexpr = false;
610  KVString det_label;
611  while (!var.End()) {
612  KVString t = var.Next();
613  // check if we have an explicit reference to a detector
614  if (t.Contains("::")) {
615  t.Begin("::");
616  det_label = t.Next();
617  auto _det = (KVDetector*)GetDetectors()->FindObjectByLabel(det_label);
618  if (_det) {
619  explicit_det_ref = true;
620  sig_type = t.Next();
621  // check signal is defined for detector
622  if (_det->GetDetectorSignal(sig_type)) {
623  if (det_labels.Length()) {
624  if (!det_labels.Contains(det_label)) det_labels += Form(",%s", det_label.Data());
625  }
626  else
627  det_labels = det_label;
628  }
629  else {
630  Warning("GetSignalFromGridVar", "signal '%s' not found in det '%s' -> replaced by '0'", sig_type.Data(), _det->GetName());
631  dum.ReplaceAll(Form("%s::%s", _det->GetLabel(), sig_type.Data()), "0");
632  }
633  }
634  if (det && (_det != det)) multidetexpr = true; // more than 1 detector is referenced
635  det = _det;
636  }
637  }
638  // treat all cases with explicit detector references
639  if (explicit_det_ref) {
640  if (!multidetexpr) {
641  // only 1 detector is directly referenced
642  if (!is_expression) {
643  // it is not a mathematical expression: this is just the name of a detector signal
644  return det->GetDetectorSignal(sig_type);
645  }
646  else {
647  // math expression with explicit reference to signal(s) of 1 detector
648  sig_type = var;
649  // remove all explicit references to detector: 'det' pointer has been set
650  sig_type.ReplaceAll(Form("%s::", det_label.Data()), "");
651  // expression will be handled below
652  }
653  }
654  else {
655  // use specific constructor for explicit detector references
656  ds = new KVDetectorSignalExpression(var, dum, GetDetectors());
657  if (!ds->IsValid()) {
658  delete ds;
659  ds = nullptr;
660  Error("GetSignalFromGridVar",
661  "Problem initialising ID-grid %s coordinate for telescope %s."
662  " Check definition of VAR%s for grid (=%s)",
663  axe.Data(), GetName(), axe.Data(), var.Data());
664  return ds;
665  }
667  return ds;
668  }
669  }
670  else {
671  // no explicit reference to detectors - by default, use detector (1) for Y axis,
672  // and the last detector for X axis
673  if (axe == "Y" || GetSize() == 1) det = GetDetector(1);
674  else det = GetDetector(GetSize());
675  sig_type = var;
676  det_labels = det->GetLabel();
677  }
678  ds = det->GetDetectorSignal(sig_type);
679  if (!ds) {
680  // sig_type is not the name of a known signal: assume it is an expression using known signal names
681  if (!det->AddDetectorSignalExpression(sig_type, sig_type)) {
682  Error("GetSignalFromGridVar",
683  "Problem initialising ID-grid %s coordinate for telescope %s. Request for unknown signal %s for detector %s. Check definition of VAR%s for grid (=%s)",
684  axe.Data(), GetName(), sig_type.Data(), det->GetName(), axe.Data(), var.Data());
685  ds = nullptr;
686  }
687  else
688  ds = det->GetDetectorSignal(sig_type);
689  }
690  return ds;
691 }
692 
693 
694 
696 
698 {
700  if (!cl || !cl->InheritsFrom("KVIDZAGrid")) cl = TClass::GetClass("KVIDZAGrid");
701  KVIDGrid* idgrid = (KVIDGrid*)cl->New();
702 
703  idgrid->AddIDTelescope(this);
704  idgrid->SetOnlyZId(onlyZ);
705  idgrid->SetRunList("1-10000");
706 
707  KVIDCutLine* B_line = (KVIDCutLine*)idgrid->Add("OK", "KVIDCutLine");
708  Int_t npoi_bragg = 0;
709  B_line->SetName("Bragg_line");
710  B_line->SetAcceptedDirection("right");
711 
712  return idgrid;
713 }
714 
715 
716 
724 
725 void KVIDTelescope::addLineToGrid(KVIDGrid* idgrid, int zz, int aa, int npoints)
726 {
727 
728  //loop over energy
729  //first find :
730  // ****E1 = energy at which particle passes 1st detector and starts to enter in the 2nd one****
731  // E2 = energy at which particle passes the 2nd detector
732  //then perform npoints calculations between these two energies and use these
733  //to construct a KVIDZALine
734 
735  double xfactor = 1.;
736 
737  KVNucleus part;
738 
739  KVDetector* det_de = GetDetector(1);
740  KVDetector* det_eres = GetDetector(2);
741 
742  Double_t SeuilE = 0.1;
743 
744  part.SetZ(zz);
745  part.SetA(aa);
746 
747  Double_t E1, E2;
748  //find E1
749  //go from SeuilE MeV to det_de->GetEIncOfMaxDeltaE(part.GetZ(),part.GetA()))
750  Double_t E1min = SeuilE, E1max = det_de->GetEIncOfMaxDeltaE(zz, aa);
751  E1 = (E1min + E1max) / 2.;
752 
753  while ((E1max - E1min) > SeuilE) {
754 
755  part.SetEnergy(E1);
756  det_de->Clear();
757  det_eres->Clear();
758 
759  det_de->DetectParticle(&part);
760  det_eres->DetectParticle(&part);
761  if (det_eres->GetEnergy() > SeuilE) {
762  //particle got through - decrease energy
763  E1max = E1;
764  E1 = (E1max + E1min) / 2.;
765  }
766  else {
767  //particle stopped - increase energy
768  E1min = E1;
769  E1 = (E1max + E1min) / 2.;
770  }
771  }
772 
773  //add point to Bragg line
774  Double_t dE_B = det_de->GetMaxDeltaE(zz, aa);
775  Double_t E_B = det_de->GetEIncOfMaxDeltaE(zz, aa);
776  Double_t Eres_B = det_de->GetERes(zz, aa, E_B);
777 
778  KVIDCutLine* B_line = (KVIDCutLine*)idgrid->GetCut("Bragg_line");
779  if (B_line) B_line->SetPoint(B_line->GetN(), Eres_B, dE_B);
780 
781  //find E2
782  //go from E1 MeV to maximum value where the energy loss formula is valid
783  Double_t E2min = E1, E2max = det_eres->GetEmaxValid(part.GetZ(), part.GetA());
784  E2 = (E2min + E2max) / 2.;
785 
786  while ((E2max - E2min > SeuilE)) {
787 
788  part.SetEnergy(E2);
789  det_de->Clear();
790  det_eres->Clear();
791 
792  det_de->DetectParticle(&part);
793  det_eres->DetectParticle(&part);
794  if (part.GetEnergy() > SeuilE) {
795  //particle got through - decrease energy
796  E2max = E2;
797  E2 = (E2max + E2min) / 2.;
798  }
799  else {
800  //particle stopped - increase energy
801  E2min = E2;
802  E2 = (E2max + E2min) / 2.;
803  }
804  }
805  E2 *= xfactor;
806  if ((!strcmp(det_eres->GetType(), "CSI")) && (E2 > 5000)) E2 = 5000;
807  // printf("z=%d a=%d E1=%lf E2=%lf\n",zz,aa,E1,E2);
808  KVIDZALine* line = (KVIDZALine*)idgrid->Add("ID", "KVIDZALine");
809  if (TMath::Even(zz)) line->SetLineColor(4);
810  line->SetZ(zz);
811  line->SetA(aa);
812 
813  Double_t logE1 = TMath::Log(E1);
814  Double_t logE2 = TMath::Log(E2);
815  Double_t dLog = (logE2 - logE1) / (npoints - 1.);
816 
817  for (Int_t i = 0; i < npoints; i++) {
818  // Double_t E = E1 + i*(E2-E1)/(npoints-1.);
819  Double_t E = TMath::Exp(logE1 + i * dLog);
820 
821  Double_t Eres = 0.;
822  Int_t niter = 0;
823  while (Eres < SeuilE && niter <= 20) {
824  det_de->Clear();
825  det_eres->Clear();
826 
827  part.SetEnergy(E);
828 
829  det_de->DetectParticle(&part);
830  det_eres->DetectParticle(&part);
831 
832  Eres = det_eres->GetEnergy();
833  E += SeuilE;
834  niter += 1;
835  }
836  if (!(niter > 20)) {
837  Double_t dE = det_de->GetEnergy();
838  Double_t gEres, gdE;
839  line->GetPoint(i - 1, gEres, gdE);
840  line->SetPoint(i, Eres, dE);
841 
842  }
843  }
844  //printf("sort de boucle points");
845 
846 
847 }
848 
849 
850 
859 
861 {
862  // Returns a comma-separated list of the labels of the detectors used to determine
863  // the "x" or "y" coordinates of the identification grid(s)
864  //
865  // If there is more than 1 grid and the list is not the same for all grids,
866  // prints a warning message.
867  //
868  // \param[in] axis name of grid axis i.e. "x", "X", "y" or "Y" (case insensitive)
869 
870  KVString _axis = axis;
871  _axis.ToUpper();
872 
873  if (_axis != "X" && _axis != "Y") {
874  Error("GetDetectorLabelsForGridCoord", "Called with illegal axis name '%s'", axis.Data());
875  return "";
876  }
877  KVString lab_list;
878 
879  for (auto& __gc : fGraphCoords) {
880  KVString labs;
881  if (_axis == "X") labs = __gc.second.fDetLabelsX;
882  else labs = __gc.second.fDetLabelsY;
883  if (lab_list.Length() && lab_list != labs) {
884  Error("GetDetectorLabelsForGridCoord", "Grids for telescope %s use different detector types for %s-coordinate: "
885  "%s and %s", GetName(), _axis.Data(), lab_list.Data(), labs.Data());
886  return "";
887  }
888  lab_list = labs;
889  }
890  return lab_list;
891 }
892 
893 
894 
895 
899 
901 {
902  //Return the first in the list of identification grids used by this telescope
903  //(this is for backwards compatibility with ID telescopes which had only one grid).
904  return (KVIDGraph*)GetListOfIDGrids()->First();
905 }
906 
907 
908 
909 
912 
914 {
915  //Return pointer to grid using position in list. First grid has index = 1.
916  if (index < 1) {
917  Error("GetIDGrid(int)", "Index must be >=1!");
918  return nullptr;
919  }
920  return (KVIDGraph*)GetListOfIDGrids()->At(index - 1);
921 }
922 
923 
924 
925 
929 
931 {
932  //Return pointer to grid using "label" to search in list of grids associated
933  //to this telescope.
934  return (KVIDGraph*)GetListOfIDGrids()->FindObjectByLabel(label);
935 }
936 
937 
938 
939 
941 
943 {
944  AbstractMethod("GetIDMapX");
945  return -1.;
946 }
947 
948 
949 
954 
956 {
957  // Returns the pedestal associated with the 2nd detector of the telescope,
958  // optionally depending on the given option string.
959  // By default this returns 0, and should be overridden in specific implementations.
960 
961  return 0.;
962 }
963 
964 
965 
970 
972 {
973  // Returns the pedestal associated with the 1st detector of the telescope,
974  // optionally depending on the given option string.
975  // By default this returns 0, and should be overridden in specific implementations.
976 
977  return 0.;
978 }
979 
980 
981 
985 
987 {
988  // Return value of X coordinate to be used with the given ID grid
989  // This corresponds to whatever was given as parameter "VARX" for the grid
990 
991  KVDetectorSignal* ds = fGraphCoords[g].fVarX;
992  if (ds) return ds->GetValue();
993  return -1;
994 }
995 
996 
997 
1001 
1003 {
1004  // Return value of Y coordinate to be used with the given ID grid
1005  // This corresponds to whatever was given as parameter "VARY" for the grid
1006 
1007  KVDetectorSignal* ds = fGraphCoords[g].fVarY;
1008  if (ds) return ds->GetValue();
1009  return -1;
1010 }
1011 
1012 
1013 
1014 
1016 
1018 {
1019  AbstractMethod("GetIDMapY");
1020  return -1.;
1021 }
1022 
1023 
1024 
1025 
1029 
1031 {
1032  //Remove all identification grids for this ID telescope
1033  //Grids are not deleted as this is handled by gIDGridManager
1034  fIDGrids.Clear();
1035  fGraphCoords.clear();
1036 }
1037 
1038 
1039 
1040 
1059 
1061 {
1062  //Static function which will create an instance of the KVIDTelescope-derived class
1063  //corresponding to 'name'
1064  //These are defined as 'Plugin' objects in the file $KVROOT/KVFiles/.kvrootrc :
1065  //~~~~~~~
1066  // # The KVMultiDetArray::GetIDTelescopes(KVDetector*de, KVDetector*e) method uses these plugins to
1067  // # create KVIDTelescope instances adapted to the specific array geometry and detector types.
1068  // # For each pair of detectors we look for a plugin with one of the following names:
1069  // # [name_of_dataset].de_detector_type[de detector thickness]-e_detector_type[de detector thickness]
1070  // # Each characteristic in [] brackets may or may not be present in the name; first we test for names
1071  // # with these characteristics, then all combinations where one or other of the characteristics is not present.
1072  // # In addition, we first test all combinations which begin with [name_of_dataset].
1073  // # The first plugin found in this order will be used.
1074  // # In addition, if for one of the two detectors there is a plugin called
1075  // # [name_of_dataset].de_detector_type[de detector thickness]
1076  // # [name_of_dataset].e_detector_type[e detector thickness]
1077  // # then we add also an instance of this 1-detector identification telescope.
1078  //~~~~~~~
1079 
1080  TPluginHandler* ph;
1081  //check and load plugin library
1082  if (!(ph = LoadPlugin("KVIDTelescope", uri)))
1083  return 0;
1084 
1085  //execute constructor/macro for identification telescope
1086  KVIDTelescope* mda = (KVIDTelescope*) ph->ExecPlugin(0);
1087  if (mda) {
1088  //set label of telescope with URI used to find plugin (minus dataset name)
1089  mda->SetLabelFromURI(uri);
1090  }
1091 
1092  return mda;
1093 }
1094 
1095 
1096 
1097 
1101 
1103 {
1104  //PRIVATE METHOD
1105  //Sets label of telescope based on URI of plugin describing child class for this telescope
1106 
1107  TString _uri(uri);
1108  if (gDataSet && _uri.BeginsWith(gDataSet->GetName())) _uri.Remove(0, strlen(gDataSet->GetName()) + 1);
1109  SetLabel(_uri.Data());
1110 }
1111 
1112 
1113 
1114 
1147 
1149 {
1150  // Initialise the identification parameters (grids, etc.) of ALL identification telescopes of this
1151  // kind (label) in the multidetector array. Therefore this method need only be called once, and not
1152  // called for each telescope. The kind/label (returned by GetLabel) of the telescope corresponds
1153  // to the URI used to find the plugin class in $KVROOT/KVFiles/.kvrootrc.
1154  //
1155  // By default, this method first looks for a file containing the list of all files containing the
1156  // identification grids. The default name of this file is:
1157  //
1158  //~~~
1159  // [label].IdentificationFiles.dat
1160  //~~~
1161  //
1162  // or any name given by the environment variable
1163  //
1164  //~~~
1165  // [dataset].IdentificationParameterList.[label]: [filename]
1166  //~~~
1167  //
1168  // If neither file is found, we look for a file with name:
1169  //
1170  //~~~
1171  // IDGrids.[label].dat
1172  //~~~
1173  //
1174  // or any name given by the environment variable
1175  //
1176  //~~~
1177  // [dataset].IdentificationParameterFile.[label]: [filename]
1178  //~~~
1179  //
1180  // which is assumed to contain all identification grids for this type of telescope.
1181 
1183  if(!ok) ok = LoadIdentificationParameters(multidet);
1184 
1185  if (!ok)
1186  Warning("SetIdentificationParameters", "No files found. See doc for this method.");
1187 
1188  return ok;
1189 }
1190 
1191 
1192 
1193 
1211 
1213 {
1214  // In the case where the identification grids are stored in several files, this method looks
1215  // for a file with name
1216  //
1217  //~~~
1218  // [label].IdentificationFiles.dat
1219  //~~~
1220  //
1221  // or any name given by the environment variable
1222  //
1223  //~~~
1224  // [dataset].IdentificationParameterList.[label]: [filename]
1225  //~~~
1226  //
1227  // and if found reads from it the list of files containing the identification grids.
1228  //
1229  // \returns kTRUE if file found and read, kFALSE if no file found
1230 
1231  TString filename = gDataSet->GetDataSetEnv(Form("IdentificationParameterList.%s", GetLabel()),
1232  Form("%s.IdentificationFiles.dat", GetLabel()));
1233 
1234  if(gDataSet->FindDataSetFile(filename))
1235  {
1236  KVFileReader fr;
1237  auto path = gDataSet->GetFullPathToDataSetFile(filename);
1238  Info("ReadIdentificationParameterFiles", "Using file %s", path.Data());
1239  fr.OpenFileToRead(path);
1240  Bool_t ok = kTRUE;
1241  while (fr.IsOK()) {
1242  fr.ReadLine(0);
1243 
1244  if (fr.GetCurrentLine() != "")
1245  ok &= LoadIdentificationParameters(multidet, fr.GetCurrentLine().Data());
1246  }
1247 
1248  fr.CloseFile();
1249  return ok;
1250  }
1251 
1252  return kFALSE;
1253 }
1254 
1255 
1256 
1257 
1275 
1277 {
1278  // \param[in] _filename optional name of dataset file containing grids for this telescope type
1279  // \returns kTRUE if a file is found and read in
1280  //
1281  // If no value for _filename is given, this method looks for a dataset file named
1282  //
1283  //~~~
1284  // IDGrids.[label].dat
1285  //~~~
1286  //
1287  // or any name given by the environment variable
1288  //
1289  //~~~
1290  // [dataset].IdentificationParameterFile.[label]: [filename]
1291  //~~~
1292  //
1293  // which is assumed to contain all identification grids for this type of telescope.
1294 
1295  TString filename;
1296  if(_filename.IsNull())
1297  {
1298  filename = gDataSet->GetDataSetEnv(Form("IdentificationParameterFile.%s", GetLabel()),
1299  Form("IDGrids.%s.dat",GetLabel()));
1300  }
1301  else
1302  filename = _filename;
1303 
1304  if(gDataSet->FindDataSetFile(filename))
1305  {
1306  auto path = gDataSet->GetFullPathToDataSetFile(filename);
1307  //Read grids from file
1308  Info("LoadIdentificationParameters", "Using file %s", path.Data());
1309  multidet->ReadGridsFromAsciiFile(path);
1310  return kTRUE;
1311  }
1312  else if (!_filename.IsNull())
1313  {
1314  // if the method was called with a definite filename to look for, we assume that this
1315  // file is supposed to exist, therefore if it is not found, there is a problem!
1316  Error("LoadIdentificationParameters",
1317  "File %s not found. Should be in %s",
1318  _filename.Data(), gDataSet->GetDataSetDir());
1319  return kFALSE;
1320  }
1321  return kFALSE;
1322 }
1323 
1324 
1325 
1326 
1338 
1340 {
1341  //Remove identification parameters from telescope in such a way that they
1342  //can subsequently be reset e.g. with a new version.
1343  //This is used by KVMultiDetArray::UpdateIdentifications.
1344  //Child classes with specific SetIdentificationParameters methods should
1345  //also redefine this method in order to remove (destroy) cleanly the objects
1346  //created in SetIdentificationParameters.
1347  //
1348  //This default method takes the list of grids associated to the telescope,
1349  //and for each one: 1) checks if it is still in the gIDGridManager's list
1350  //2) if yes, delete the grid and remove it from gIDGridManager
1351 
1352  TIter next_grid(GetListOfIDGrids());
1353  KVIDGrid* grid;
1354  while ((grid = (KVIDGrid*)next_grid())) {
1355 
1356  if (gIDGridManager->GetGrids()->FindObject(grid)) { //this only works if KVIDTelescope uses TObject:IsEqual method (i.e. compares pointers)
1357 
1358  gIDGridManager->DeleteGrid(grid);
1359 
1360  }
1361  }
1362  //clear list of grids
1363  fIDGrids.Clear();
1365 }
1366 
1367 
1368 
1369 //void KVIDTelescope::CalculateParticleEnergy(KVReconstructedNucleus* nuc)
1370 //{
1371 // // The energy of each particle is calculated as follows:
1372 // //
1373 // // E = dE_1 + dE_2 + ... + dE_N
1374 // //
1375 // // dE_1, dE_2, ... = energy losses measured in each detector through which
1376 // // the particle has passed (or stopped, in the case of dE_N).
1377 // // These energy losses are corrected for (Z,A)-dependent effects
1378 // // such as pulse-heigth defect in silicon detectors, losses in
1379 // // windows of gas detectors, etc.
1380 // //
1381 // // Whenever possible, the energy loss for fired detectors which are uncalibrated
1382 // // or not functioning is calculated. In this case the status returned by GetCalibStatus()
1383 // // will be KVIDTelescope::kCalibStatus_Calculated.
1384 // // If none of the detectors is calibrated, the particle's energy cannot be calculated &
1385 // // the status will be KVIDTelescope::kCalibStatus_NoCalibrations.
1386 // // Otherwise, the status code will be KVIDTelescope::kCalibStatus_OK.
1387 
1388 // //status code
1389 // fCalibStatus = kCalibStatus_NoCalibrations;
1390 
1391 // UInt_t z = nuc->GetZ();
1392 // //uncharged particles
1393 // if (z == 0) return;
1394 
1395 // KVDetector* d1 = GetDetector(1);
1396 // KVDetector* d2 = (GetSize() > 1 ? GetDetector(2) : 0);
1397 // Bool_t d1_cal = d1->IsCalibrated();
1398 // Bool_t d2_cal = (d2 ? d2->IsCalibrated() : kFALSE);
1399 
1400 // //no calibrations
1401 // if (!d1_cal && !d2)
1402 // return;
1403 // if ((d1 && d2) && !d1_cal && !d2_cal)
1404 // return;
1405 
1406 // //status code
1407 // fCalibStatus = kCalibStatus_OK;
1408 
1409 // UInt_t a = nuc->GetA();
1410 
1411 // // particles stopped in first member of telescope
1412 // if (nuc->GetStatus() == 3) {
1413 // if (d1_cal) {
1414 // nuc->SetEnergy(d1->GetCorrectedEnergy(nuc, -1, kFALSE)); //N.B.: transmission=kFALSE because particle stop in d1
1415 // }
1416 // return;
1417 // }
1418 
1419 // Double_t e1, e2, einc;
1420 // e1 = e2 = einc = 0.0;
1421 
1422 // if (!d1_cal) {//1st detector not calibrated - calculate from residual energy in 2nd detector
1423 
1424 // //second detector must exist and have all acquisition parameters fired with above-pedestal value
1425 // if (d2 && d2->Fired("Pall")) e2 = d2->GetCorrectedEnergy(nuc, -1, kFALSE); //N.B.: transmission=kFALSE because particle stop in d2
1426 // if (e2 <= 0.0) {
1427 // // zero energy loss in 2nd detector ? can't do anything...
1428 // fCalibStatus = kCalibStatus_NoCalibrations;
1429 // return;
1430 // }
1431 // //calculate & set energy loss in dE detector
1432 // //N.B. using e2 for the residual energy after detector 1 means
1433 // //that we are assuming the particle stops in detector 2.
1434 // //if this is not true, we will underestimate the energy of the particle.
1435 // e1 = d1->GetDeltaEFromERes(z, a, e2);
1436 // if (e1 < 0.0) e1 = 0.0;
1437 // else {
1438 // d1->SetEnergyLoss(e1);
1439 // d1->SetEResAfterDetector(e2);
1440 // e1 = d1->GetCorrectedEnergy(nuc);
1441 // //status code
1442 // fCalibStatus = kCalibStatus_Calculated;
1443 // }
1444 // }
1445 // else { //1st detector is calibrated too: get corrected energy loss
1446 
1447 // e1 = d1->GetCorrectedEnergy(nuc);
1448 
1449 // }
1450 
1451 // if (d2 && !d2_cal) {//2nd detector not calibrated - calculate from energy loss in 1st detector
1452 
1453 // e1 = d1->GetCorrectedEnergy(nuc);
1454 // if (e1 <= 0.0) {
1455 // // zero energy loss in 1st detector ? can't do anything...
1456 // fCalibStatus = kCalibStatus_NoCalibrations;
1457 // return;
1458 // }
1459 // //calculate & set energy loss in 2nd detector
1460 // e2 = d1->GetEResFromDeltaE(z, a);
1461 // if (e2 < 0.0) e2 = 0.0;
1462 // else {
1463 // e2 = d2->GetDeltaE(z, a, e2);
1464 // d2->SetEnergyLoss(e2);
1465 // e2 = d2->GetCorrectedEnergy(nuc);
1466 // //status code
1467 // fCalibStatus = kCalibStatus_Calculated;
1468 // }
1469 // }
1470 // else if (d2) { //2nd detector is calibrated too: get corrected energy loss
1471 
1472 // e2 = d2->GetCorrectedEnergy(nuc, -1, kFALSE);//N.B.: transmission=kFALSE because particle assumed to stop in d2
1473 // // recalculate corrected energy in first stage using info on Eres
1474 // d1->SetEResAfterDetector(e2);
1475 // e1 = d1->GetCorrectedEnergy(nuc);
1476 // }
1477 
1478 // //incident energy of particle (before 1st member of telescope)
1479 // einc = e1 + e2;
1480 
1481 // Double_t coherence_tolerance = gEnv->GetValue("KVIDTelescope.CoherencyTolerance", 1.05);
1482 // if (coherence_tolerance < 1) coherence_tolerance += 1.00;
1483 
1484 // //Now we have to work our way up the list of detectors from which the particle was
1485 // //reconstructed. For each fired & calibrated detector which is only associated with
1486 // //one particle in the events, we add the corrected measured energy loss
1487 // //to the particle. For uncalibrated, unfired detectors and detectors through which
1488 // //more than one particle has passed, we calculate the corrected energy loss and add it
1489 // //to the particle.
1490 // int ndets = nuc->GetNumDet();
1491 // if (ndets > (int)GetSize()) { //particle passed through other detectors before this idtelesocpe
1492 // //look at detectors not in this id telescope
1493 // int idet = GetSize();//next detector after delta-e member of IDTelescope (stopping detector = 0)
1494 // while (idet < ndets) {
1495 
1496 // KVDetector* det = nuc->GetDetector(idet);
1497 // if (det->Fired() && det->IsCalibrated() && det->GetNHits() == 1) {
1498 // Double_t dE = det->GetEnergy();
1499 // //in order to check if particle was really the only one to
1500 // //hit each detector, we calculate the particle's energy loss
1501 // //from its residual energy. if the measured energy loss is
1502 // //significantly larger, there may be a second particle.
1503 // e1 = det->GetDeltaEFromERes(z, a, einc);
1504 // if (e1 < 0.0) e1 = 0.0;
1505 // det->SetEResAfterDetector(einc);
1506 // dE = det->GetCorrectedEnergy(nuc);
1507 // einc += dE;
1508 // }
1509 // else {
1510 // // Uncalibrated/unfired/multihit detector. Calculate energy loss.
1511 // //calculate energy of particle before detector from energy after detector
1512 // e1 = det->GetDeltaEFromERes(z, a, einc);
1513 // if (e1 < 0.0) e1 = 0.0;
1514 // if (det->GetNHits() > 1) {
1515 // //Info("CalculateParticleEnergy",
1516 // // "Detector %s was hit by %d particles. Calculated energy loss for particle %f MeV",
1517 // // det->GetName(), det->GetNHits(), e1);
1518 // if (!(det->Fired() && det->IsCalibrated())) {
1519 // det->SetEnergyLoss(e1 + det->GetEnergy());// sum up calculated energy losses in uncalibrated detector
1520 // }
1521 // //status code
1522 // fCalibStatus = kCalibStatus_Multihit;
1523 // }
1524 // else if (!det->Fired() || !det->IsCalibrated()) {
1525 // //Info("CalculateParticleEnergy",
1526 // // "Detector %s uncalibrated/not fired. Calculated energy loss for particle %f MeV",
1527 // // det->GetName(), e1);
1528 // det->SetEnergyLoss(e1);
1529 // //status code
1530 // fCalibStatus = kCalibStatus_Calculated;
1531 // }
1532 // det->SetEResAfterDetector(einc);
1533 // e1 = det->GetCorrectedEnergy(nuc, e1);
1534 // einc += e1;
1535 // }
1536 // idet++;
1537 // }
1538 // }
1539 // //einc is now the energy of the particle before crossing the first detector
1540 // nuc->SetEnergy(einc);
1541 //}
1542 
1543 
1544 
1554 
1556 {
1557  // Returns name of default ID grid class for this ID telescope.
1558  // This is defined in a .kvrootrc configuration file by one of the following:
1559  // KVIDTelescope.DefaultGrid:
1560  // KVIDTelescope.DefaultGrid.[type]:
1561  // where [type] is the type of this identification telescope (which is given
1562  // by the character string returned by method GetLabel()... sorry :( )
1563  // If no default grid is defined for the specific type of this telescope,
1564  // the default defined by KVIDTelescope.DefaultGrid is used.
1565 
1566  TString specific;
1567  specific.Form("KVIDTelescope.DefaultGrid.%s", GetLabel());
1568  return gEnv->GetValue(specific.Data(), gEnv->GetValue("KVIDTelescope.DefaultGrid", "KVIDGraph"));
1569 }
1570 
1571 
1572 
1578 
1580 {
1581  //Create a dE-E grid (energy loss in detector 1 versus residual energy in detector 2) for a given list of isotopes
1582  // - AperZ : list of A for each Z. (ex: "1=1-3,2=3-4 5,3=6-8,4=7 9-12"...)
1583  // - npoints : number of points in each generated line
1584  // - xfactor : scales the detector 2 thickness to prolongate the lines
1585 
1586  if (GetSize() <= 1) return 0;
1587  if (!GetDetector(1) || !GetDetector(2)) return 0;
1588  double thickness = GetDetector(2)->GetThickness();
1589  GetDetector(2)->SetThickness(thickness * xfactor);
1590 
1591  Info("CalculateDeltaE_EGrid", "called with KVNameValueList");
1592 
1593  KVIDGrid* idgrid = newGrid(0);
1594 
1595  for (auto& par : AperZ) {
1596  KVString tmp = par.GetName();
1597  int zz = tmp.Atoi();
1598  KVNumberList alist = par.GetString();
1599  for (auto aa : alist) {
1600  addLineToGrid(idgrid, zz, aa, npoints);
1601  }
1602  }
1603 
1604  GetDetector(2)->SetThickness(thickness);
1605 
1606  return idgrid;
1607 
1608 }
1609 
1610 
1611 
1619 
1620 KVIDGrid* KVIDTelescope::CalculateDeltaE_EGrid(const KVNumberList& Zrange, Int_t deltaA, Int_t npoints, Double_t lifetime, UChar_t massformula, Double_t xfactor)
1621 {
1622  //Create a dE-E grid (energy loss in detector 1 versus residual energy in detector 2) for a given list of isotopes
1623  // - Zrange : list of element for which we create lines
1624  // - deltaA : number of isotopes generated for each Z around massformula (ex: deltaA=1, Aref-1 Aref Aref+1)
1625  // - npoints : number of points in each generated line
1626  // - lifetime: remove isotopes with lifetime lower than this value
1627  // - xfactor : scales the detector 2 thickness to prolongate the lines
1628 
1629  if (GetSize() <= 1) return 0;
1630  if (!GetDetector(1) || !GetDetector(2)) return 0;
1631  double thickness = GetDetector(2)->GetThickness();
1632  GetDetector(2)->SetThickness(thickness * xfactor);
1633 
1634  Info("CalculateDeltaE_EGrid", "called with KVNumberList");
1635 
1636  KVIDGrid* idgrid = newGrid(0);
1637  KVNucleus part;
1638 
1639  Zrange.Begin();
1640  while (!Zrange.End()) {
1641  Int_t zz = Zrange.Next();
1642  part.SetZ(zz, massformula);
1643  Int_t aref = part.GetA();
1644  for (Int_t aa = aref - deltaA; aa <= aref + deltaA; aa += 1) {
1645  part.SetA(aa);
1646  if (part.IsKnown() && (part.GetLifeTime() > lifetime)) {
1647  addLineToGrid(idgrid, zz, aa, npoints);
1648 
1649  }
1650  }
1651  }
1652  GetDetector(2)->SetThickness(thickness);
1653  return idgrid;
1654 }
1655 
1656 
1657 
1663 
1665 {
1666  //Create a dE-E grid (energy loss in detector 1 versus residual energy in detector 2) for a given list of isotopes
1667  // - haa_zz : lines will be generated for A,Z filled with 1 in this histogram
1668  // - Zonly : if true, generate only one line per Z with the <A>(Z) of the histogram
1669  // - npoints : number of points in each generated line
1670 
1671  if (GetSize() <= 1) return 0;
1672  if (!GetDetector(1) || !GetDetector(2)) return 0;
1673  double thickness = GetDetector(2)->GetThickness();
1674 
1675  Info("CalculateDeltaE_EGrid", "called with TH2");
1676 
1677  KVIDGrid* idgrid = newGrid(0);
1678  KVNucleus part;
1679 
1680  for (Int_t nx = 1; nx <= haa_zz->GetNbinsX(); nx += 1) {
1681 
1682  Int_t zz = TMath::Nint(haa_zz->GetXaxis()->GetBinCenter(nx));
1683  KVNumberList nlA;
1684  Double_t sumA = 0, sum = 0;
1685  for (Int_t ny = 1; ny <= haa_zz->GetNbinsY(); ny += 1) {
1686  Double_t stat = haa_zz->GetBinContent(nx, ny);
1687  if (stat > 0) {
1688  Double_t val = haa_zz->GetYaxis()->GetBinCenter(ny);
1689  nlA.Add(TMath::Nint(val));
1690  sumA += val * stat;
1691  sum += stat;
1692  }
1693  }
1694  sumA /= sum;
1695  Int_t nA = nlA.GetNValues();
1696  if (nA == 0) {
1697  Warning("CalculateDeltaE_EGrid", "no count for Z=%d", zz);
1698  }
1699  else {
1700  if (Zonly) {
1701  nlA.Clear();
1702  nlA.Add(TMath::Nint(sumA));
1703  }
1704  else {
1705  if (nA == 1) {
1706  Int_t aref = nlA.Last();
1707  nlA.Add(aref - 1);
1708  nlA.Add(aref + 1);
1709  }
1710  }
1711  part.SetZ(zz);
1712  nlA.Begin();
1713  while (!nlA.End()) {
1714  Int_t aa = nlA.Next();
1715  part.SetA(aa);
1716  if (part.IsKnown()) {
1717  addLineToGrid(idgrid, zz, aa, npoints);
1718  }
1719  }
1720  }
1721  }
1722  return idgrid;
1723 }
1724 
1725 
1726 
1742 
1744 {
1745  // Returns the Y-axis value in the 2D identification map containing isotope (Z,A)
1746  // corresponding to either the given X-axis/Eres value or the current X-axis value given by GetIDGridXCoord()
1747  // If no mass information is available, just give Z.
1748  //
1749  // In this (default) implementation this means scanning the ID grids associated with
1750  // this telescope until we find an identification line Z or (Z,A), and then interpolating
1751  // the Y-coordinate for the current X-coordinate value.
1752  //
1753  // Status variable can take one of following values:
1754  //
1755  // KVIDTelescope::kMeanDE_OK all OK
1756  // KVIDTelescope::kMeanDE_XtooSmall X-coordinate is smaller than smallest X-coordinate of ID line
1757  // KVIDTelescope::kMeanDE_XtooLarge X-coordinate is larger than largest X-coordinate of ID line
1758  // KVIDTelescope::kMeanDE_NoIdentifie No identifier found for Z or (Z,A)
1759 
1760  status = kMeanDE_OK;
1761  // loop over grids
1762  TIter next(GetListOfIDGrids());
1763  KVIDGrid* grid;
1764  KVIDLine* idline = 0;
1765  while ((grid = (KVIDGrid*)next())) {
1766  idline = (KVIDLine*)grid->GetIdentifier(Z, A);
1767  if (idline) break;
1768  }
1769  if (!idline) {
1770  status = kMeanDE_NoIdentifier;
1771  return -1.;
1772  }
1773  Double_t x, x1, y1, x2, y2;
1774  x = (Eres < 0 ? GetIDGridXCoord(grid) : Eres);
1775  idline->GetEndPoint(x2, y2);
1776  if (x > x2) {
1777  status = kMeanDE_XtooLarge;
1778  return -1;
1779  }
1780  idline->GetStartPoint(x1, y1);
1781  if (x < x1) {
1782  status = kMeanDE_XtooSmall;
1783  return -1.;
1784  }
1785  return idline->Eval(x);
1786 }
1787 
1788 
1789 
1790 
1799 
1801 {
1802  // Return kTRUE if energy of ION is > minimum incident energy required for identification
1803  // This theoretical limit is defined here to be the incident energy for which the
1804  // dE in the first detector of a dE-E telescope is maximum.
1805  // If EINC>0 it is assumed to be the energy of the ion just before the first detector
1806  // (case where ion would have to pass other detectors before reaching this telescope).
1807  //
1808  // If this is not a dE-E telescope, we return kTRUE by default.
1809 
1810  if (GetSize() < 2) return kTRUE;
1811 
1812  KVDetector* dEdet = GetDetector(1);
1813  Double_t emin = dEdet->GetEIncOfMaxDeltaE(ION->GetZ(), ION->GetA());
1814  if (EINC > 0.0) return (EINC > emin);
1815  return (ION->GetEnergy() > emin);
1816 }
1817 
1818 
1819 
1848 
1850 {
1851  // For filtering simulations
1852  //
1853  // \param IDR identification result object for this telescope
1854  // \param n the simulated particle currently being considered
1855  //
1856  // Set the IDR->Zident and IDR->Aident status of the identification of n.
1857  // In principle this depends on whether this telescope provides mass
1858  // identification or not, but this may depend on the particle's energy.
1859  //
1860  // In order to enable mass identification for certain telescopes without a dedicated
1861  // implementation (e.g. for simulating array response), put the following lines
1862  // in your .kvrootrc:
1863  //
1864  // [dataset].[telescope label].MassID: yes
1865  //
1866  // If you want to limit mass identification to certain values of Z and/or A,
1867  // add the following line:
1868  //
1869  // [dataset].[telescope label].MassID.Validity: [expression]
1870  //
1871  // where [expression] is some valid C++ boolean expression involving Z and/or A,
1872  // for example
1873  //
1874  // [dataset].[telescope label].MassID.Validity: (Z>3)&&(A<20)
1875  //
1876  // Then this expression will be tested here in order to determine particle
1877  // identification status
1878 
1879  KVNucleus test(IDR->Z, IDR->A);
1880  if (!test.IsKnown()) { // reject weird identifications as pile-up (in between lines)
1881  IDR->Zident = false;
1882  IDR->IDquality = 4;
1883  IDR->IDOK = false;
1884  return;
1885  }
1886  IDR->Zident = true;
1887  if (!HasMassID()) {
1888  IDR->Aident = false;
1889  }
1890  else {
1891  if (fMassIDValidity) IDR->Aident = fMassIDValidity->Test(n); // test expression for mass ID validity
1892  else IDR->Aident = true; // no expression set; all nuclei are identified in mass
1893  }
1894 }
1895 
1896 
1897 
1900 
1902 {
1903  // Open IdentificationBilan.dat file with given path
1904 
1906  fgIdentificationBilan = new TEnv(path);
1907 }
1908 
1909 
1910 
1913 
1915 {
1916  // Set status of ID Telescope for given system
1917  if (!(fgIdentificationBilan->GetValue(Form("%s.%s", system.Data(), GetName()), kTRUE))) ResetBit(kReadyForID);
1918 }
1919 
1920 
int Int_t
unsigned int UInt_t
#define d(i)
#define e(i)
bool Bool_t
unsigned char UChar_t
char Char_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
R__EXTERN TEnv * gEnv
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 Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t g
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void gc
Option_t Option_t TPoint TPoint const char y1
char name[80]
char * Form(const char *fmt,...)
void SetLabel(const Char_t *lab)
Definition: KVBase.h:194
virtual const Char_t * GetType() const
Definition: KVBase.h:176
void Error(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1650
const Char_t * GetLabel() const
Definition: KVBase.h:198
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:772
void Warning(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1637
UInt_t GetNumber() const
Definition: KVBase.h:219
const Char_t * GetDataSetDir() const
Definition: KVDataSet.cpp:744
Bool_t HasCalibIdentInfos() const
Definition: KVDataSet.h:229
KVString GetDataSetEnv(const Char_t *type, const Char_t *defval="") const
Definition: KVDataSet.cpp:782
TString GetFullPathToDataSetFile(const Char_t *filename)
Definition: KVDataSet.cpp:1924
static Bool_t FindDataSetFile(const TString &dataset, const Char_t *filename)
Definition: KVDataSet.cpp:1958
Signal output from a mathematical combination of other signals.
Base class for output signal data produced by a detector.
virtual Bool_t IsValid() const
virtual Double_t GetValue(const KVNameValueList &params="") const
Handle reading columns of numeric data in text files.
Definition: KVFileReader.h:121
KVString GetCurrentLine()
Definition: KVFileReader.h:320
void CloseFile()
Definition: KVFileReader.h:237
ReadStatus ReadLine(const KVString &pattern="")
Definition: KVFileReader.h:243
Bool_t IsOK()
Definition: KVFileReader.h:231
Bool_t OpenFileToRead(const KVString &filename)
Definition: KVFileReader.h:210
Group of detectors which can be treated independently of all others in array.
Definition: KVGroup.h:20
Line in ID grid used to delimit regions where no identification is possible.
Definition: KVIDCutLine.h:23
void SetName(const char *name) override
This is redeclared to make it appear in context menus for KVIDCutLines.
Definition: KVIDCutLine.h:84
virtual void SetAcceptedDirection(const Char_t *dir)
Definition: KVIDCutLine.cpp:74
Base class for particle identification in a 2D map.
Definition: KVIDGraph.h:32
void Add(TString, KVIDentifier *)
Definition: KVIDGraph.cpp:838
void SetRunList(const char *runlist)
Definition: KVIDGraph.h:156
void SetName(const char *name) override
Definition: KVIDGraph.h:140
void AddIDTelescope(KVBase *t)
Definition: KVIDGraph.h:407
KVIDentifier * GetCut(const Char_t *name) const
Definition: KVIDGraph.h:280
virtual void Identify(Double_t, Double_t, KVIdentificationResult *) const =0
const Char_t * GetName() const override
Definition: KVIDGraph.cpp:1332
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 SetOnlyZId(Bool_t yes=kTRUE)
Definition: KVIDGraph.cpp:1496
KVList * GetGrids()
void DeleteGrid(KVIDGraph *, Bool_t update=kTRUE)
Abstract base class for 2D identification grids in e.g. (dE,E) maps.
Definition: KVIDGrid.h:74
Base class for lines/cuts used for particle identification in 2D data maps.
Definition: KVIDLine.h:143
void GetStartPoint(Double_t &x, Double_t &y) const
void GetEndPoint(Double_t &x, Double_t &y) const
Base class for all detectors or associations of detectors in array which can identify charged particl...
Definition: KVIDTelescope.h:84
KVIDGrid * newGrid(bool onlyZ)
virtual Double_t GetIDMapY(Option_t *opt="")
void init()
default init
KVIDGrid * CalculateDeltaE_EGrid(const KVNameValueList &AperZ, Int_t npoints=30, Double_t xfactor=1.)
void SetGroup(KVGroup *kvg)
KVList fMultiDetExpressions
used to clean up any multi-detector signal expressions generated to calculate X/Y coordinates
void SetLabelFromURI(const Char_t *uri)
Double_t GetIDGridYCoord(KVIDGraph *) const
virtual Bool_t Identify(KVIdentificationResult *, Double_t x=-1., Double_t y=-1.)
virtual Double_t GetIDMapX(Option_t *opt="")
static KVIDTelescope * MakeIDTelescope(const Char_t *name)
virtual Double_t GetPedestalY(Option_t *opt="")
virtual Bool_t SetIdentificationParameters(const KVMultiDetArray *)
void SetIDGrid(KVIDGraph *)
Bool_t HasMassID() const
KVDetector * GetDetector(UInt_t n) const
virtual Bool_t CheckTheoreticalIdentificationThreshold(KVNucleus *, Double_t=0.0)
UInt_t GetGroupNumber()
void addLineToGrid(KVIDGrid *gg, int zz, int aa, int npoints)
KVGroup * GetGroup() const
virtual Double_t GetPedestalX(Option_t *opt="")
virtual void AddDetector(KVDetector *d)
virtual Double_t GetMeanDEFromID(Int_t &status, Int_t Z, Int_t A=-1, Double_t Eres=-1.0)
Bool_t ReadIdentificationParameterFiles(const KVMultiDetArray *multidet)
const KVList * GetDetectors() const
KVDetectorSignal * GetSignalFromGridVar(const KVString &var, const KVString &axe, KVString &det_labels)
KVIDGraph * GetIDGrid()
void Print(Option_t *opt="") const override
std::unordered_map< KVIDGraph *, GraphCoords > fGraphCoords
X/Y coordinates from detector signals for ID maps.
virtual UShort_t GetIDCode()
Bool_t LoadIdentificationParameters(const KVMultiDetArray *multidet, const TString &_filename="")
void CheckIdentificationBilan(const TString &system)
Set status of ID Telescope for given system.
virtual void RemoveIdentificationParameters()
static void OpenIdentificationBilan(const TString &path)
Open IdentificationBilan.dat file with given path.
void SetHasMassID(Bool_t yes=kTRUE)
virtual void Initialize(void)
Double_t GetIDGridXCoord(KVIDGraph *) const
const Char_t * GetDefaultIDGridClass()
void GetIDGridCoords(Double_t &X, Double_t &Y, KVIDGraph *grid, Double_t x=-1, Double_t y=-1)
UInt_t GetSize() const
KVUnownedList fIDGrids
identification grid(s)
Definition: KVIDTelescope.h:93
KVUnownedList fDetectors
list of detectors in telescope
Definition: KVIDTelescope.h:91
KVString GetDetectorLabelsForGridCoord(const KVString &axis) const
const KVList * GetListOfIDGrids() const
std::unique_ptr< KVParticleCondition > fMassIDValidity
may be used to limit mass identification to certain Z and/or A range
KVGroup * fGroup
group to which telescope belongs
Definition: KVIDTelescope.h:92
static TEnv * fgIdentificationBilan
Definition: KVIDTelescope.h:87
virtual void SetIdentificationStatus(KVIdentificationResult *IDR, const KVNucleus *)
virtual void RemoveGrids()
Base class for identification ridge lines corresponding to different nuclear species.
Definition: KVIDZALine.h:33
Full result of one attempted particle identification.
Bool_t IDattempted
=kTRUE if identification was attempted
Bool_t IDOK
general quality of identification, =kTRUE if acceptable identification made
void SetGridName(const Char_t *n)
TString Rejecting_Cut
name of cut in grid which rejected particle for identification
Bool_t Aident
= kTRUE if A of particle established
Int_t A
A of particle found (if Aident==kTRUE)
Int_t Z
Z of particle found (if Zident==kTRUE)
Int_t IDquality
specific quality code returned by identification procedure
void Clear(Option_t *opt="") override
Reset to initial values.
Int_t IDcode
a general identification code for this type of identification
void SetIDType(const Char_t *t)
Bool_t Zident
=kTRUE if Z of particle established
Base class for describing the geometry of a detector array.
Bool_t ReadGridsFromAsciiFile(const Char_t *) const
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
Bool_t IsKnown(int z=-1, int a=-1) const
Definition: KVNucleus.cpp:1272
Int_t GetA() const
Definition: KVNucleus.cpp:792
void SetA(Int_t a)
Definition: KVNucleus.cpp:647
void SetZ(Int_t z, Char_t mt=-1)
Definition: KVNucleus.cpp:696
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:763
Double_t GetLifeTime(Int_t z=-1, Int_t a=-1) const
Definition: KVNucleus.cpp:1033
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:85
void Clear(Option_t *="") override
Empty number list, reset it to initial state.
Bool_t End(void) const
Definition: KVNumberList.h:199
Int_t GetNValues() const
void Begin(void) const
void Add(Int_t)
Add value 'n' to the list.
Int_t Last() const
Returns largest number included in list.
Int_t Next(void) const
Implement selections for KVNucleus objects.
Double_t GetEnergy() const
Definition: KVParticle.h:624
void SetEnergy(Double_t e)
Definition: KVParticle.h:602
virtual TObject * FindObjectByLabel(const Char_t *) const
TObject * First() const override
void Add(TObject *obj) override
TObject * FindObject(const char *name) const override
void Clear(Option_t *option="") override
virtual void SetCleanup(Bool_t enable=kTRUE)
TObject * At(Int_t idx) const override
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
void Begin(TString delim) const
Definition: KVString.cpp:565
Bool_t End() const
Definition: KVString.cpp:634
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:695
Int_t GetNValues(TString delim) const
Definition: KVString.cpp:886
Optimised list in which named objects can only be placed once.
void Add(TObject *obj) override
virtual void SetLineColor(Color_t lcolor)
virtual Double_t GetBinCenter(Int_t bin) const
static TClass * GetClass(Bool_t load=kTRUE, Bool_t silent=kFALSE)
const char * GetVarX() const
const char * GetVarY() const
virtual const char * GetValue(const char *name, const char *dflt) const
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Int_t GetN() const
virtual Double_t Eval(Double_t x, TSpline *spline=nullptr, Option_t *option="") const
virtual Int_t GetNbinsY() const
TAxis * GetXaxis()
virtual Int_t GetNbinsX() const
TAxis * GetYaxis()
virtual Double_t GetBinContent(Int_t bin) const
const char * GetName() const override
virtual void SetName(const char *name)
void AbstractMethod(const char *method) const
void SetBit(UInt_t f)
void ResetBit(UInt_t f)
virtual void Info(const char *method, const char *msgfmt,...) const
Longptr_t ExecPlugin(int nargs)
Ssiz_t Length() const
Int_t Atoi() const
const char * Data() const
void ToUpper()
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
Bool_t IsNull() const
void Form(const char *fmt,...)
TString & Remove(EStripType s, char c)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
TString & ReplaceAll(const char *s1, const char *s2)
TLine * line
Double_t y[n]
Double_t x[n]
const Int_t n
TGraphErrors * gr
Int_t Nint(T x)
Double_t Exp(Double_t x)
constexpr Double_t E()
Double_t Log(Double_t x)
Bool_t Even(Long_t a)
ClassImp(TPyArg)