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) {
163  ok = kFALSE;
164  Warning("Initialize",
165  "ID tel. %s: grid %s has undefined VarX(%s:%p) - WILL NOT USE",
166  GetName(), gr->GetName(), gr->GetVarX(), fGraphCoords[gr].fVarX);
167  }
168  if (!fGraphCoords[gr].fVarY) {
169  ok = kFALSE;
170  Warning("Initialize",
171  "ID tel. %s: grid %s has undefined VarY(%s:%p) - WILL NOT USE",
172  GetName(), gr->GetName(), gr->GetVarY(), fGraphCoords[gr].fVarY);
173  }
174  }
175  // set to true if at least one grid can provide mass identification
176  SetHasMassID(mass_id);
177  // if more than one grid, need to re-order them according to [Dataset].[telescope label].GridOrder
178  if (GetListOfIDGrids()->GetEntries() > 1 && gDataSet) {
179  KVString grid_list = gDataSet->GetDataSetEnv(Form("%s.GridOrder", GetLabel()));
180  ok = kFALSE;
181  if (grid_list == "")
182  Warning("Initialize", "ID telescope %s has %d grids but no %s variable defined",
183  GetName(), GetListOfIDGrids()->GetEntries(), Form("%s.GridOrder", GetLabel()));
184  else if (grid_list.GetNValues(",") != GetListOfIDGrids()->GetEntries())
185  Warning("Initialize", "ID telescope %s has %d grids but %d grids appear in variable %s",
186  GetName(), GetListOfIDGrids()->GetEntries(), grid_list.GetNValues(","), Form("%s.GridOrder", GetLabel()));
187  else {
188  fIDGrids.Clear();
189  grid_list.Begin(",");
190  while (!grid_list.End()) {
191  auto gr_name = grid_list.Next();
192  auto gr_ob = tmp_list.FindObject(gr_name);
193  if (!gr_ob) {
194  Info("Initialize", "IDtel=%s grid %s missing", GetName(), gr_name.Data());
195  }
196  else {
197  fIDGrids.Add(gr_ob);
198  }
199  }
200  ok = kTRUE;
201  }
202  }
203  if (ok) SetBit(kReadyForID);
204  }
205  }
206 
207  if (gDataSet) {
208  SetHasMassID(gDataSet->GetDataSetEnv(Form("%s.MassID", GetLabel()), kFALSE));
209  KVString valid;
210  if ((valid = gDataSet->GetDataSetEnv(Form("%s.MassID.Validity", GetLabel()), "")) != "") {
211  valid.ReplaceAll("Z", "_NUC_->GetZ()");
212  valid.ReplaceAll("A", "_NUC_->GetA()");
213  fMassIDValidity.reset(new KVParticleCondition(valid));
214  }
215  }
216 }
217 
218 
219 
228 
229 void KVIDTelescope::AddDetector(KVDetector* d)
230 {
231  // Add a detector to the telescope.
232  //
233  // Detectors must be added in the order they will be hit by impinging particles,
234  // with the last detector being the one particles stopped in the telescope will stop in.
235  // i.e. dE1, dE2, ..., Eres
236  //
237  // Update name of telescope to "ID_[name of 1st detector]_[name of 2nd detector]_ ... _[name of last detector]"
238 
239  if (d) {
240  fDetectors.Add(d);
241  if (GetSize() > 1) {
242  TString old_name = GetName();
243  old_name += Form("_%s", GetDetectors()->Last()->GetName());
244  SetName(old_name);
245  }
246  else SetName(Form("ID_%s", GetDetector(1)->GetName()));
247  //d->AddIDTelescope(this); <= caused multiple copies to exist in detector's list
248  }
249  else {
250  Warning("AddDetector", "Called with null pointer");
251  }
252 }
253 
254 
255 
259 
261 {
262  // print out telescope structure
263  //if opt="fired" only fired detectors are printed
264 
265  TIter next(GetDetectors());
266  KVDetector* obj;
267 
268  if (!strcmp(opt, "fired")) {
269  while ((obj = (KVDetector*) next())) {
270 
271  if (obj->Fired() || obj->GetEnergy())
272  obj->Print("data");
273  }
274  }
275  else {
276  cout << "\n" << opt << "Structure of KVIDTelescope object: " <<
277  GetName() << " " << GetType() << endl;
278  cout << opt <<
279  "--------------------------------------------------------" <<
280  endl;
281  while ((obj = (KVDetector*) next())) {
282  cout << opt << "Detector: " << obj->GetName() << endl;
283  }
284  }
285 }
286 
287 
288 
289 
292 
293 KVDetector* KVIDTelescope::GetDetector(const Char_t* name) const
294 {
295  // Return a pointer to the detector in the telescope with the name "name".
296 
297  KVDetector* tmp = (KVDetector*) GetDetectors()->FindObject(name);
298  if (!tmp)
299  Warning("GetDetector(const Char_t *name)",
300  "Detector %s not found in telescope %s", name, GetName());
301  return tmp;
302 }
303 
304 
305 
306 
308 
310 {
311  return fGroup;
312 }
313 
314 
315 
316 
318 
320 {
321  fGroup = kvg;
322 }
323 
324 
325 
327 
329 {
330  return (GetGroup() ? GetGroup()->GetNumber() : 0);
331 }
332 
333 
334 
357 
359 {
360  //Default identification method.
361  //
362  //Works for ID telescopes for which one or more identification grids are defined, each
363  //with VARX/VARY parameters corresponding to a KVDetectorSignal or KVDetectorSignalExpression
364  //associated with one or other of the detectors constituting the telescope.
365  //
366  //For identifications using more than one grid, the default behaviour is to try identification
367  //with each grid in turn until a successful identification is obtained. The order in which
368  //the grids should be tried should be specified by a variable with the following format:
369  //
370  //~~~~~~~~~~~~~~~~
371  //[Dataset].[Array Name].[ID type].GridOrder: [Grid1],[Grid2],...
372  //~~~~~~~~~~~~~~~~
373  //
374  //where the name of each grid is given as "VARY_VARX". If no variable defining the order is found,
375  //the grids will be tried in the order they were found in the file containing the grids for this
376  //telescope.
377  //
378  // The KVIdentificationResult is first Clear()ed; then it is filled with IDtype = GetType()
379  // of this identification telescope, IDattempted = true, and the results of the identification
380  // procedure.
381 
382  idr->Clear();
383  idr->IDattempted = true;
384  idr->SetIDType(GetType());
385 
386  KVIDGraph* grid;
387  TIter it(GetListOfIDGrids());
388  while ((grid = (KVIDGraph*)it())) { //loop over grids in order given by [Dataset].[Array Name].[ID type].GridOrder:
389  Double_t de, e;
390  auto fired = GetIDGridCoords(e, de, grid, x, y);// true if all signals required by e & de fired/present in event
391  idr->SetGridName(grid->GetName());
392  if (fired && grid->IsIdentifiable(e, de, &idr->Rejecting_Cut)) {
393  grid->Identify(e, de, idr);
394  if (idr->IDOK) break; // stop on first successful identification
395  }
396  else {
397  // either all signals required to calculate coordinates in grid are not present/fired,
398  // or particle rejected by cut in grid, in which case idr->Rejecting_Cut contains its name.
399  idr->IDOK = kFALSE;
401  }
402  }
403  idr->IDcode = GetIDCode();
404 
405  return kTRUE;
406 }
407 
408 
409 
410 
446 
448 {
449  // Add an identification grid to the list of grids used by this telescope.
450  //
451  // If the grid's VARX and VARY parameters are set and contain the names of valid
452  // detector signals (see formatting rules below) they will be used by
453  // GetIDGridXCoord() and GetIDGridYCoord() to return the coordinates
454  // needed to perform particle identification using the grid.
455  //
456  // The name of the grid is set to "VARY_VARX" (just the signal names, not the detector
457  // label part - see below). This value will be stored in the
458  // KVIdentificationResult corresponding to an attempted identification of a
459  // KVReconstructedNucleus by this grid.
460  //
461  // VARX/VARY Formatting
462  //
463  // To be valid, grid VARX/Y parameters should be set as follows:
464  //
465  //~~~~~~~~~~~~~~~~~~
466  // [signal name]
467  // or [det_label]::[signal name]
468  //~~~~~~~~~~~~~~~~~~
469  //
470  // where
471  //
472  //~~~~~~~~~~~~~~~~~~
473  // [det_label] (optional) = detector label i.e. string returned by KVDetector::GetLabel()
474  // method for detector. By default, VARX is assumed to be the Eres detector
475  // or last detector and VARY the DE detector or first detector
476  // [signal_name] = name of a signal defined for the detector, possibly depending
477  // on availability of calibration
478  //
479  // To see all available signals for a detector, use
480  //
481  // KVDetector::GetListOfDetectorSignals()
482  //~~~~~~~~~~~~~~~~~~
483 
484  if (grid) {
485  fIDGrids.Add(grid);
486  KVString det_labels_x, det_labels_y;
487  KVDetectorSignal* xx = GetSignalFromGridVar(grid->GetVarX(), "X", det_labels_x);
488  KVDetectorSignal* yy = GetSignalFromGridVar(grid->GetVarY(), "Y", det_labels_y);
489  GraphCoords gc;
490  gc.fVarX = xx;
491  gc.fDetLabelsX = det_labels_x;
492  gc.fVarY = yy;
493  gc.fDetLabelsY = det_labels_y;
494  fGraphCoords[grid] = gc;
495  TString grid_name;
496  if (xx && yy) {
497  grid_name.Form("%s_%s", yy->GetName(), xx->GetName());
498  grid->SetName(grid_name);
499  }
500  }
501 }
502 
503 
504 
550 
552 {
553  // Deduce & return pointer to detector signal from grid VARX/VARY parameter
554  //
555  // To be valid, grid VARX/Y parameters should be set using mathematical expressions
556  // which use the following references to detector signals for the telescope:
557  //
558  //~~~~~~~~~~~~~~~~~~
559  // [signal name]
560  // OR [det_label]::[signal name]
561  //~~~~~~~~~~~~~~~~~~
562  //
563  // where
564  //
565  //~~~~~~~~~~~~~~~~~~
566  // [signal_name] = name of a signal defined for 1 of the detectors of the telescope
567  // [det_label] = optional detector label i.e. string returned by
568  // KVDetector::GetLabel() method for detector
569  //~~~~~~~~~~~~~~~~~~
570  //
571  // If `[det_label]` is not given, we assume for `VARX` the last (E) detector,
572  // while for `VARY` we assume the first (dE) detector. If this telescope has only
573  // one detector, we use it for both variables.
574  //
575  // To see all available signals for a detector, use
576  //
577  //~~~~~~~~~~~~~~~~~~
578  // KVDetector::GetListOfDetectorSignals()
579  //~~~~~~~~~~~~~~~~~~
580  //
581  // #### Example ####
582  // Imagine a telescope which combines 2 detectors, with labels SI and CSI (in that order, i.e. SI is the dE detector,
583  // CSI is the residual energy detector). The following cases are valid:
584  //
585  //~~~~~
586  // VARX = Energy : use 'Energy' signal of CSI detector (default for VARX)
587  // VARY = ADC : use 'ADC' signal of SI detector (default for VARY)
588  // VARY = CSI::Energy : use 'Energy' signal of CSI detector (overrides default for VARY, SI)
589  // VARX = (Q3-Q3Fast)/(0.8*Q3) : uses 'Q3' and 'Q3Fast' signals of CSI detector (default for VARX)
590  // VARY = 1.5*SI::ADC - CSI::Q3Fast/CSI::Q3 : combination of signals from both detectors
591  //~~~~~
592  //
593  // The 'det_labels' string will be filled with a comma-separated list of the labels of each detector used
594  // in the expressions.
595  //
596  // \note if any signals are not defined, they will be evaluated as zero
597 
598  if (var == "") {
599  Warning("GetSignalFromGridVar",
600  "No VAR%s defined for grid for telescope %s. KVIDTelescope-derived class handling identification must override GetIDMapX/GetIDMapY",
601  axe.Data(), GetName());
602  return nullptr;
603  }
604 
605  KVString dum = var;
606  KVDetector* det = nullptr;
607  KVDetectorSignal* ds(nullptr);
608  KVString sig_type;
609 
610  // check if VARX/Y is a mathematical expression
611  Bool_t is_expression = var.GetNValues("+-*/()") > 1;
612  // examine each term in expression (this will split the elements in a mathematical expression,
613  // or just take the whole expression if no math operators are present)
614  var.Begin("+-*/()");
615  bool explicit_det_ref = false;
616  bool multidetexpr = false;
617  KVString det_label;
618  while (!var.End()) {
619  KVString t = var.Next();
620  // check if we have an explicit reference to a detector
621  if (t.Contains("::")) {
622  t.Begin("::");
623  det_label = t.Next();
624  auto _det = (KVDetector*)GetDetectors()->FindObjectByLabel(det_label);
625  if (_det) {
626  explicit_det_ref = true;
627  sig_type = t.Next();
628  // check signal is defined for detector
629  if (_det->GetDetectorSignal(sig_type)) {
630  if (det_labels.Length()) {
631  if (!det_labels.Contains(det_label)) det_labels += Form(",%s", det_label.Data());
632  }
633  else
634  det_labels = det_label;
635  }
636  else {
637  Warning("GetSignalFromGridVar", "signal '%s' not found in det '%s' -> replaced by '0'", sig_type.Data(), _det->GetName());
638  dum.ReplaceAll(Form("%s::%s", _det->GetLabel(), sig_type.Data()), "0");
639  }
640  }
641  if (det && (_det != det)) multidetexpr = true; // more than 1 detector is referenced
642  det = _det;
643  }
644  }
645  // treat all cases with explicit detector references
646  if (explicit_det_ref) {
647  if (!multidetexpr) {
648  // only 1 detector is directly referenced
649  if (!is_expression) {
650  // it is not a mathematical expression: this is just the name of a detector signal
651  return det->GetDetectorSignal(sig_type);
652  }
653  else {
654  // math expression with explicit reference to signal(s) of 1 detector
655  sig_type = var;
656  // remove all explicit references to detector: 'det' pointer has been set
657  sig_type.ReplaceAll(Form("%s::", det_label.Data()), "");
658  // expression will be handled below
659  }
660  }
661  else {
662  // use specific constructor for explicit detector references
663  ds = new KVDetectorSignalExpression(var, dum, GetDetectors());
664  if (!ds->IsValid()) {
665  delete ds;
666  ds = nullptr;
667  Error("GetSignalFromGridVar",
668  "Problem initialising ID-grid %s coordinate for telescope %s."
669  " Check definition of VAR%s for grid (=%s)",
670  axe.Data(), GetName(), axe.Data(), var.Data());
671  return ds;
672  }
674  return ds;
675  }
676  }
677  else {
678  // no explicit reference to detectors - by default, use detector (1) for Y axis,
679  // and the last detector for X axis
680  if (axe == "Y" || GetSize() == 1) det = GetDetector(1);
681  else det = GetDetector(GetSize());
682  sig_type = var;
683  det_labels = det->GetLabel();
684  }
685  ds = det->GetDetectorSignal(sig_type);
686  if (!ds && is_expression) {
687  // sig_type is not the name of a known signal: it must be an expression using known signal names
688  if (!det->AddDetectorSignalExpression(sig_type, sig_type)) {
689  Error("GetSignalFromGridVar",
690  "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)",
691  axe.Data(), GetName(), sig_type.Data(), det->GetName(), axe.Data(), var.Data());
692  ds = nullptr;
693  }
694  else
695  ds = det->GetDetectorSignal(sig_type);
696  }
697  return ds;
698 }
699 
700 
701 
703 
705 {
707  if (!cl || !cl->InheritsFrom("KVIDZAGrid")) cl = TClass::GetClass("KVIDZAGrid");
708  KVIDGrid* idgrid = (KVIDGrid*)cl->New();
709 
710  idgrid->AddIDTelescope(this);
711  idgrid->SetOnlyZId(onlyZ);
712  idgrid->SetRunList("1-10000");
713 
714  KVIDCutLine* B_line = (KVIDCutLine*)idgrid->Add("OK", "KVIDCutLine");
715  Int_t npoi_bragg = 0;
716  B_line->SetName("Bragg_line");
717  B_line->SetAcceptedDirection("right");
718 
719  return idgrid;
720 }
721 
722 
723 
733 
734 void KVIDTelescope::addLineToGrid(KVIDGrid* idgrid, int zz, int aa, int npoints, std::optional<Double_t> einc_max)
735 {
736 
737  //loop over energy
738  //first find :
739  // ****E1 = energy at which particle passes 1st detector and starts to enter in the 2nd one****
740  // E2 = energy at which particle passes the 2nd detector
741  //then perform npoints calculations between these two energies and use these
742  //to construct a KVIDZALine
743  //
744  // if einc_max is given, limit incident energy of particles
745 
746  double xfactor = 1.;
747 
748  KVNucleus part;
749 
750  KVDetector* det_de = GetDetector(1);
751  KVDetector* det_eres = GetDetector(2);
752 
753  Double_t SeuilE = 0.1;
754 
755  part.SetZ(zz);
756  part.SetA(aa);
757 
758  Double_t E1, E2;
759  //find E1
760  //go from SeuilE MeV to det_de->GetEIncOfMaxDeltaE(part.GetZ(),part.GetA()))
761  Double_t E1min = SeuilE, E1max = det_de->GetEIncOfMaxDeltaE(zz, aa);
762  E1 = (E1min + E1max) / 2.;
763 
764  while ((E1max - E1min) > SeuilE) {
765 
766  part.SetEnergy(E1);
767  det_de->Clear();
768  det_eres->Clear();
769 
770  det_de->DetectParticle(&part);
771  det_eres->DetectParticle(&part);
772  if (det_eres->GetEnergy() > SeuilE) {
773  //particle got through - decrease energy
774  E1max = E1;
775  E1 = (E1max + E1min) / 2.;
776  }
777  else {
778  //particle stopped - increase energy
779  E1min = E1;
780  E1 = (E1max + E1min) / 2.;
781  }
782  }
783 
784  //add point to Bragg line
785  Double_t dE_B = det_de->GetMaxDeltaE(zz, aa);
786  Double_t E_B = det_de->GetEIncOfMaxDeltaE(zz, aa);
787  Double_t Eres_B = det_de->GetERes(zz, aa, E_B);
788 
789  KVIDCutLine* B_line = (KVIDCutLine*)idgrid->GetCut("Bragg_line");
790  if (B_line) B_line->SetPoint(B_line->GetN(), Eres_B, dE_B);
791 
792  //find E2
793  //go from E1 MeV to maximum value where the energy loss formula is valid
794  Double_t E2min = E1, E2max = det_eres->GetEmaxValid(part.GetZ(), part.GetA());
795  E2 = (E2min + E2max) / 2.;
796 
797  while ((E2max - E2min > SeuilE)) {
798 
799  part.SetEnergy(E2);
800  det_de->Clear();
801  det_eres->Clear();
802 
803  det_de->DetectParticle(&part);
804  det_eres->DetectParticle(&part);
805  if (part.GetEnergy() > SeuilE) {
806  //particle got through - decrease energy
807  E2max = E2;
808  E2 = (E2max + E2min) / 2.;
809  }
810  else {
811  //particle stopped - increase energy
812  E2min = E2;
813  E2 = (E2max + E2min) / 2.;
814  }
815  }
816  E2 *= xfactor;
817  if ((!strcmp(det_eres->GetType(), "CSI")) && (E2 > 5000)) E2 = 5000;
818  if(einc_max) E2 = std::min(E2, einc_max.value());
819  // printf("z=%d a=%d E1=%lf E2=%lf\n",zz,aa,E1,E2);
820  KVIDZALine* line = (KVIDZALine*)idgrid->Add("ID", "KVIDZALine");
821  if (TMath::Even(zz)) line->SetLineColor(4);
822  line->SetZ(zz);
823  line->SetA(aa);
824 
825  Double_t logE1 = TMath::Log(E1);
826  Double_t logE2 = TMath::Log(E2);
827  Double_t dLog = (logE2 - logE1) / (npoints - 1.);
828 
829  for (Int_t i = 0; i < npoints; i++) {
830  // Double_t E = E1 + i*(E2-E1)/(npoints-1.);
831  Double_t E = TMath::Exp(logE1 + i * dLog);
832 
833  Double_t Eres = 0.;
834  Int_t niter = 0;
835  while (Eres < SeuilE && niter <= 20) {
836  det_de->Clear();
837  det_eres->Clear();
838 
839  part.SetEnergy(E);
840 
841  det_de->DetectParticle(&part);
842  det_eres->DetectParticle(&part);
843 
844  Eres = det_eres->GetEnergy();
845  E += SeuilE;
846  niter += 1;
847  }
848  if (!(niter > 20)) {
849  Double_t dE = det_de->GetEnergy();
850  Double_t gEres, gdE;
851  line->GetPoint(i - 1, gEres, gdE);
852  line->SetPoint(i, Eres, dE);
853 
854  }
855  }
856  //printf("sort de boucle points");
857 
858 
859 }
860 
861 
862 
871 
873 {
874  // Returns a comma-separated list of the labels of the detectors used to determine
875  // the "x" or "y" coordinates of the identification grid(s)
876  //
877  // If there is more than 1 grid and the list is not the same for all grids,
878  // prints a warning message.
879  //
880  // \param[in] axis name of grid axis i.e. "x", "X", "y" or "Y" (case insensitive)
881 
882  KVString _axis = axis;
883  _axis.ToUpper();
884 
885  if (_axis != "X" && _axis != "Y") {
886  Error("GetDetectorLabelsForGridCoord", "Called with illegal axis name '%s'", axis.Data());
887  return "";
888  }
889  KVString lab_list;
890 
891  for (auto& __gc : fGraphCoords) {
892  KVString labs;
893  if (_axis == "X") labs = __gc.second.fDetLabelsX;
894  else labs = __gc.second.fDetLabelsY;
895  if (lab_list.Length() && lab_list != labs) {
896  Error("GetDetectorLabelsForGridCoord", "Grids for telescope %s use different detector types for %s-coordinate: "
897  "%s and %s", GetName(), _axis.Data(), lab_list.Data(), labs.Data());
898  return "";
899  }
900  lab_list = labs;
901  }
902  return lab_list;
903 }
904 
905 
906 
907 
911 
913 {
914  //Return the first in the list of identification grids used by this telescope
915  //(this is for backwards compatibility with ID telescopes which had only one grid).
916  return (KVIDGraph*)GetListOfIDGrids()->First();
917 }
918 
919 
920 
921 
924 
926 {
927  //Return pointer to grid using position in list. First grid has index = 1.
928  if (index < 1) {
929  Error("GetIDGrid(int)", "Index must be >=1!");
930  return nullptr;
931  }
932  return (KVIDGraph*)GetListOfIDGrids()->At(index - 1);
933 }
934 
935 
936 
937 
941 
943 {
944  //Return pointer to grid using "label" to search in list of grids associated
945  //to this telescope.
946  return (KVIDGraph*)GetListOfIDGrids()->FindObjectByLabel(label);
947 }
948 
949 
950 
951 
953 
955 {
956  AbstractMethod("GetIDMapX");
957  return -1.;
958 }
959 
960 
961 
966 
968 {
969  // Returns the pedestal associated with the 2nd detector of the telescope,
970  // optionally depending on the given option string.
971  // By default this returns 0, and should be overridden in specific implementations.
972 
973  return 0.;
974 }
975 
976 
977 
982 
984 {
985  // Returns the pedestal associated with the 1st detector of the telescope,
986  // optionally depending on the given option string.
987  // By default this returns 0, and should be overridden in specific implementations.
988 
989  return 0.;
990 }
991 
992 
993 
1000 
1002 {
1003  // \return value of X coordinate to be used with the given ID grid
1004  //
1005  // This corresponds to whatever was given as parameter "VARX" for the grid
1006  //
1007  // \param[out] fired set to true or false depending on whether detector signal(s) used for X were fired (read from last raw event)
1008 
1009  KVDetectorSignal* ds = fGraphCoords[g].fVarX;
1010  if (ds) {
1011  fired = ds->IsFired();
1012  return ds->GetValue();
1013  }
1014  fired = false;
1015  return -1;
1016 }
1017 
1018 
1019 
1026 
1028 {
1029  // \return value of Y coordinate to be used with the given ID grid
1030  //
1031  // This corresponds to whatever was given as parameter "VARY" for the grid
1032  //
1033  // \param[out] fired set to true or false depending on whether detector signal(s) used for Y were fired (read from last raw event)
1034 
1035  KVDetectorSignal* ds = fGraphCoords[g].fVarY;
1036  if (ds) {
1037  fired = ds->IsFired();
1038  return ds->GetValue();
1039  }
1040  fired = false;
1041  return -1;
1042 }
1043 
1044 
1045 
1046 
1048 
1050 {
1051  AbstractMethod("GetIDMapY");
1052  return -1.;
1053 }
1054 
1055 
1056 
1057 
1061 
1063 {
1064  //Remove all identification grids for this ID telescope
1065  //Grids are not deleted as this is handled by gIDGridManager
1066  fIDGrids.Clear();
1067  fGraphCoords.clear();
1068 }
1069 
1070 
1071 
1072 
1091 
1093 {
1094  //Static function which will create an instance of the KVIDTelescope-derived class
1095  //corresponding to 'name'
1096  //These are defined as 'Plugin' objects in the file $KVROOT/KVFiles/.kvrootrc :
1097  //~~~~~~~
1098  // # The KVMultiDetArray::GetIDTelescopes(KVDetector*de, KVDetector*e) method uses these plugins to
1099  // # create KVIDTelescope instances adapted to the specific array geometry and detector types.
1100  // # For each pair of detectors we look for a plugin with one of the following names:
1101  // # [name_of_dataset].de_detector_type[de detector thickness]-e_detector_type[de detector thickness]
1102  // # Each characteristic in [] brackets may or may not be present in the name; first we test for names
1103  // # with these characteristics, then all combinations where one or other of the characteristics is not present.
1104  // # In addition, we first test all combinations which begin with [name_of_dataset].
1105  // # The first plugin found in this order will be used.
1106  // # In addition, if for one of the two detectors there is a plugin called
1107  // # [name_of_dataset].de_detector_type[de detector thickness]
1108  // # [name_of_dataset].e_detector_type[e detector thickness]
1109  // # then we add also an instance of this 1-detector identification telescope.
1110  //~~~~~~~
1111 
1112  TPluginHandler* ph;
1113  //check and load plugin library
1114  if (!(ph = LoadPlugin("KVIDTelescope", uri)))
1115  return 0;
1116 
1117  //execute constructor/macro for identification telescope
1118  KVIDTelescope* mda = (KVIDTelescope*) ph->ExecPlugin(0);
1119  if (mda) {
1120  //set label of telescope with URI used to find plugin (minus dataset name)
1121  mda->SetLabelFromURI(uri);
1122  }
1123 
1124  return mda;
1125 }
1126 
1127 
1128 
1129 
1133 
1135 {
1136  //PRIVATE METHOD
1137  //Sets label of telescope based on URI of plugin describing child class for this telescope
1138 
1139  TString _uri(uri);
1140  if (gDataSet && _uri.BeginsWith(gDataSet->GetName())) _uri.Remove(0, strlen(gDataSet->GetName()) + 1);
1141  SetLabel(_uri.Data());
1142 }
1143 
1144 
1145 
1146 
1179 
1181 {
1182  // Initialise the identification parameters (grids, etc.) of ALL identification telescopes of this
1183  // kind (label) in the multidetector array. Therefore this method need only be called once, and not
1184  // called for each telescope. The kind/label (returned by GetLabel) of the telescope corresponds
1185  // to the URI used to find the plugin class in $KVROOT/KVFiles/.kvrootrc.
1186  //
1187  // By default, this method first looks for a file containing the list of all files containing the
1188  // identification grids. The default name of this file is:
1189  //
1190  //~~~
1191  // [label].IdentificationFiles.dat
1192  //~~~
1193  //
1194  // or any name given by the environment variable
1195  //
1196  //~~~
1197  // [dataset].IdentificationParameterList.[label]: [filename]
1198  //~~~
1199  //
1200  // If neither file is found, we look for a file with name:
1201  //
1202  //~~~
1203  // IDGrids.[label].dat
1204  //~~~
1205  //
1206  // or any name given by the environment variable
1207  //
1208  //~~~
1209  // [dataset].IdentificationParameterFile.[label]: [filename]
1210  //~~~
1211  //
1212  // which is assumed to contain all identification grids for this type of telescope.
1213 
1215  if(!ok) ok = LoadIdentificationParameters(multidet);
1216 
1217  if (!ok)
1218  Warning("SetIdentificationParameters", "No files found. See doc for this method.");
1219 
1220  return ok;
1221 }
1222 
1223 
1224 
1225 
1243 
1245 {
1246  // In the case where the identification grids are stored in several files, this method looks
1247  // for a file with name
1248  //
1249  //~~~
1250  // [label].IdentificationFiles.dat
1251  //~~~
1252  //
1253  // or any name given by the environment variable
1254  //
1255  //~~~
1256  // [dataset].IdentificationParameterList.[label]: [filename]
1257  //~~~
1258  //
1259  // and if found reads from it the list of files containing the identification grids.
1260  //
1261  // \returns kTRUE if file found and read, kFALSE if no file found
1262 
1263  TString filename = gDataSet->GetDataSetEnv(Form("IdentificationParameterList.%s", GetLabel()),
1264  Form("%s.IdentificationFiles.dat", GetLabel()));
1265 
1266  if(gDataSet->FindDataSetFile(filename))
1267  {
1268  KVFileReader fr;
1269  auto path = gDataSet->GetFullPathToDataSetFile(filename);
1270  Info("ReadIdentificationParameterFiles", "Using file %s", path.Data());
1271  fr.OpenFileToRead(path);
1272  Bool_t ok = kTRUE;
1273  while (fr.IsOK()) {
1274  fr.ReadLine(0);
1275 
1276  if (fr.GetCurrentLine() != "")
1277  ok &= LoadIdentificationParameters(multidet, fr.GetCurrentLine().Data());
1278  }
1279 
1280  fr.CloseFile();
1281  return ok;
1282  }
1283 
1284  return kFALSE;
1285 }
1286 
1287 
1288 
1289 
1307 
1309 {
1310  // \param[in] _filename optional name of dataset file containing grids for this telescope type
1311  // \returns kTRUE if a file is found and read in
1312  //
1313  // If no value for _filename is given, this method looks for a dataset file named
1314  //
1315  //~~~
1316  // IDGrids.[label].dat
1317  //~~~
1318  //
1319  // or any name given by the environment variable
1320  //
1321  //~~~
1322  // [dataset].IdentificationParameterFile.[label]: [filename]
1323  //~~~
1324  //
1325  // which is assumed to contain all identification grids for this type of telescope.
1326 
1327  TString filename;
1328  if(_filename.IsNull())
1329  {
1330  filename = gDataSet->GetDataSetEnv(Form("IdentificationParameterFile.%s", GetLabel()),
1331  Form("IDGrids.%s.dat",GetLabel()));
1332  }
1333  else
1334  filename = _filename;
1335 
1336  if(gDataSet->FindDataSetFile(filename))
1337  {
1338  auto path = gDataSet->GetFullPathToDataSetFile(filename);
1339  //Read grids from file
1340  Info("LoadIdentificationParameters", "Using file %s", path.Data());
1341  multidet->ReadGridsFromAsciiFile(path);
1342  return kTRUE;
1343  }
1344  else if (!_filename.IsNull())
1345  {
1346  // if the method was called with a definite filename to look for, we assume that this
1347  // file is supposed to exist, therefore if it is not found, there is a problem!
1348  Error("LoadIdentificationParameters",
1349  "File %s not found. Should be in %s",
1350  _filename.Data(), gDataSet->GetDataSetDir());
1351  return kFALSE;
1352  }
1353  return kFALSE;
1354 }
1355 
1356 
1357 
1358 
1370 
1372 {
1373  //Remove identification parameters from telescope in such a way that they
1374  //can subsequently be reset e.g. with a new version.
1375  //This is used by KVMultiDetArray::UpdateIdentifications.
1376  //Child classes with specific SetIdentificationParameters methods should
1377  //also redefine this method in order to remove (destroy) cleanly the objects
1378  //created in SetIdentificationParameters.
1379  //
1380  //This default method takes the list of grids associated to the telescope,
1381  //and for each one: 1) checks if it is still in the gIDGridManager's list
1382  //2) if yes, delete the grid and remove it from gIDGridManager
1383 
1384  TIter next_grid(GetListOfIDGrids());
1385  KVIDGrid* grid;
1386  while ((grid = (KVIDGrid*)next_grid())) {
1387 
1388  if (gIDGridManager->GetGrids()->FindObject(grid)) { //this only works if KVIDTelescope uses TObject:IsEqual method (i.e. compares pointers)
1389 
1390  gIDGridManager->DeleteGrid(grid);
1391 
1392  }
1393  }
1394  //clear list of grids
1395  fIDGrids.Clear();
1397 }
1398 
1399 
1400 
1401 //void KVIDTelescope::CalculateParticleEnergy(KVReconstructedNucleus* nuc)
1402 //{
1403 // // The energy of each particle is calculated as follows:
1404 // //
1405 // // E = dE_1 + dE_2 + ... + dE_N
1406 // //
1407 // // dE_1, dE_2, ... = energy losses measured in each detector through which
1408 // // the particle has passed (or stopped, in the case of dE_N).
1409 // // These energy losses are corrected for (Z,A)-dependent effects
1410 // // such as pulse-heigth defect in silicon detectors, losses in
1411 // // windows of gas detectors, etc.
1412 // //
1413 // // Whenever possible, the energy loss for fired detectors which are uncalibrated
1414 // // or not functioning is calculated. In this case the status returned by GetCalibStatus()
1415 // // will be KVIDTelescope::kCalibStatus_Calculated.
1416 // // If none of the detectors is calibrated, the particle's energy cannot be calculated &
1417 // // the status will be KVIDTelescope::kCalibStatus_NoCalibrations.
1418 // // Otherwise, the status code will be KVIDTelescope::kCalibStatus_OK.
1419 
1420 // //status code
1421 // fCalibStatus = kCalibStatus_NoCalibrations;
1422 
1423 // UInt_t z = nuc->GetZ();
1424 // //uncharged particles
1425 // if (z == 0) return;
1426 
1427 // KVDetector* d1 = GetDetector(1);
1428 // KVDetector* d2 = (GetSize() > 1 ? GetDetector(2) : 0);
1429 // Bool_t d1_cal = d1->IsCalibrated();
1430 // Bool_t d2_cal = (d2 ? d2->IsCalibrated() : kFALSE);
1431 
1432 // //no calibrations
1433 // if (!d1_cal && !d2)
1434 // return;
1435 // if ((d1 && d2) && !d1_cal && !d2_cal)
1436 // return;
1437 
1438 // //status code
1439 // fCalibStatus = kCalibStatus_OK;
1440 
1441 // UInt_t a = nuc->GetA();
1442 
1443 // // particles stopped in first member of telescope
1444 // if (nuc->GetStatus() == 3) {
1445 // if (d1_cal) {
1446 // nuc->SetEnergy(d1->GetCorrectedEnergy(nuc, -1, kFALSE)); //N.B.: transmission=kFALSE because particle stop in d1
1447 // }
1448 // return;
1449 // }
1450 
1451 // Double_t e1, e2, einc;
1452 // e1 = e2 = einc = 0.0;
1453 
1454 // if (!d1_cal) {//1st detector not calibrated - calculate from residual energy in 2nd detector
1455 
1456 // //second detector must exist and have all acquisition parameters fired with above-pedestal value
1457 // if (d2 && d2->Fired("Pall")) e2 = d2->GetCorrectedEnergy(nuc, -1, kFALSE); //N.B.: transmission=kFALSE because particle stop in d2
1458 // if (e2 <= 0.0) {
1459 // // zero energy loss in 2nd detector ? can't do anything...
1460 // fCalibStatus = kCalibStatus_NoCalibrations;
1461 // return;
1462 // }
1463 // //calculate & set energy loss in dE detector
1464 // //N.B. using e2 for the residual energy after detector 1 means
1465 // //that we are assuming the particle stops in detector 2.
1466 // //if this is not true, we will underestimate the energy of the particle.
1467 // e1 = d1->GetDeltaEFromERes(z, a, e2);
1468 // if (e1 < 0.0) e1 = 0.0;
1469 // else {
1470 // d1->SetEnergyLoss(e1);
1471 // d1->SetEResAfterDetector(e2);
1472 // e1 = d1->GetCorrectedEnergy(nuc);
1473 // //status code
1474 // fCalibStatus = kCalibStatus_Calculated;
1475 // }
1476 // }
1477 // else { //1st detector is calibrated too: get corrected energy loss
1478 
1479 // e1 = d1->GetCorrectedEnergy(nuc);
1480 
1481 // }
1482 
1483 // if (d2 && !d2_cal) {//2nd detector not calibrated - calculate from energy loss in 1st detector
1484 
1485 // e1 = d1->GetCorrectedEnergy(nuc);
1486 // if (e1 <= 0.0) {
1487 // // zero energy loss in 1st detector ? can't do anything...
1488 // fCalibStatus = kCalibStatus_NoCalibrations;
1489 // return;
1490 // }
1491 // //calculate & set energy loss in 2nd detector
1492 // e2 = d1->GetEResFromDeltaE(z, a);
1493 // if (e2 < 0.0) e2 = 0.0;
1494 // else {
1495 // e2 = d2->GetDeltaE(z, a, e2);
1496 // d2->SetEnergyLoss(e2);
1497 // e2 = d2->GetCorrectedEnergy(nuc);
1498 // //status code
1499 // fCalibStatus = kCalibStatus_Calculated;
1500 // }
1501 // }
1502 // else if (d2) { //2nd detector is calibrated too: get corrected energy loss
1503 
1504 // e2 = d2->GetCorrectedEnergy(nuc, -1, kFALSE);//N.B.: transmission=kFALSE because particle assumed to stop in d2
1505 // // recalculate corrected energy in first stage using info on Eres
1506 // d1->SetEResAfterDetector(e2);
1507 // e1 = d1->GetCorrectedEnergy(nuc);
1508 // }
1509 
1510 // //incident energy of particle (before 1st member of telescope)
1511 // einc = e1 + e2;
1512 
1513 // Double_t coherence_tolerance = gEnv->GetValue("KVIDTelescope.CoherencyTolerance", 1.05);
1514 // if (coherence_tolerance < 1) coherence_tolerance += 1.00;
1515 
1516 // //Now we have to work our way up the list of detectors from which the particle was
1517 // //reconstructed. For each fired & calibrated detector which is only associated with
1518 // //one particle in the events, we add the corrected measured energy loss
1519 // //to the particle. For uncalibrated, unfired detectors and detectors through which
1520 // //more than one particle has passed, we calculate the corrected energy loss and add it
1521 // //to the particle.
1522 // int ndets = nuc->GetNumDet();
1523 // if (ndets > (int)GetSize()) { //particle passed through other detectors before this idtelesocpe
1524 // //look at detectors not in this id telescope
1525 // int idet = GetSize();//next detector after delta-e member of IDTelescope (stopping detector = 0)
1526 // while (idet < ndets) {
1527 
1528 // KVDetector* det = nuc->GetDetector(idet);
1529 // if (det->Fired() && det->IsCalibrated() && det->GetNHits() == 1) {
1530 // Double_t dE = det->GetEnergy();
1531 // //in order to check if particle was really the only one to
1532 // //hit each detector, we calculate the particle's energy loss
1533 // //from its residual energy. if the measured energy loss is
1534 // //significantly larger, there may be a second particle.
1535 // e1 = det->GetDeltaEFromERes(z, a, einc);
1536 // if (e1 < 0.0) e1 = 0.0;
1537 // det->SetEResAfterDetector(einc);
1538 // dE = det->GetCorrectedEnergy(nuc);
1539 // einc += dE;
1540 // }
1541 // else {
1542 // // Uncalibrated/unfired/multihit detector. Calculate energy loss.
1543 // //calculate energy of particle before detector from energy after detector
1544 // e1 = det->GetDeltaEFromERes(z, a, einc);
1545 // if (e1 < 0.0) e1 = 0.0;
1546 // if (det->GetNHits() > 1) {
1547 // //Info("CalculateParticleEnergy",
1548 // // "Detector %s was hit by %d particles. Calculated energy loss for particle %f MeV",
1549 // // det->GetName(), det->GetNHits(), e1);
1550 // if (!(det->Fired() && det->IsCalibrated())) {
1551 // det->SetEnergyLoss(e1 + det->GetEnergy());// sum up calculated energy losses in uncalibrated detector
1552 // }
1553 // //status code
1554 // fCalibStatus = kCalibStatus_Multihit;
1555 // }
1556 // else if (!det->Fired() || !det->IsCalibrated()) {
1557 // //Info("CalculateParticleEnergy",
1558 // // "Detector %s uncalibrated/not fired. Calculated energy loss for particle %f MeV",
1559 // // det->GetName(), e1);
1560 // det->SetEnergyLoss(e1);
1561 // //status code
1562 // fCalibStatus = kCalibStatus_Calculated;
1563 // }
1564 // det->SetEResAfterDetector(einc);
1565 // e1 = det->GetCorrectedEnergy(nuc, e1);
1566 // einc += e1;
1567 // }
1568 // idet++;
1569 // }
1570 // }
1571 // //einc is now the energy of the particle before crossing the first detector
1572 // nuc->SetEnergy(einc);
1573 //}
1574 
1575 
1576 
1586 
1588 {
1589  // Returns name of default ID grid class for this ID telescope.
1590  // This is defined in a .kvrootrc configuration file by one of the following:
1591  // KVIDTelescope.DefaultGrid:
1592  // KVIDTelescope.DefaultGrid.[type]:
1593  // where [type] is the type of this identification telescope (which is given
1594  // by the character string returned by method GetLabel()... sorry :( )
1595  // If no default grid is defined for the specific type of this telescope,
1596  // the default defined by KVIDTelescope.DefaultGrid is used.
1597 
1598  TString specific;
1599  specific.Form("KVIDTelescope.DefaultGrid.%s", GetLabel());
1600  return gEnv->GetValue(specific.Data(), gEnv->GetValue("KVIDTelescope.DefaultGrid", "KVIDGraph"));
1601 }
1602 
1603 
1604 
1610 
1612 {
1613  //Create a dE-E grid (energy loss in detector 1 versus residual energy in detector 2) for a given list of isotopes
1614  // - AperZ : list of A for each Z. (ex: "1=1-3,2=3-4 5,3=6-8,4=7 9-12"...)
1615  // - npoints : number of points in each generated line
1616  // - xfactor : scales the detector 2 thickness to prolongate the lines
1617 
1618  if (GetSize() <= 1) return 0;
1619  if (!GetDetector(1) || !GetDetector(2)) return 0;
1620  double thickness = GetDetector(2)->GetThickness();
1621  GetDetector(2)->SetThickness(thickness * xfactor);
1622 
1623  Info("CalculateDeltaE_EGrid", "called with KVNameValueList");
1624 
1625  KVIDGrid* idgrid = newGrid(0);
1626 
1627  for (auto& par : AperZ) {
1628  KVString tmp = par.GetName();
1629  int zz = tmp.Atoi();
1630  KVNumberList alist = par.GetString();
1631  for (auto aa : alist) {
1632  addLineToGrid(idgrid, zz, aa, npoints);
1633  }
1634  }
1635 
1636  GetDetector(2)->SetThickness(thickness);
1637 
1638  return idgrid;
1639 
1640 }
1641 
1642 
1643 
1651 
1652 KVIDGrid* KVIDTelescope::CalculateDeltaE_EGrid(const KVNumberList& Zrange, Int_t deltaA, Int_t npoints, Double_t lifetime, UChar_t massformula, Double_t xfactor)
1653 {
1654  //Create a dE-E grid (energy loss in detector 1 versus residual energy in detector 2) for a given list of isotopes
1655  // - Zrange : list of element for which we create lines
1656  // - deltaA : number of isotopes generated for each Z around massformula (ex: deltaA=1, Aref-1 Aref Aref+1)
1657  // - npoints : number of points in each generated line
1658  // - lifetime: remove isotopes with lifetime lower than this value
1659  // - xfactor : scales the detector 2 thickness to prolongate the lines
1660 
1661  if (GetSize() <= 1) return 0;
1662  if (!GetDetector(1) || !GetDetector(2)) return 0;
1663  double thickness = GetDetector(2)->GetThickness();
1664  GetDetector(2)->SetThickness(thickness * xfactor);
1665 
1666  Info("CalculateDeltaE_EGrid", "called with KVNumberList");
1667 
1668  KVIDGrid* idgrid = newGrid(0);
1669  KVNucleus part;
1670 
1671  Zrange.Begin();
1672  while (!Zrange.End()) {
1673  Int_t zz = Zrange.Next();
1674  part.SetZ(zz, massformula);
1675  Int_t aref = part.GetA();
1676  for (Int_t aa = aref - deltaA; aa <= aref + deltaA; aa += 1) {
1677  part.SetA(aa);
1678  if (part.IsKnown() && (part.GetLifeTime() > lifetime)) {
1679  addLineToGrid(idgrid, zz, aa, npoints);
1680 
1681  }
1682  }
1683  }
1684  GetDetector(2)->SetThickness(thickness);
1685  return idgrid;
1686 }
1687 
1688 
1689 
1692 
1693 KVIDGrid* KVIDTelescope::CalculateDeltaE_EGrid(const KVNumberList& Zrange, Int_t npoints, Double_t Einc_max, UChar_t massformula)
1694 {
1695  // Generate ID grid for given list of Z values (1 line per Z), up to the given maximum incident energy
1696  if (GetSize() <= 1) return 0;
1697  if (!GetDetector(1) || !GetDetector(2)) return 0;
1698 
1699  Info("CalculateDeltaE_EGrid", "called with Einc_max=%f", Einc_max);
1700 
1701  KVIDGrid* idgrid = newGrid(0);
1702  KVNucleus part;
1703 
1704  Zrange.Begin();
1705  while (!Zrange.End()) {
1706  Int_t zz = Zrange.Next();
1707  part.SetZ(zz, massformula);
1708  addLineToGrid(idgrid, zz, part.GetA(), npoints, Einc_max);
1709  }
1710  return idgrid;
1711 }
1712 
1713 
1714 
1720 
1722 {
1723  //Create a dE-E grid (energy loss in detector 1 versus residual energy in detector 2) for a given list of isotopes
1724  // - haa_zz : lines will be generated for A,Z filled with 1 in this histogram
1725  // - Zonly : if true, generate only one line per Z with the <A>(Z) of the histogram
1726  // - npoints : number of points in each generated line
1727 
1728  if (GetSize() <= 1) return 0;
1729  if (!GetDetector(1) || !GetDetector(2)) return 0;
1730  double thickness = GetDetector(2)->GetThickness();
1731 
1732  Info("CalculateDeltaE_EGrid", "called with TH2");
1733 
1734  KVIDGrid* idgrid = newGrid(0);
1735  KVNucleus part;
1736 
1737  for (Int_t nx = 1; nx <= haa_zz->GetNbinsX(); nx += 1) {
1738 
1739  Int_t zz = TMath::Nint(haa_zz->GetXaxis()->GetBinCenter(nx));
1740  KVNumberList nlA;
1741  Double_t sumA = 0, sum = 0;
1742  for (Int_t ny = 1; ny <= haa_zz->GetNbinsY(); ny += 1) {
1743  Double_t stat = haa_zz->GetBinContent(nx, ny);
1744  if (stat > 0) {
1745  Double_t val = haa_zz->GetYaxis()->GetBinCenter(ny);
1746  nlA.Add(TMath::Nint(val));
1747  sumA += val * stat;
1748  sum += stat;
1749  }
1750  }
1751  sumA /= sum;
1752  Int_t nA = nlA.GetNValues();
1753  if (nA == 0) {
1754  Warning("CalculateDeltaE_EGrid", "no count for Z=%d", zz);
1755  }
1756  else {
1757  if (Zonly) {
1758  nlA.Clear();
1759  nlA.Add(TMath::Nint(sumA));
1760  }
1761  else {
1762  if (nA == 1) {
1763  Int_t aref = nlA.Last();
1764  nlA.Add(aref - 1);
1765  nlA.Add(aref + 1);
1766  }
1767  }
1768  part.SetZ(zz);
1769  nlA.Begin();
1770  while (!nlA.End()) {
1771  Int_t aa = nlA.Next();
1772  part.SetA(aa);
1773  if (part.IsKnown()) {
1774  addLineToGrid(idgrid, zz, aa, npoints);
1775  }
1776  }
1777  }
1778  }
1779  return idgrid;
1780 }
1781 
1782 
1783 
1799 
1801 {
1802  // Returns the Y-axis value in the 2D identification map containing isotope (Z,A)
1803  // corresponding to either the given X-axis/Eres value or the current X-axis value given by GetIDGridXCoord()
1804  // If no mass information is available, just give Z.
1805  //
1806  // In this (default) implementation this means scanning the ID grids associated with
1807  // this telescope until we find an identification line Z or (Z,A), and then interpolating
1808  // the Y-coordinate for the current X-coordinate value.
1809  //
1810  // Status variable can take one of following values:
1811  //
1812  // KVIDTelescope::kMeanDE_OK all OK
1813  // KVIDTelescope::kMeanDE_XtooSmall X-coordinate is smaller than smallest X-coordinate of ID line
1814  // KVIDTelescope::kMeanDE_XtooLarge X-coordinate is larger than largest X-coordinate of ID line
1815  // KVIDTelescope::kMeanDE_NoIdentifie No identifier found for Z or (Z,A)
1816 
1817  status = kMeanDE_OK;
1818  // loop over grids
1819  TIter next(GetListOfIDGrids());
1820  KVIDGrid* grid;
1821  KVIDLine* idline = 0;
1822  while ((grid = (KVIDGrid*)next())) {
1823  idline = (KVIDLine*)grid->GetIdentifier(Z, A);
1824  if (idline) break;
1825  }
1826  if (!idline) {
1827  status = kMeanDE_NoIdentifier;
1828  return -1.;
1829  }
1830  Double_t x, x1, y1, x2, y2; bool f;
1831  x = (Eres < 0 ? GetIDGridXCoord(grid,f) : Eres);
1832  idline->GetEndPoint(x2, y2);
1833  if (x > x2) {
1834  status = kMeanDE_XtooLarge;
1835  return -1;
1836  }
1837  idline->GetStartPoint(x1, y1);
1838  if (x < x1) {
1839  status = kMeanDE_XtooSmall;
1840  return -1.;
1841  }
1842  return idline->Eval(x);
1843 }
1844 
1845 
1846 
1847 
1856 
1858 {
1859  // Return kTRUE if energy of ION is > minimum incident energy required for identification
1860  // This theoretical limit is defined here to be the incident energy for which the
1861  // dE in the first detector of a dE-E telescope is maximum.
1862  // If EINC>0 it is assumed to be the energy of the ion just before the first detector
1863  // (case where ion would have to pass other detectors before reaching this telescope).
1864  //
1865  // If this is not a dE-E telescope, we return kTRUE by default.
1866 
1867  if (GetSize() < 2) return kTRUE;
1868 
1869  KVDetector* dEdet = GetDetector(1);
1870  Double_t emin = dEdet->GetEIncOfMaxDeltaE(ION->GetZ(), ION->GetA());
1871  if (EINC > 0.0) return (EINC > emin);
1872  return (ION->GetEnergy() > emin);
1873 }
1874 
1875 
1876 
1905 
1907 {
1908  // For filtering simulations
1909  //
1910  // \param IDR identification result object for this telescope
1911  // \param n the simulated particle currently being considered
1912  //
1913  // Set the IDR->Zident and IDR->Aident status of the identification of n.
1914  // In principle this depends on whether this telescope provides mass
1915  // identification or not, but this may depend on the particle's energy.
1916  //
1917  // In order to enable mass identification for certain telescopes without a dedicated
1918  // implementation (e.g. for simulating array response), put the following lines
1919  // in your .kvrootrc:
1920  //
1921  // [dataset].[telescope label].MassID: yes
1922  //
1923  // If you want to limit mass identification to certain values of Z and/or A,
1924  // add the following line:
1925  //
1926  // [dataset].[telescope label].MassID.Validity: [expression]
1927  //
1928  // where [expression] is some valid C++ boolean expression involving Z and/or A,
1929  // for example
1930  //
1931  // [dataset].[telescope label].MassID.Validity: (Z>3)&&(A<20)
1932  //
1933  // Then this expression will be tested here in order to determine particle
1934  // identification status
1935 
1936  KVNucleus test(IDR->Z, IDR->A);
1937  if (!test.IsKnown()) { // reject weird identifications as pile-up (in between lines)
1938  IDR->Zident = false;
1939  IDR->IDquality = 4;
1940  IDR->IDOK = false;
1941  return;
1942  }
1943  IDR->Zident = true;
1944  if (!HasMassID()) {
1945  IDR->Aident = false;
1946  }
1947  else {
1948  if (fMassIDValidity) IDR->Aident = fMassIDValidity->Test(n); // test expression for mass ID validity
1949  else IDR->Aident = true; // no expression set; all nuclei are identified in mass
1950  }
1951 }
1952 
1953 
1954 
1957 
1959 {
1960  // Open IdentificationBilan.dat file with given path
1961 
1963  fgIdentificationBilan = new TEnv(path);
1964 }
1965 
1966 
1967 
1970 
1972 {
1973  // Set status of ID Telescope for given system
1974  if (!(fgIdentificationBilan->GetValue(Form("%s.%s", system.Data(), GetName()), kTRUE))) ResetBit(kReadyForID);
1975 }
1976 
1977 
int Int_t
unsigned int UInt_t
#define d(i)
#define f(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:670
ValType GetDataSetEnv(const Char_t *type, const ValType &defval) const
Definition: KVDataSet.h:268
Bool_t HasCalibIdentInfos() const
Definition: KVDataSet.h:393
TString GetFullPathToDataSetFile(const Char_t *filename)
Definition: KVDataSet.cpp:1875
static Bool_t FindDataSetFile(const TString &dataset, const Char_t *filename)
Definition: KVDataSet.cpp:1909
Signal output from a mathematical combination of other signals.
Base class for output signal data produced by a detector.
virtual Bool_t IsFired() const
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:844
void SetRunList(const char *runlist)
Definition: KVIDGraph.h:184
void SetName(const char *name) override
Definition: KVIDGraph.h:168
void AddIDTelescope(KVBase *t)
Definition: KVIDGraph.h:435
KVIDentifier * GetCut(const Char_t *name) const
Definition: KVIDGraph.h:308
virtual void Identify(Double_t, Double_t, KVIdentificationResult *) const =0
const Char_t * GetName() const override
Definition: KVIDGraph.cpp:1344
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:1281
virtual void SetOnlyZId(Bool_t yes=kTRUE)
Definition: KVIDGraph.cpp:1508
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:85
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)
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()
KVGroup * GetGroup() const
virtual Double_t GetPedestalX(Option_t *opt="")
Double_t GetIDGridXCoord(KVIDGraph *, bool &) const
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)
bool GetIDGridCoords(Double_t &X, Double_t &Y, KVIDGraph *grid, Double_t x=-1, Double_t y=-1)
KVIDGraph * GetIDGrid()
void Print(Option_t *opt="") const override
Double_t GetIDGridYCoord(KVIDGraph *, bool &) const
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)
const Char_t * GetDefaultIDGridClass()
UInt_t GetSize() const
KVUnownedList fIDGrids
identification grid(s)
Definition: KVIDTelescope.h:94
KVUnownedList fDetectors
list of detectors in telescope
Definition: KVIDTelescope.h:92
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:93
static TEnv * fgIdentificationBilan
Definition: KVIDTelescope.h:88
virtual void SetIdentificationStatus(KVIdentificationResult *IDR, const KVNucleus *)
virtual void RemoveGrids()
void addLineToGrid(KVIDGrid *gg, int zz, int aa, int npoints, std::optional< Double_t > einc_max=std::nullopt)
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)