StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
PeakWindow.h
1 /*
2 Author: David Kapukchyan
3 @[October 20, 2022]
4 > Added doxygen style comments
5 
6 @[June 27, 2022]
7 > Added a TMarker for the peak position in drawing PeakWindow objects. Also added options to the Draw() and Paint() methods for this object.
8 
9 @[June 26, 2022]
10 > Fixed SetPeak to work for the extreme cases where the peak point is the first or last point
11 
12 @[June 12, 2022]
13 > Implemented TObject's "Clone" and "Copy" methods in line with *PeakAna* for ease of copying these objects.
14 
15 @[May 10, 2022]
16 > Added TLine objects here for drawing. As well they're accessor functions and styling options.
17 
18 @[February 28, 2022]
19 > Added Sinh(*mse*) into chirality function as both a multiplicative factor and addition to the "hyperbola" as described in Feb. 25, 2022. It did not work as expected since the "hyperbola" form goes to zero but so does the Sinh(*mse*) and so as a multiplicative factor that means whether the slope is zero or the "hyperbola" is zero you get low chirality since it should not be a one or the either case but a mixture. The sum also did not work as this meant the probablity function maxes at 0.5 instead of 1. The bigger issue is that the "hyperbola" term was dominating the chirality and thus the Sinh(*mse*) was not contributing at all to the chirality. The main issue that arises when one takes away the "hyperbola" part, is that the slope is still not sufficient as noisy data can easily cause the slope to be zero since the start and end points will match up. This is a "false" matching in the sense that the overall behavior of the curve could still be increasing but a large noise level causes the next point to be as low as the start point. In order to solve this one may need to take into account the noise level (e.g. "sigma" in tunnel scale) and treat it as some kind of error on the start and end points and then propagate that to the uncertainty in slope and use that as a part of the chirality ??
20 
21 @[February 25, 2022]
22 > Testing of "peak" chirality function of the form (mpe - mse)^2 - (mps - mse)^2/(scale)^2; essentially a hyperbola where *mpe* is slope of line containing points (mPeakX,mPeakY) and (mEndX,mEndY), *mps* is slope of line containing points (mPeakX,mPeakY) and (mStartX,mStartY), *mse* is slope of line containing (mEndY,mEndX) and (mStartX,mStartY); in all cases order of points mattter; only one scale is needed since scale>1 shifts to the "start" and 0<scale<1 shifts to the "end". This form was okay but only measures how far away the peak position is from the center of the window. In future want to also take into account *mse* into the chirality ??
23 
24 @[February 22, 2022]
25 > Added a "peak" chirality function and functions to check if a given PeakWindow should be merged based on chirality
26 
27 @[February 18, 2022]
28 > Moved *PeakWindow* to it's implementation file as it's become a more composite structure
29 */
36 #ifndef PEAKWINDOW_H
37 #define PEAKWINDOW_H
38 
39 //C++ Headers
40 #include <iostream>
41 #include <vector>
42 #include <utility> //std::move (Needs C++ 11)
43 
44 //ROOT Headers
45 #include "TObject.h"
46 #include "TKey.h"
47 #include "TH1.h"
48 #include "TMarker.h"
49 #include "TGraph.h"
50 #include "TLine.h"
51 #include "TMath.h"
52 
53 #include "St_base/StMessMgr.h"
54 
55 class PeakWindow : public TObject{
56 public:
61  PeakWindow();
62  PeakWindow(Double_t start, Double_t end);
63  PeakWindow(const PeakWindow& oldpeak);
64  PeakWindow& operator=(const PeakWindow& rhs);
65  virtual ~PeakWindow();
66 
67  virtual void Copy(TObject& obj) const;
68  virtual TObject* Clone(const char* newname="") const;
69 
70  //To ease computation time I store the x,y values of the peak window/range rather than the graph point
71  Double_t mStartX;
72  Double_t mEndX;
73  Double_t mStartY;
74  Double_t mEndY;
75  Int_t mP_Peak;
76  Double_t mPeakX;
77  Double_t mPeakY;
78 
79  void SetWindow(Double_t s, Double_t e);
80  void GetWindow(Double_t &s, Double_t &e) const;
81 
86  void SetPeak(TGraph* gdata);
87 
95  static PeakWindow Combine(const PeakWindow &leftpeak, const PeakWindow &rightpeak, bool keepleft=true);
102  virtual void Combine( const PeakWindow &other, bool keepthis=true );
103 
115  virtual UShort_t CompareTo( const PeakWindow& other ) const;
116  Double_t StartEndLineSlope() const;
117  Double_t StartEndSlopeUncertainty(Double_t sigma) const;
118  Double_t StartEndLineYint() const;
119 
126  Double_t MidPoint(TGraph* graph=0) const;
127  virtual Double_t SlopeChirality(Double_t scale) const;//function to compute chirality factor for startendslope
128  virtual Double_t PeakChirality(Double_t slopescale, Double_t peakscale) const;//for peakscale 1 means peak is centered in window, peakscale>1 peak center shifted towards start, peakscale<1 peak center shifted towards end (peakscale>0)
129  virtual Double_t PeakChiralityProb(Double_t probscale, Double_t chirality) const;
130  virtual Double_t PeakChiralityProb(Double_t probscale, Double_t peakscale, Double_t chirscale) const;//Returns probabilty not to merge or probabilty is a "real" peak
131 
145  virtual Double_t PeakTunnelProb(TGraph* graph, Double_t scale=1.0, Double_t sigma=1.0 ) const;
146  //bool ChirMerge(const Rtools::PeakWindow& other, Double_t scalechir=1) const;//returns true if "this" should be merged with "other" PeakWindow according to PeakChirality
147  //bool ChirMergeLeft(const Rtools::PeakWindow& left, const Rtools::PeakWindow& right, Double_t scalechir=1) const;//based on chirality - returns true if should be merged with "left" PeakWindow; false if should be merged with "right" PeakWindow
148  virtual void Reset(Double_t start, Double_t end);
149  virtual void Print(Option_t* opt="") const;
150 
159  virtual void Draw(Option_t* opt="");
160  virtual void Paint(Option_t* opt="");
161 
162  //Options for drawing
163  TLine* GetStartLine(Double_t ymin=0, Double_t ymax=0);
164  Color_t GetStartLineColor() const;
165  Style_t GetStartLineStyle() const;
166  Width_t GetStartLineWidth() const;
167  void SetStartLineColor(Color_t color);
168  void SetStartLineColorAlpha(Color_t color,Float_t alpha);
169  void SetStartLineStyle(Style_t style);
170  void SetStartLineWidth(Width_t width);
171 
172  TMarker* GetPeakMarker();
173  Color_t GetPeakMarkerColor() const;
174  Style_t GetPeakMarkerStyle() const;
175  Size_t GetPeakMarkerSize() const;
176  void SetPeakMarkerColor(Color_t color);
177  void SetPeakMarkerColorAlpha(Color_t color, Float_t alpha);
178  void SetPeakMarkerStyle(Style_t style);
179  void SetPeakMarkerSize(Size_t size);
180 
181  TLine* GetEndLine(Double_t ymin=0, Double_t ymax=0);
182  Color_t GetEndLineColor() const;
183  Style_t GetEndLineStyle() const;
184  Width_t GetEndLineWidth() const;
185  void SetEndLineColor(Color_t color);
186  void SetEndLineColorAlpha(Color_t color,Float_t alpha);
187  void SetEndLineStyle(Style_t style);
188  void SetEndLineWidth(Width_t width);
189 
190 protected:
191  TLine* mStartLine;
192  TMarker* mPeakMarker;
193  TLine* mEndLine;
194 
195  ClassDef(PeakWindow,3);
196 };
197 
198 #endif
static PeakWindow Combine(const PeakWindow &leftpeak, const PeakWindow &rightpeak, bool keepleft=true)
combine two PeakWindow objects
Definition: PeakWindow.cxx:156
PeakWindow & operator=(const PeakWindow &rhs)
Assignment operator.
Definition: PeakWindow.cxx:23
TLine * GetStartLine(Double_t ymin=0, Double_t ymax=0)
Create and return a TLine for the start of the peak window.
Definition: PeakWindow.cxx:317
TLine * mStartLine
TLine for drawing the start of the peak window.
Definition: PeakWindow.h:191
void SetWindow(Double_t s, Double_t e)
Definition: PeakWindow.cxx:70
PeakWindow()
Default Constructor.
Definition: PeakWindow.cxx:5
Double_t StartEndLineYint() const
Computes the y-intercept of the line formed by the points (mStartX,mStartY) and (mEndX,mEndY)
Definition: PeakWindow.cxx:189
virtual void Copy(TObject &obj) const
Only copies variables, to copy TLines use Clone()
Definition: PeakWindow.cxx:33
Double_t StartEndLineSlope() const
Computes the slope of the line formed by the points (mStartX,mStartY) and (mEndX,mEndY) ...
Definition: PeakWindow.cxx:174
TMarker * GetPeakMarker()
Create and return a TMarker to mark the location of the peak.
Definition: PeakWindow.cxx:349
Double_t mPeakX
x-value of peak position as determined by SetPeak()
Definition: PeakWindow.h:76
virtual ~PeakWindow()
Destructor.
Definition: PeakWindow.cxx:63
virtual TObject * Clone(const char *newname="") const
Clone whole object, name is irrelevant.
Definition: PeakWindow.cxx:51
Int_t mP_Peak
Point Number of peak in a TGraph object (P for point), point is such that slope with previous point w...
Definition: PeakWindow.h:75
virtual void Print(Option_t *opt="") const
Prints information about PeakWindow.
Definition: PeakWindow.cxx:96
void GetWindow(Double_t &s, Double_t &e) const
Definition: PeakWindow.cxx:76
Double_t StartEndSlopeUncertainty(Double_t sigma) const
Uncertainty int the slope of the line formed by the points (mStartX,mStartY) and (mEndX,mEndY)
Definition: PeakWindow.cxx:181
TLine * mEndLine
TLine for drawing the end of the peak window.
Definition: PeakWindow.h:193
virtual void Reset(Double_t start, Double_t end)
Reset PeakWindow to constructor state.
Definition: PeakWindow.cxx:79
Double_t mStartY
y value associated with mStartX
Definition: PeakWindow.h:73
Double_t MidPoint(TGraph *graph=0) const
Computes the the line formed by the points (mStartX,mStartY) and (mEndX,mEndY) and evaluates that lin...
Definition: PeakWindow.cxx:194
Double_t mStartX
x value for start of the peak window
Definition: PeakWindow.h:71
void SetPeak(TGraph *gdata)
sets mPeakX based on mP_Peak using line of slopes from points (mP_Peak-1,mP_Peak) and (mP_Peak...
Definition: PeakWindow.cxx:112
virtual void Paint(Option_t *opt="")
paint method see Draw() for options
Definition: PeakWindow.cxx:293
virtual UShort_t CompareTo(const PeakWindow &other) const
compare two PeakWindow objects
Definition: PeakWindow.cxx:271
Double_t mEndX
x value for end of the peak window
Definition: PeakWindow.h:72
virtual void Draw(Option_t *opt="")
Draw the PeakWindow.
Definition: PeakWindow.cxx:288
Double_t mEndY
y value associated with mEndX
Definition: PeakWindow.h:74
Double_t mPeakY
y-value at mP_Peak
Definition: PeakWindow.h:77
TLine * GetEndLine(Double_t ymin=0, Double_t ymax=0)
Create and return a TLine for the end of the peak window.
Definition: PeakWindow.cxx:373
virtual Double_t PeakTunnelProb(TGraph *graph, Double_t scale=1.0, Double_t sigma=1.0) const
Compute probablity that a given PeakWindow is a real peak.
Definition: PeakWindow.cxx:247
TMarker * mPeakMarker
TMarker for drawing the peak location.
Definition: PeakWindow.h:192