KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVIDTelescope.h
Go to the documentation of this file.
1 /***************************************************************************
2 $Id: KVIDTelescope.h,v 1.33 2009/04/01 15:58:10 ebonnet Exp $
3  KVIDTelescope.h - description
4  -------------------
5  begin : Wed Jun 18 2003
6  copyright : (C) 2003 by J.D Frankland
7  email : frankland@ganil.fr
8  ***************************************************************************/
9 
10 /***************************************************************************
11  * *
12  * This program is free software; you can redistribute it and/or modify *
13  * it under the terms of the GNU General Public License as published by *
14  * the Free Software Foundation; either version 2 of the License, or *
15  * (at your option) any later version. *
16  * *
17  ***************************************************************************/
18 
19 #ifndef KVIDTELESCOPE_H
20 #define KVIDTELESCOPE_H
21 
22 #include "KVBase.h"
23 #include "KVDetector.h"
24 #include "TGraph.h"
25 #include "KVParticleCondition.h"
27 class KVGroup;
28 class KVIDGraph;
29 class KVIDGrid;
30 class KVMultiDetArray;
32 class TH2;
33 
34 #include <KVUnownedList.h>
35 #include <unordered_map>
36 
83 class KVIDTelescope: public KVBase {
84 
85 
88 
89 protected:
93  enum {
94  kMassID = BIT(15), //set if telescope is capable of mass identification i.e. isotopic resolution
95  kReadyForID = BIT(16) //set if telescope is ready and able for identification. set in Initialize()
96  };
98 
99  void SetLabelFromURI(const Char_t* uri);
100  std::unique_ptr<KVParticleCondition> fMassIDValidity;
101  KVDetectorSignal* GetSignalFromGridVar(const KVString& var, const KVString& axe, KVString& det_labels);
102  struct GraphCoords {
107  };
108  mutable std::unordered_map<KVIDGraph*, GraphCoords> fGraphCoords;
110 
111  KVIDGrid* newGrid(bool onlyZ);
112  void addLineToGrid(KVIDGrid* gg, int zz, int aa, int npoints);
113 
114 public:
115 
117 
119  enum {
120  kCalibStatus_OK, // fine, OK, all detectors calculated and functioning properly
121  kCalibStatus_Calculated, // one or more detectors not calibrated/functioning, energy loss calculated
122  kCalibStatus_NoCalibrations, // no calibrations available for any detectors
123  kCalibStatus_Multihit, // one or more detectors hit by more than one particle, energy loss calculated for this detector
124  kCalibStatus_Coherency // comparison of calculated and measured energy loss in one or more detectors reveals presence of other particles
125  };
126 
127  KVIDTelescope();
128  void init();
129  virtual void AddDetector(KVDetector* d);
130  const KVList* GetDetectors() const
131  {
132  return &fDetectors;
133  }
135  {
137  if (n > GetSize() || n < 1) {
138  Error("GetDetector(UInt_t n)", "n must be between 1 and %u",
139  GetSize());
140  return 0;
141  }
142  KVDetector* d = (KVDetector*) GetDetectors()->At((n - 1));
143  return d;
144  };
145  KVDetector* GetDetector(const Char_t* name) const;
146  Bool_t HasDetector(const KVDetector* d) const
147  {
149  return GetDetectorRank(d) > 0;
150  }
152  {
155 
156  return GetDetectors()->IndexOf(kvd) + 1;
157  }
158  UInt_t GetSize() const
159  {
161  return GetDetectors()->GetSize();
162  }
163 
164  KVGroup* GetGroup() const;
165  void SetGroup(KVGroup* kvg);
167 
168  virtual TGraph* MakeIDLine(KVNucleus* nuc, Double_t Emin, Double_t Emax,
169  Double_t Estep = 0.0);
170 
171  virtual void Initialize(void);
172 
173  virtual Bool_t Identify(KVIdentificationResult*, Double_t x = -1., Double_t y = -1.);
174 
176  virtual Int_t GetCalibStatus() const
177  {
183  return fCalibStatus;
184  }
185 
186  virtual void Print(Option_t* opt = "") const;
187 
188  void SetIDGrid(KVIDGraph*);
189  KVIDGraph* GetIDGrid();
191  KVIDGraph* GetIDGrid(const Char_t*);
192  const KVList* GetListOfIDGrids() const
193  {
194  return &fIDGrids;
195  }
196 
197  virtual void RemoveGrids();
198 
199  virtual Double_t GetIDMapX(Option_t* opt = "");
200  virtual Double_t GetIDMapY(Option_t* opt = "");
201  virtual Double_t GetPedestalX(Option_t* opt = "");
202  virtual Double_t GetPedestalY(Option_t* opt = "");
203 
206  void GetIDGridCoords(Double_t& X, Double_t& Y, KVIDGraph* grid, Double_t x = -1, Double_t y = -1)
207  {
212 
213  Y = (y < 0. ? GetIDGridYCoord(grid) : y);
214  X = (x < 0. ? GetIDGridXCoord(grid) : x);
215  }
216  void SetHasMassID(Bool_t yes = kTRUE)
217  {
218  SetBit(kMassID, yes);
219  }
220 
224  {
225  return TestBit(kMassID);
226  }
227 
228  static KVIDTelescope* MakeIDTelescope(const Char_t* name);
229 
231  virtual void RemoveIdentificationParameters();
232 
233  void LoadIdentificationParameters(const Char_t* filename, const KVMultiDetArray* multidet);
234  void ReadIdentificationParameterFiles(const Char_t* filename, const KVMultiDetArray* multidet);
235 
239  {
240  return TestBit(kReadyForID);
241  };
242 
243  const Char_t* GetDefaultIDGridClass();
244  KVIDGrid* CalculateDeltaE_EGrid(const KVNameValueList& AperZ, Int_t npoints = 30, Double_t xfactor = 1.);
245  KVIDGrid* CalculateDeltaE_EGrid(const KVNumberList& Zrange, Int_t deltaMasse, Int_t npoints, Double_t lifetime = -10/*s*/, UChar_t massformula = 0, Double_t xfactor = 1.);
246  KVIDGrid* CalculateDeltaE_EGrid(TH2* haa_zz, Bool_t Zonly, Int_t npoints);
248  enum {
249  kMeanDE_OK, // all OK
250  kMeanDE_XtooSmall, // X-coordinate is smaller than smallest X-coordinate of ID line
251  kMeanDE_XtooLarge, // X-coordinate is larger than largest X-coordinate of ID line
252  kMeanDE_NoIdentifier // No identifier found for Z or (Z,A)
253  };
254  virtual Double_t GetMeanDEFromID(Int_t& status, Int_t Z, Int_t A = -1, Double_t Eres = -1.0);
255  virtual UShort_t GetIDCode()
256  {
259  return fIDCode;
260  }
261 
262  virtual void SetIDCode(UShort_t c)
263  {
265  fIDCode = c;
266  }
267 
268  virtual Bool_t CheckTheoreticalIdentificationThreshold(KVNucleus* /*ION*/, Double_t /*EINC*/ = 0.0);
269  virtual Bool_t CanIdentify(Int_t Z, Int_t /*A*/)
270  {
276  return (Z > 0);
277  }
280  {
286 
287  return (GetSize() == 1 || (GetSize() == 2 && GetDetector(1)->GetNode()->GetNTraj() == 1
288  && GetDetector(2)->GetNode()->GetNTraj() == 1));
289  }
290 
291  static void OpenIdentificationBilan(const TString& path);
292  void CheckIdentificationBilan(const TString& system);
293 
294  ClassDef(KVIDTelescope, 5) //A delta-E - E identification telescope
295 };
296 
297 #endif
int Int_t
unsigned int UInt_t
#define d(i)
#define c(i)
unsigned short UShort_t
unsigned char UChar_t
char Char_t
bool Bool_t
double Double_t
const char Option_t
#define ClassDef(name, id)
#define BIT(n)
Base class for KaliVeda framework.
Definition: KVBase.h:141
Base class for output signal data produced by a detector.
Base class for detector geometry description.
Definition: KVDetector.h:159
Group of detectors which can be treated independently of all others in array.
Definition: KVGroup.h:19
Base class for particle identification in a 2D map.
Definition: KVIDGraph.h:31
Abstract base class for 2D identification grids in e.g. (dE,E) maps.
Definition: KVIDGrid.h:73
Base class for all detectors or associations of detectors in array which can identify charged particl...
Definition: KVIDTelescope.h:83
KVIDGrid * newGrid(bool onlyZ)
void LoadIdentificationParameters(const Char_t *filename, const KVMultiDetArray *multidet)
This method add to the gIDGridManager list the identification grids.
virtual Double_t GetIDMapY(Option_t *opt="")
virtual Bool_t IsReadyForID()
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)
@ kCalibStatus_NoCalibrations
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
void ReadIdentificationParameterFiles(const Char_t *filename, const KVMultiDetArray *multidet)
virtual Bool_t CheckTheoreticalIdentificationThreshold(KVNucleus *, Double_t=0.0)
UInt_t GetGroupNumber()
void addLineToGrid(KVIDGrid *gg, int zz, int aa, int npoints)
virtual void SetIDCode(UShort_t c)
KVGroup * GetGroup() const
virtual Bool_t CanIdentify(Int_t Z, Int_t)
Int_t fCalibStatus
temporary variable holding status code for last call to Calibrate(KVReconstructedNucleus*)
Definition: KVIDTelescope.h:97
virtual Double_t GetPedestalX(Option_t *opt="")
virtual void CalculateParticleEnergy(KVReconstructedNucleus *nuc)
virtual void AddDetector(KVDetector *d)
virtual Double_t GetMeanDEFromID(Int_t &status, Int_t Z, Int_t A=-1, Double_t Eres=-1.0)
const KVList * GetDetectors() const
KVDetectorSignal * GetSignalFromGridVar(const KVString &var, const KVString &axe, KVString &det_labels)
KVIDGraph * GetIDGrid()
virtual void Print(Option_t *opt="") const
std::unordered_map< KVIDGraph *, GraphCoords > fGraphCoords
X/Y coordinates from detector signals for ID maps.
UShort_t fIDCode
general id code corresponding to correct identification by this type of telescope
Definition: KVIDTelescope.h:87
virtual Int_t GetCalibStatus() const
virtual UShort_t GetIDCode()
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()
virtual void SetIdentificationStatus(KVReconstructedNucleus *)
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:92
Bool_t HasDetector(const KVDetector *d) const
KVUnownedList fDetectors
list of detectors in telescope
Definition: KVIDTelescope.h:90
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
UInt_t GetDetectorRank(const KVDetector *kvd) const
KVGroup * fGroup
group to which telescope belongs
Definition: KVIDTelescope.h:91
Bool_t IsIndependent() const
static TEnv * fgIdentificationBilan
Definition: KVIDTelescope.h:86
virtual void RemoveGrids()
virtual TGraph * MakeIDLine(KVNucleus *nuc, Double_t Emin, Double_t Emax, Double_t Estep=0.0)
Full result of one attempted particle identification.
Extended TList class which owns its objects by default.
Definition: KVList.h:27
Base class for describing the geometry of a detector array.
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
Nuclei reconstructed from data measured by a detector array .
virtual Int_t GetSize() const
virtual TObject * At(Int_t idx) const
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
Extended TList class which does not own its objects by default.
Definition: KVUnownedList.h:16
void SetBit(UInt_t f)
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual Int_t IndexOf(const TObject *obj) const
Double_t y[n]
Double_t x[n]
const Int_t n
KVDetectorSignal * fVarY
KVDetectorSignal * fVarX