![]() |
KaliVeda
Toolkit for HIC analysis
|
Identification grid with lines corresponding to different nuclear isotopes (KVIDZALine)
Such a grid can be used either to identify simultaneously both the mass and charge of detected particles (if lines for several isotopes of each atomic number are drawn), or solely the charge (if only one isotope per Z is drawn, and if SetOnlyZId(kTRUE) is called).
After each identification attempt, the value returned by GetQualityCode() indicates whether the identification was successful or not. The meaning of the different codes depends on the type of identification.
In both cases, an acceptable identification is achieved if the quality code is kICODE0, kICODE1, kICODE2, or kICODE3.
Points with codes kICODE4 or kICODE5 are normally considered as "noise" (pile-up ?) and should be rejected.
Points which are (vertically) out of range for this grid have code kICODE6 (point too far below) or kICODE7 (point too far above).
Points with code kICODE8 are totally out of range.
Definition at line 66 of file KVIDZAGrid.h.
#include <KVIDZAGrid.h>
Public Types | |
enum | { kICODE0 , kICODE1 , kICODE2 , kICODE3 , kICODE4 , kICODE5 , kICODE6 , kICODE7 , kICODE8 , kICODE9 , kICODE10 } |
Public Member Functions | |
KVIDZAGrid () | |
default ctor. More... | |
KVIDZAGrid (const KVIDZAGrid &) | |
Copy constructor. More... | |
virtual | ~ KVIDZAGrid () |
virtual void | CalculateLineWidths () |
virtual void | Copy (TObject &) const |
Copy this to 'obj'. More... | |
virtual TClass * | DefaultIDLineClass () |
void | DrawLinesWithWidth () |
KVIDLine * | GetClosestLine () const |
Double_t | GetDistanceClosestLine () const |
UChar_t | GetIndexClosest () const |
Int_t | GetQualityCode () const |
virtual KVIDZALine * | GetZALine (Int_t z, Int_t a, Int_t &) const |
virtual KVIDZALine * | GetZLine (Int_t z, Int_t &) const |
Int_t | GetZmax () const |
KVIDZALine * | GetZmaxLine () const |
virtual void | Identify (Double_t x, Double_t y, KVIdentificationResult *) const |
virtual void | IdentZ (Double_t x, Double_t y, Double_t &Z) |
virtual void | IdentZA (Double_t x, Double_t y, Int_t &Z, Double_t &A) |
virtual void | Initialize () |
KVIDGraph * | MakeSubsetGraph (Int_t Zmin, Int_t Zmax, const Char_t *="") |
virtual void MakeEDeltaEZGrid(Int_t Zmin, Int_t Zmax, Int_t npoints=20, Double_t gamma = 2);//*MENU* More... | |
KVIDGraph * | MakeSubsetGraph (TList *, TClass *=0) |
void | ReCheckQuality (Int_t &Z, Double_t &A) |
void | RemoveLine (Int_t Z, Int_t A=-1) |
Remove and destroy identifier. More... | |
void | RemoveZLines (const Char_t *ZList) |
Remove and destroy identifiers. More... | |
void | SetManualWidth (Double_t manual_width=.3, Double_t manual_width_scaling=0.05) |
void | SetVarXVarY (const char *VarX, const char *VarY) |
![]() | |
KVIDGrid () | |
Default constructor. More... | |
virtual | ~ KVIDGrid () |
virtual TClass * | DefaultOKLineClass () |
KVIDLine * | FindNearestEmbracingIDLine (Double_t x, Double_t y, const Char_t *position, const Char_t *axis, Int_t &idx, Int_t &idx_min, Int_t &idx_max, Double_t &dist, Double_t &dist_min, Double_t &dist_max) const |
KVIDLine * | FindNearestIDLineFast (Double_t x, Double_t y, const Char_t *position, Int_t &idx, Int_t &idx_min, Int_t &idx_max, Double_t &dist, Double_t &dist_min, Double_t &dist_max) const |
KVIDLine * | FindNextEmbracingLine (Int_t &index, Int_t inc_index, Double_t x, Double_t y, const Char_t *axis) const |
Int_t | GetIDLinesEmbracingPoint (const Char_t *direction, Double_t x, Double_t y) const |
KVIDLine * | NewLine (const Char_t *idline_class="") |
![]() | |
KVIDGraph () | |
KVIDGraph (const KVIDGraph &) | |
Copy constructor. More... | |
virtual | ~KVIDGraph () |
virtual Bool_t | AcceptIDForTest (const KVIdentificationResult &idr) |
void | Add (TString, KVIDentifier *) |
KVIDentifier * | Add (TString, TString) |
virtual void | AddCut (KVIDentifier *cut) |
virtual void | AddIdentifier (KVIDentifier *id) |
void | AddIDTelescope (KVBase *t) |
void | AddIDTelescopes (const TList *) |
Associate this graph with all ID telescopes in list. More... | |
virtual void | AddInfo (KVIDentifier *info) |
void | AddParameter (char *Name, char *Value) |
virtual Int_t | CheckVersion (Int_t version) |
virtual void | Clear (Option_t *opt="") |
void | ClearListOfTelescopes () |
void | ClearPad (TVirtualPad *) |
virtual void | Delete (Option_t *option="") |
void | Draw (Option_t *opt="") |
virtual void | DrawAndAdd (const Char_t *type="ID", const Char_t *classname="KVIDentifier") |
virtual void | DrawClass () const |
virtual TObject * | DrawClone (Option_t *option="") const |
virtual void | DrawPanel () |
virtual void | Dump () const |
virtual Bool_t | ExistVersion (Int_t version) |
void | ExtendBeginningAllIdentLines (Double_t, Option_t *="") |
void | ExtendEndAllIdentLines (Double_t, Option_t *="") |
void | FindAxisLimits () |
Calculate X/Y min/max of all objects in graph. More... | |
virtual TFitResultPtr | Fit (const char *formula, Option_t *option="", Option_t *goption="", Axis_t xmin=0, Axis_t xmax=0) |
virtual TFitResultPtr | Fit (TF1 *f1, Option_t *option="", Option_t *goption="", Axis_t xmin=0, Axis_t xmax=0) |
virtual void | FitPanel () |
virtual void SetTitle(const char *title="") {TGraph::SetTitle(title);}; More... | |
KVIDentifier * | GetCut (const Char_t *name) const |
KVList * | GetCuts () |
const KVList * | GetCuts () const |
KVIDentifier * | GetIdentifier (const Char_t *name) const |
KVIDentifier * | GetIdentifier (Int_t Z, Int_t A) const |
KVIDentifier * | GetIdentifierAt (Int_t index) const |
KVList * | GetIdentifiers () |
const KVList * | GetIdentifiers () const |
const Char_t * | GetIDTelescopeLabel () const |
const TList * | GetIDTelescopes () const |
KVIDentifier * | GetInfo (const Char_t *name) const |
KVList * | GetInfos () |
const KVList * | GetInfos () const |
KVIDGraph * | GetLastSavedVersion () const |
Int_t | GetMassFormula () const |
const Char_t * | GetName () const |
const Char_t * | GetNamesOfIDTelescopes () const |
Int_t | GetNumberOfCuts () const |
Int_t | GetNumberOfIdentifiers () const |
Int_t | GetNumberOfInfos () const |
TVirtualPad * | GetPad () const |
KVNameValueList * | GetParameters () |
const KVNameValueList * | GetParameters () const |
TString | GetPattern (void) |
const Char_t * | GetRunList () const |
const KVNumberList & | GetRuns () const |
Axis_t | GetXmax () const |
Axis_t | GetXmin () const |
Double_t | GetXScaleFactor () |
Axis_t | GetYmax () const |
Axis_t | GetYmin () const |
Double_t | GetYScaleFactor () |
Bool_t | HandlesIDTelescope (KVBase *t) const |
virtual Bool_t | HasMassIDCapability () const |
void | Increment (Float_t x) |
virtual Int_t | InsertPoint () |
virtual void | Inspect () const |
TVirtualPad * | IsDrawn () const |
virtual Bool_t | IsIdentifiable (Double_t, Double_t, TString *rejected_by=nullptr) const |
Bool_t | IsSorted () const |
void | Modified () |
void | Print (Option_t *opt="") const |
Print out all objects in graph. More... | |
void | ReadAsciiFile (const Char_t *filename) |
virtual void | ReadAsciiFile_WP (Int_t version=-1) |
lecture des grilles avec version suivant un modele de fichier More... | |
virtual void | ReadFromAsciiFile (std::ifstream &gridfile) |
void | RemoveCut (KVIDentifier *) |
Remove and destroy cut. More... | |
void | RemoveIdentifier (KVIDentifier *) |
Remove and destroy identifier. More... | |
void | RemoveIDTelescope (KVBase *t) |
void | RemoveInfo (KVIDentifier *) |
Remove and destroy cut. More... | |
virtual Int_t | RemovePoint () |
virtual Int_t | RemovePoint (Int_t i) |
void | ResetDraw () |
void | ResetPad () |
void | RevertToLastSavedVersion () |
virtual void | SaveAs (const char *filename="", Option_t *option="") const |
void | Scale (TF1 *sx, TF1 *sy) |
virtual void | SetDrawOption (Option_t *option="") |
virtual void | SetEditable (Bool_t editable=kTRUE) |
virtual void | SetFillAttributes () |
virtual void | SetInfos (Double_t, Double_t, KVIdentificationResult *) const |
loop over KVIDGraph::fInfoZones to set flags in KVIdentificationResult More... | |
virtual void | SetLineAttributes () |
void | SetLineColor (Color_t lcolor) |
void | SetLineStyle (Style_t lstyle) |
void | SetLineWidth (Width_t lwidth) |
virtual void | SetMarkerAttributes () |
void | SetMassFormula (Int_t) |
virtual void | SetMassIDCapability (Bool_t yes=kTRUE) |
virtual void | SetMaximum (Double_t maximum=-1111) |
virtual void | SetMinimum (Double_t minimum=-1111) |
virtual void | SetName (const char *name) |
virtual void | SetOnlyZId (Bool_t yes=kTRUE) |
void | SetPattern (TString pattern) |
void | SetRunList (const char *runlist) |
void | SetRuns (const KVNumberList &nl) |
Set list of runs for which grid is valid. More... | |
virtual void | SetVarX (const char *v) |
virtual void | SetVarY (const char *v) |
void | SetXScaleFactor (Double_t=0) |
void | SetXVariable (const char *v) |
void SetIDTelescopes();// MENU={Hierarchy="Set.../ID Telescopes"} More... | |
void | SetYScaleFactor (Double_t=0) |
void | SetYVariable (const char *v) |
void | SortIdentifiers () |
virtual void | TestIdentification (TH2F *data, KVHashList &histos, KVNameValueList &histo_names) |
void | UnDraw () |
void | UpdateLastSavedVersion () |
update last saved version. mkae copy of current state of graph. More... | |
void | WriteAsciiFile (const Char_t *filename) |
Open, write and close ascii file containing this grid. More... | |
virtual void | WriteAsciiFile_WP (Int_t version=-1) |
virtual void | WriteToAsciiFile (std::ofstream &gridfile) |
Private Member Functions | |
virtual Bool_t | FindFourEmbracingLines (Double_t x, Double_t y, const Char_t *position) |
void | init () |
initialisation More... | |
void | SetZmax (Int_t z) |
Private Attributes | |
Int_t | Ainf |
Int_t | Ainfi |
Int_t | Aint |
mass of line used to identify particle More... | |
Int_t | Asup |
Int_t | Asups |
Double_t | dinf |
Double_t | dinfi |
Double_t | dsup |
Double_t | dsups |
KVIDLine * | fClosest |
closest line to last-identified point More... | |
Double_t | fDistanceClosest |
distance from point to closest line More... | |
Int_t | fICode |
code de retour More... | |
Int_t | fIdxClosest |
index of closest line in main list fIdentifiers More... | |
KVIDLine * | fLinf |
KVIDLine * | fLinfi |
KVIDLine * | fLsup |
KVIDLine * | fLsups |
UShort_t | fZMax |
largest Z of lines in grid More... | |
KVIDZALine * | fZMaxLine |
line with biggest Z and A More... | |
Int_t | kinf |
Int_t | kinfi |
Int_t | ksup |
Int_t | ksups |
used by IdentZA and IdentZ More... | |
Double_t | winf |
Double_t | winfi |
Double_t | wsup |
Double_t | wsups |
Int_t | Zinf |
Int_t | Zinfi |
Int_t | Zint |
Z of line used to identify particle. More... | |
Int_t | Zsup |
Int_t | Zsups |
Additional Inherited Members | |
![]() | |
static KVIDGraph * | AddGraphs (KVIDGraph *g1, Int_t id1_min, Int_t id1_max, KVIDGraph *g2, Int_t id2_min, Int_t id2_max) |
static Bool_t | GetAutoAdd () |
static KVIDGraph * | MakeIDGraph (const Char_t *) |
static void | SetAutoAdd (Bool_t yes=kTRUE) |
anonymous enum |
Enumerator | |
---|---|
kICODE0 | |
kICODE1 | |
kICODE2 | |
kICODE3 | |
kICODE4 | |
kICODE5 | |
kICODE6 | |
kICODE7 | |
kICODE8 | |
kICODE9 | |
kICODE10 |
Definition at line 101 of file KVIDZAGrid.h.
KVIDZAGrid::KVIDZAGrid | ( | ) |
default ctor.
Definition at line 33 of file KVIDZAGrid.cpp.
KVIDZAGrid::KVIDZAGrid | ( | const KVIDZAGrid & | grid | ) |
Copy constructor.
Definition at line 54 of file KVIDZAGrid.cpp.
|
virtual |
|
virtual |
Calculate natural "width" of each line in the grid. The lines in the grid are first sorted so that they are in order of ascending 'Y' i.e. first line is 1H, last line is the heaviest isotope (highest line).
Then, for a given line :
**** if the grid is to be used for A & Z identification (default)**** :
**** if the grid is to be used for Z identification (fOnlyZid=kTRUE)**** :
In each case we find D_L (the separation between the two lines at their extreme left) and D_R (their separation at extreme right). The width of the line is then calculated from these two using the method KVIDZALine::SetAsymWidth (which may be overridden in child classes). Info("CalculateLineWidths", "For grid %s (%s vs. %s, runs %s).", GetName(), GetVarY(), GetVarX(), GetRunList());
Reimplemented from KVIDGrid.
Definition at line 287 of file KVIDZAGrid.cpp.
|
virtual |
Copy this to 'obj'.
Reimplemented from KVIDGraph.
Reimplemented in KVIDZAFromZGrid.
Definition at line 66 of file KVIDZAGrid.cpp.
|
inlinevirtual |
Reimplemented from KVIDGrid.
Reimplemented in KVIDGCsI.
Definition at line 145 of file KVIDZAGrid.h.
void KVIDZAGrid::DrawLinesWithWidth | ( | ) |
This method displays the grid as in KVIDGrid::Draw, but the natural line widths are shown as error bars
Definition at line 459 of file KVIDZAGrid.cpp.
|
privatevirtual |
This method will locate (at most) four lines close to the point (x,y), the point must lie within the endpoints (in X) of each line (the lines "embrace" the point). Returns kTRUE if at least one line is found. Identification can then be carried out with either IdentZA or IdentZ (see Identify). Returns kFALSE if no lines are found (not even a closest embracing line).
We look for two lines above the point and two lines below the point, as in one of the following two cases:
---------------------— ksups ------------------—
closest —> ---------------------— ksup ------------------— X
X ---------------------— kinf ------------------— <— closest
---------------------— kinfi ------------------—
First we find the closest embracing line to the point, using FindNearestEmbracingIDLine. Then we search above and below for the other 'embracing' lines. Note that no condition is applied regarding the distances to these lines: the lines must have been sorted in order of increasing ordinate before hand in Initialize(), we simply use the order of lines in the list of identifiers. The Z, A, width and distance to each of these lines are stored in the variables Zsups, Asups, wsups, dsups etc. etc. to be used by IdentZA or IdentZ.
Definition at line 518 of file KVIDZAGrid.cpp.
|
inline |
Definition at line 159 of file KVIDZAGrid.h.
|
inline |
Definition at line 163 of file KVIDZAGrid.h.
|
inline |
Definition at line 167 of file KVIDZAGrid.h.
|
inlinevirtual |
Return quality code for previously-attempted identification Meanings of code values are given in class description
Reimplemented from KVIDGraph.
Definition at line 150 of file KVIDZAGrid.h.
|
virtual |
Returns ID line corresponding to nucleus (Z,A) and its index in fIDLines. First we use GetLine(z) to find a line with the right Z (in principal there are several in the grid), then search for the correct isotope among these.
Reimplemented in KVIDGCsI.
Definition at line 212 of file KVIDZAGrid.cpp.
|
virtual |
Returns ID line for which GetZ() returns 'z'. index=index of line found in fIDLines list (-1 if not found). To increase speed, this is done by dichotomy, rather than by looping over all the lines in the list.
Reimplemented in KVIDGCsI.
Definition at line 151 of file KVIDZAGrid.cpp.
|
inline |
Definition at line 133 of file KVIDZAGrid.h.
|
inline |
Definition at line 138 of file KVIDZAGrid.h.
|
virtual |
Fill the KVIdentificationResult object with the results of identification for point (x,y) corresponding to some physically measured quantities related to a reconstructed nucleus.
By default (OnlyZId()=kFALSE) this means identifying the Z & A of the nucleus. In this case, we consider that the nucleus' Z & A have been correctly measured if the 'quality code' returned by IdentZA() is < kICODE4: we set idr->Zident and idr->Aident to kTRUE if fICode<kICODE4
If OnlyZId()=kTRUE, only the Z of the nucleus is established. In this case, we consider that the nucleus' Z has been correctly measured if the 'quality code' returned by IdentZ() is < kICODE4, thus: we set idr->Zident to kTRUE if fICode<kICODE4 The mass idr->A is set to the mass of the nearest line.
For points lying between two lines of same Z and different A (fICode<kIDCode4) the "real" mass is given by interpolation between the two masses. The integer mass is the A of the line closest to the point. This means that the integer A is not always = nint("real" A), as for example if a grid is drawn with lines for 7Be & 9Be but not 8Be (usual case), then particles between the two lines can have "real" masses between 7.5 and 8.5, but their integer A will be =7 or =9, never 8.
Implements KVIDGraph.
Reimplemented in KVIDZAFromZGrid, KVIDGCsI, and KVIDGChIoSi.
Definition at line 1367 of file KVIDZAGrid.cpp.
|
virtual |
Finds Z & 'real Z' for point (x,y) once closest lines to point have been found (see GetNearestIDLine). This is is based on the algorithm developed by L. Tassan-Got in IdnCsOr, even the same variable names and comments have been used (as much as possible).
Definition at line 1076 of file KVIDZAGrid.cpp.
|
virtual |
Finds Z, A and 'real A' for point (x,y) once closest lines to point have been found by calling method FindFourEmbracingLines beforehand. This is a line-for-line copy of the latter part of IdnCsOr, even the same variable names and comments have been used (as much as possible).
Reimplemented in KVIDGCsI.
Definition at line 643 of file KVIDZAGrid.cpp.
|
private |
initialisation
Definition at line 80 of file KVIDZAGrid.cpp.
|
virtual |
General initialisation method for identification grid. This method MUST be called once before using the grid for identifications. The ID lines are sorted. The natural line widths of all ID lines are calculated. The line with the largest Z (Zmax line) is found.
Reimplemented from KVIDGrid.
Reimplemented in KVIDZAFromZGrid, KVIDGCsI, and KVIDGChIoSi.
Definition at line 1507 of file KVIDZAGrid.cpp.
KVIDGraph * KVIDZAGrid::MakeSubsetGraph | ( | Int_t | Zmin, |
Int_t | Zmax, | ||
const Char_t * | graph_class = "" |
||
) |
virtual void MakeEDeltaEZGrid(Int_t Zmin, Int_t Zmax, Int_t npoints=20, Double_t gamma = 2);//*MENU*
Create a new graph/grid using the subset of lines of this grid with Zmin <= Z <= Zmax. By default the new graph/grid will be of the same class as this one, unless graph_class !="", in which case it must contain the name of a class which derives from KVIDGraph. A clone of each line will be made and added to the new graph, which will have the same name and be associated with the same ID telescopes as this one.
Definition at line 1716 of file KVIDZAGrid.cpp.
KVIDGraph * KVIDZAGrid::MakeSubsetGraph | ( | TList * | lines, |
TClass * | graph_class = 0 |
||
) |
Create a new graph/grid using the subset of lines of this grid contained in TList 'lines'. By default the new graph/grid will be of the same class as this one, unless graph_class !=0, in which case it must contain the address of a TClass object representing a class which derives from KVIDGraph. A clone of each line will be made and added to the new graph, which will have the same name and be associated with the same ID telescopes as this one.
Definition at line 1671 of file KVIDZAGrid.cpp.
void KVIDZAGrid::ReCheckQuality | ( | Int_t & | Z, |
Double_t & | A | ||
) |
Recheck the identification quality using the 'manual' width set for each Z in the parameter list
If abs(Aint-Areal)> ManualWidth + ManualWidthScaling*sqrt(Z), kIDCode5 is returned. To be activated, a parameter "ManualWidth" should be added to the grid's parameter list.
Acceptable values for FAZIA are :
Example of use :
Definition at line 1022 of file KVIDZAGrid.cpp.
void KVIDZAGrid::RemoveLine | ( | Int_t | Z, |
Int_t | A = -1 |
||
) |
Remove and destroy identifier.
Definition at line 101 of file KVIDZAGrid.cpp.
void KVIDZAGrid::RemoveZLines | ( | const Char_t * | ZList | ) |
Remove and destroy identifiers.
Definition at line 131 of file KVIDZAGrid.cpp.
void KVIDZAGrid::SetManualWidth | ( | Double_t | manual_width = .3 , |
Double_t | manual_width_scaling = 0.05 |
||
) |
Definition at line 1062 of file KVIDZAGrid.cpp.
|
inline |
Definition at line 124 of file KVIDZAGrid.h.
|
inlineprivate |
Definition at line 73 of file KVIDZAGrid.h.
|
private |
Definition at line 91 of file KVIDZAGrid.h.
|
private |
Definition at line 91 of file KVIDZAGrid.h.
|
private |
mass of line used to identify particle
Definition at line 93 of file KVIDZAGrid.h.
|
private |
Definition at line 91 of file KVIDZAGrid.h.
|
private |
Definition at line 91 of file KVIDZAGrid.h.
|
private |
Definition at line 88 of file KVIDZAGrid.h.
|
private |
Definition at line 88 of file KVIDZAGrid.h.
|
private |
Definition at line 88 of file KVIDZAGrid.h.
|
private |
Definition at line 88 of file KVIDZAGrid.h.
|
private |
closest line to last-identified point
Definition at line 78 of file KVIDZAGrid.h.
|
private |
distance from point to closest line
Definition at line 83 of file KVIDZAGrid.h.
|
private |
code de retour
Definition at line 85 of file KVIDZAGrid.h.
|
private |
index of closest line in main list fIdentifiers
Definition at line 84 of file KVIDZAGrid.h.
|
private |
Definition at line 81 of file KVIDZAGrid.h.
|
private |
Definition at line 82 of file KVIDZAGrid.h.
|
private |
Definition at line 80 of file KVIDZAGrid.h.
|
private |
Definition at line 79 of file KVIDZAGrid.h.
|
private |
largest Z of lines in grid
Definition at line 70 of file KVIDZAGrid.h.
|
private |
line with biggest Z and A
Definition at line 71 of file KVIDZAGrid.h.
|
private |
Definition at line 87 of file KVIDZAGrid.h.
|
private |
Definition at line 87 of file KVIDZAGrid.h.
|
private |
Definition at line 87 of file KVIDZAGrid.h.
|
private |
used by IdentZA and IdentZ
Definition at line 87 of file KVIDZAGrid.h.
|
private |
Definition at line 89 of file KVIDZAGrid.h.
|
private |
Definition at line 89 of file KVIDZAGrid.h.
|
private |
Definition at line 89 of file KVIDZAGrid.h.
|
private |
Definition at line 89 of file KVIDZAGrid.h.
|
private |
Definition at line 90 of file KVIDZAGrid.h.
|
private |
Definition at line 90 of file KVIDZAGrid.h.
|
private |
Z of line used to identify particle.
Definition at line 94 of file KVIDZAGrid.h.
|
private |
Definition at line 90 of file KVIDZAGrid.h.
|
private |
Definition at line 90 of file KVIDZAGrid.h.