KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVIDTelescope.h
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"
27class KVGroup;
28class KVIDGraph;
29class KVIDGrid;
30class KVMultiDetArray;
32class TH2;
33
34#include <KVUnownedList.h>
35#include <unordered_map>
36
84class KVIDTelescope: public KVBase {
85
86
89
90protected:
94 enum {
95 kMassID = BIT(15), //set if telescope is capable of mass identification i.e. isotopic resolution
96 kReadyForID = BIT(16) //set if telescope is ready and able for identification. set in Initialize()
97 };
99
100 void SetLabelFromURI(const Char_t* uri);
101 std::unique_ptr<KVParticleCondition> fMassIDValidity;
102 KVDetectorSignal* GetSignalFromGridVar(const KVString& var, const KVString& axe, KVString& det_labels);
109 mutable std::unordered_map<KVIDGraph*, GraphCoords> fGraphCoords;
111
112 KVIDGrid* newGrid(bool onlyZ);
113 void addLineToGrid(KVIDGrid* gg, int zz, int aa, int npoints);
114
115public:
116
118
120 enum {
121 kCalibStatus_OK, // fine, OK, all detectors calculated and functioning properly
122 kCalibStatus_Calculated, // one or more detectors not calibrated/functioning, energy loss calculated
123 kCalibStatus_NoCalibrations, // no calibrations available for any detectors
124 kCalibStatus_Multihit, // one or more detectors hit by more than one particle, energy loss calculated for this detector
125 kCalibStatus_Coherency // comparison of calculated and measured energy loss in one or more detectors reveals presence of other particles
126 };
127
129 void init();
130 virtual void AddDetector(KVDetector* d);
131 const KVList* GetDetectors() const
132 {
133 return &fDetectors;
134 }
136 {
138 if (n > GetSize() || n < 1) {
139 Error("GetDetector(UInt_t n)", "n must be between 1 and %u",
140 GetSize());
141 return 0;
142 }
143 KVDetector* d = (KVDetector*) GetDetectors()->At((n - 1));
144 return d;
145 };
146 KVDetector* GetDetector(const Char_t* name) const;
148 {
150 return GetDetectorRank(d) > 0;
151 }
153 {
156
157 return GetDetectors()->IndexOf(kvd) + 1;
158 }
160 {
162 return GetDetectors()->GetSize();
163 }
164
165 KVGroup* GetGroup() const;
166 void SetGroup(KVGroup* kvg);
168
169 virtual TGraph* MakeIDLine(KVNucleus* nuc, Double_t Emin, Double_t Emax,
170 Double_t Estep = 0.0);
171
172 virtual void Initialize(void);
173
174 virtual Bool_t Identify(KVIdentificationResult*, Double_t x = -1., Double_t y = -1.);
175
177 virtual Int_t GetCalibStatus() const
178 {
184 return fCalibStatus;
185 }
186
187 virtual void Print(Option_t* opt = "") const;
188
189 void SetIDGrid(KVIDGraph*);
192 KVIDGraph* GetIDGrid(const Char_t*);
194 {
195 return &fIDGrids;
196 }
197
198 virtual void RemoveGrids();
199
200 virtual Double_t GetIDMapX(Option_t* opt = "");
201 virtual Double_t GetIDMapY(Option_t* opt = "");
202 virtual Double_t GetPedestalX(Option_t* opt = "");
203 virtual Double_t GetPedestalY(Option_t* opt = "");
204
207 void GetIDGridCoords(Double_t& X, Double_t& Y, KVIDGraph* grid, Double_t x = -1, Double_t y = -1)
208 {
213
214 Y = (y < 0. ? GetIDGridYCoord(grid) : y);
215 X = (x < 0. ? GetIDGridXCoord(grid) : x);
216 }
217 void SetHasMassID(Bool_t yes = kTRUE)
218 {
219 SetBit(kMassID, yes);
220 }
221
225 {
226 return TestBit(kMassID);
227 }
228
229 static KVIDTelescope* MakeIDTelescope(const Char_t* name);
230
232 virtual void RemoveIdentificationParameters();
233
234 void LoadIdentificationParameters(const Char_t* filename, const KVMultiDetArray* multidet);
235 void ReadIdentificationParameterFiles(const Char_t* filename, const KVMultiDetArray* multidet);
236
240 {
241 return TestBit(kReadyForID);
242 };
243
245 KVIDGrid* CalculateDeltaE_EGrid(const KVNameValueList& AperZ, Int_t npoints = 30, Double_t xfactor = 1.);
246 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.);
247 KVIDGrid* CalculateDeltaE_EGrid(TH2* haa_zz, Bool_t Zonly, Int_t npoints);
249 enum {
250 kMeanDE_OK, // all OK
251 kMeanDE_XtooSmall, // X-coordinate is smaller than smallest X-coordinate of ID line
252 kMeanDE_XtooLarge, // X-coordinate is larger than largest X-coordinate of ID line
253 kMeanDE_NoIdentifier // No identifier found for Z or (Z,A)
254 };
255 virtual Double_t GetMeanDEFromID(Int_t& status, Int_t Z, Int_t A = -1, Double_t Eres = -1.0);
257 {
260 return fIDCode;
261 }
262
263 virtual void SetIDCode(UShort_t c)
264 {
266 fIDCode = c;
267 }
268
269 virtual Bool_t CheckTheoreticalIdentificationThreshold(KVNucleus* /*ION*/, Double_t /*EINC*/ = 0.0);
270 virtual Bool_t CanIdentify(Int_t Z, Int_t /*A*/)
271 {
277 return (Z > 0);
278 }
281 {
287
288 return (GetSize() == 1 || (GetSize() == 2 && GetDetector(1)->GetNode()->GetNTraj() == 1
289 && GetDetector(2)->GetNode()->GetNTraj() == 1));
290 }
291
292 static void OpenIdentificationBilan(const TString& path);
293 void CheckIdentificationBilan(const TString& system);
294
295 ClassDef(KVIDTelescope, 5) //A delta-E - E identification telescope
296};
297
298#endif
int Int_t
unsigned int UInt_t
#define d(i)
#define c(i)
bool Bool_t
unsigned short UShort_t
unsigned char UChar_t
char Char_t
double Double_t
const char Option_t
#define ClassDef(name, id)
#define BIT(n)
#define X(type, name)
Base class for KaliVeda framework.
Definition KVBase.h:142
Base class for output signal data produced by a detector.
Base class for detector geometry description.
Definition KVDetector.h:160
Group of detectors which can be treated independently of all others in array.
Definition KVGroup.h:20
Base class for particle identification in a 2D map.
Definition KVIDGraph.h:32
Abstract base class for 2D identification grids in e.g. (dE,E) maps.
Definition KVIDGrid.h:74
Base class for all detectors or associations of detectors in array which can identify charged particl...
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)
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
void ReadIdentificationParameterFiles(const Char_t *filename, const KVMultiDetArray *multidet)
virtual Bool_t CheckTheoreticalIdentificationThreshold(KVNucleus *, Double_t=0.0)
const KVList * GetListOfIDGrids() const
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*)
virtual Double_t GetPedestalX(Option_t *opt="")
KVDetector * GetDetector(UInt_t n) const
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)
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
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)
Bool_t HasDetector(const KVDetector *d) const
KVUnownedList fDetectors
list of detectors in telescope
KVString GetDetectorLabelsForGridCoord(const KVString &axis) const
const KVList * GetDetectors() 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
Bool_t IsIndependent() const
static TEnv * fgIdentificationBilan
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:28
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:126
Strings used to represent a set of ranges of values.
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:73
Extended TList class which does not own its objects by default.
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