]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSNeuralTracker.h
Fixed some warnings
[u/mrichter/AliRoot.git] / ITS / AliITSNeuralTracker.h
1 #ifndef ALIITSNEURALTRACKER_H
2 #define ALIITSNEURALTRACKER_H
3
4 class TObjArray;
5 class TCanvas;
6
7 ///////////////////////////////////////////////////////////////////////
8 //
9 // AliITSneuralTracker:
10 //
11 // neural network MFT algorithm
12 // for track finding in ITS stand alone
13 // according to the Denby-Peterson model with adaptments to the
14 // ALICE multiplicity
15 //
16 ///////////////////////////////////////////////////////////////////////
17
18 class AliITSNeuralTracker : public TObject {
19
20 public:
21
22                  AliITSNeuralTracker();
23                                 AliITSNeuralTracker(AliITSNeuralTracker &t) : TObject((TObject&)t)
24                                 { Fatal("AliITSNeuralTracker", "No copy constructor allowed!");exit(0);}
25         virtual ~AliITSNeuralTracker();
26
27         // ******************************************************************************
28         // * Embedded utility class --> >>> NODE <<<
29         // ******************************************************************************
30         // * This class inherits from AliITSNeuralPoint and adds some
31         // * utility pointers for quick path-finding among neurons.
32         // ******************************************************************************
33         class AliITSNode : public AliITSNeuralPoint {
34         public:
35                 AliITSNode() 
36                         {fInnerOf = fOuterOf = fMatches = 0; fNext = fPrev = 0;}
37                         
38                 AliITSNode(AliITSNeuralPoint *p, Bool_t init = kTRUE) // declared inline
39                         : AliITSNeuralPoint(p)
40                         {
41                                 fInnerOf = fOuterOf = fMatches = 0;
42                                 fNext = fPrev = 0;
43                                 if (init) {
44                                         fInnerOf = new TObjArray;
45                                         fOuterOf = new TObjArray;
46                                         fMatches = new TObjArray;
47                                 }
48                         }
49                         
50                 AliITSNode(AliITSRecPoint *p, AliITSgeomMatrix *gm)
51                         : AliITSNeuralPoint(p,gm) 
52                         {fInnerOf = fOuterOf = fMatches = 0; fNext = fPrev = 0;}
53
54                 virtual  ~AliITSNode() 
55                         {fInnerOf = fOuterOf = fMatches = 0; fNext = fPrev = 0;}
56
57                 Double_t  ThetaDeg()                    {return GetTheta()*180.0/TMath::Pi();}
58
59                 Int_t     GetSector(Double_t sec_width) {return (Int_t)(GetPhi()/sec_width);}
60                 Int_t     GetThetaCell()                {return (Int_t)(ThetaDeg());}
61                 
62                 Int_t        fPosInTree;  // position in tree of converted points
63
64                 TObjArray   *fInnerOf; //!
65                 TObjArray   *fOuterOf; //! 
66                 TObjArray   *fMatches; //!
67
68                 AliITSNode  *fNext; //!
69                 AliITSNode  *fPrev; //!
70         };
71         // ******************************************************************************
72
73
74
75         // ******************************************************************************
76         // * Embedded utility class --> >>> NEURON <<<
77         // ******************************************************************************
78         // * Simple class implementing the neural unit.
79         // * Contains the activation and some other pointers
80         // * for purposes similar to the ones in AliITSnode.
81         // ******************************************************************************
82         class AliITSneuron : public TObject {
83         public:
84                             AliITSneuron():fUsed(0),fActivation(0.),fInner(0),fOuter(0),fGain(0) { }
85                 virtual    ~AliITSneuron() {fInner=fOuter=0;fGain=0;}
86
87                 Double_t    Weight(AliITSneuron *n);
88                 void        Add2Gain(AliITSneuron *n, Double_t mult_const, Double_t exponent);
89
90                 Int_t             fUsed;        //  utility flag
91                 Double_t          fActivation;  //  Activation value
92                 AliITSNode       *fInner;       //! inner point
93                 AliITSNode       *fOuter;       //! outer point
94                 TObjArray        *fGain;        //! list of sequenced units
95         };
96         // ******************************************************************************
97
98
99
100         // ******************************************************************************
101         // * Embedded utility class --> >>> CONNECTION <<<
102         // ******************************************************************************
103         // * Used to implement the neural weighted connection
104         // * in such a way to speed up the retrieval of the
105         // * links among neuron, for a fast update procedure.
106         // ******************************************************************************
107         class AliITSlink : public TObject {
108         public:
109                          AliITSlink() : fWeight(0.), fLinked(0) { }
110                 virtual ~AliITSlink()   {fLinked = 0;}
111
112                 Double_t      fWeight;  //  Weight value
113                 AliITSneuron *fLinked;  //! the connected neuron
114         };
115         // ******************************************************************************
116
117
118         // Cut related setters
119
120         void     SetHelixMatchCuts(Double_t *min, Double_t *max);
121         void     SetThetaCuts2D(Double_t *min, Double_t *max);
122         void     SetThetaCuts3D(Double_t *min, Double_t *max);
123         void     SetCurvatureCuts(Int_t n, Double_t *cuts);
124         void     SetVertex(Double_t x, Double_t y, Double_t z)  {fVX=x; fVY=y; fVZ=z;}
125         void     SetPolarInterval(Double_t dtheta) {fPolarInterval=dtheta;}
126
127         // Neural work-flow related setters
128
129         void     SetActThreshold(Double_t val)            {fActMinimum = val;}
130         void     SetWeightExponent(Double_t a)            {fAlignExponent = a;}
131         void     SetGainToCostRatio(Double_t a)           {fGain2CostRatio = a;}
132         void     SetInitInterval(Double_t a, Double_t b)  {fEdge1 = a; fEdge2 = b;}
133         void     SetTemperature(Double_t a)               {fTemperature = a;}
134         void     SetVariationLimit(Double_t a)            {fStabThreshold = a;}
135
136         // Points array arrangement & control
137
138         void     CreateArrayStructure(Int_t nsecs);
139         Int_t    ArrangePoints(TTree *pts_tree);
140         void     StoreAbsoluteMatches();
141         Bool_t   PassCurvCut(AliITSNode *p1, AliITSNode *p2, Int_t curv_idx, Double_t vx, Double_t vy, Double_t vz);
142         Int_t    PassAllCuts(AliITSNode *p1, AliITSNode *p2, Int_t curv_idx, Double_t vx, Double_t vy, Double_t vz);
143         void     PrintPoints();
144         void     PrintMatches(Bool_t stop = kTRUE);
145
146         // Neural tracker work-flow
147
148         void     NeuralTracking(const char* filesave, TCanvas*& display);
149         void     Display(TCanvas*& canvas) const;
150         void     ResetNodes(Int_t isector);
151         Int_t    CreateNeurons(Int_t sector, Int_t curv);  // create neurons
152         Int_t    LinkNeurons() const;          // create neural connections
153         Double_t Activate(AliITSneuron* &n);   // calculates the new neural activation
154         Bool_t   Update();                     // an updating cycle
155         void     CleanNetwork();               // removes deactivated units and resolves competitions
156         Int_t    Save(Int_t sector_idx);       // save found tracks for # sector
157         TTree*   GetChains()                   {return fChains;}
158         void     WriteChains()                 {fChains->Write();}
159
160 private:
161
162         Bool_t       CheckOccupation() const; 
163
164         Int_t        fSectorNum;            //  number of azymuthal sectors
165         Double_t     fSectorWidth;          //  width of an azymuthal sector (in RADIANS) [used internally]
166         Double_t     fPolarInterval;        //  width of a polar sector (in DEGREES)
167
168         Double_t     fThetaCut2DMin[5];     //  lower edge of theta cut range (in DEGREES)
169         Double_t     fThetaCut2DMax[5];     //  upper edge of theta cut range (in DEGREES)
170         Double_t     fThetaCut3DMin[5];     //  lower edge of theta cut range (in DEGREES)
171         Double_t     fThetaCut3DMax[5];     //  upper edge of theta cut range (in DEGREES)
172         Double_t     fHelixMatchCutMin[5];  //  lower edge of helix matching cut range
173         Double_t     fHelixMatchCutMax[5];  //  lower edge of helix matching cut range
174         Int_t        fCurvNum;              //  # of curvature cut steps
175         Double_t    *fCurvCut;              //! value of all curvature cuts
176
177         Bool_t       fStructureOK;          // flag to avoid MANY segfault errors
178
179         Double_t     fVX, fVY, fVZ;         //  estimated vertex coords (for helix matching cut)
180
181         Double_t     fActMinimum;           //  minimum activation to turn 'on' the unit at the end
182         Double_t     fEdge1, fEdge2;        //  initialization interval for activations
183
184         Double_t     fTemperature;          //  logistic function temperature parameter
185         Double_t     fStabThreshold;        //  stability threshold
186         Double_t     fGain2CostRatio;       //  ratio between inhibitory and excitory contributions
187         Double_t     fAlignExponent;        //  alignment-dependent weight term
188
189         Int_t        fPoint[6];             //  Track point in layers
190         TTree       *fChains;               //! Output tree
191
192         TObjArray   *fPoints[6][180];       //! recpoints arranged into sectors for processing
193         TObjArray   *fNeurons;              //! neurons
194
195         ClassDef(AliITSNeuralTracker, 1)
196 };
197
198
199 ////////////////////////////////////////////////////////////////////////////////
200
201 #endif