StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
TCFit.h
1 // @(#)root/base:$Name: $:$Id: TCFit.h,v 1.4 2010/01/27 21:36:37 perev Exp $
2 // Author: Victor Perev 05/07/2007
3 // Class for Fit with constrains.
4 // TCFit - fitter
5 // TCFitData - base class to define data
6 // TCFitV0 - Example of user class, inherited from TCFitData. This class
7 // Fit V0's
8 
9 
10 #ifndef ROOT_TCFit
11 #define ROOT_TCFit
12 #include "Rtypes.h"
13 #include "TNamed.h"
14 #include "TString.h"
15 #include "TMatrixD.h"
16 class TCFitData;
17 
18 
19 class TCFit : public TNamed {
20 public:
21 enum FitResult {
22  kFitOK = 0,
23  kBadAprx = 0x001, //Bad approximation, do not continue
24  kBadFcn = 0x002, //Fcn too big, do not continue
25  kBadCon = 0x004, //Cnr too big, do not continue
26  kTooItr = 0x008, //Too many iterations
27  kFatal = 0xfff //Unknown fatal error
28 };
29 enum FitAction {
30  kEndFit = 0,
31  kNextStep = 1,
32  kNextCut = 2,
33  kBadFit = 3,
34  kTooIter = 4
35 };
36 
37 public:
38  TCFit(const char *name,TCFitData *dat=0);
39 virtual ~TCFit();
40  void Reset();
41  Int_t SetData(TCFitData *dat);
42  void SetMaxIter(int maxiter) { fMaxIter=maxiter;}
43  void SetMaxCuts(int maxcuts) { fMaxCuts=maxcuts;}
44  void SetDebug(int deb=1) { fDebug=deb;}
45  int Fit();
46  double ErMx(int jcol,int jrow) const;
47 static void Test0();
48 private:
49  int PriStep(const char *tit="");
50 
51 private:
52 
53  int CheckIn();
54  int CheckOut();
55  int FitStep();
56  int CheckStep();
57  int CutStep();
58  int EndStep();
59 
60 
61  TCFitData *fDat;
62  TMatrixD *fBigM;
63  TMatrixD *fBigI;
64  TMatrixD *fBigB;
65  TMatrixD *fOldP;
66  TMatrixD *fAddP;
67  int fDebug;
68  int fUPars;
69  int fUMeas;
70  int fUCons;
71  int fIter;
72  int fMaxIter;
73  int fCuts;
74  int fMaxCuts;
75  int fAkt;
76  int fFitRes;
77  double fFcnQA[2];
78  double fConQA[2];
79  double fAddQA[2];
80 ClassDef(TCFit,0)
81 
82 };
83 
84 
85 class Deriv1st;
86 class Deriv2nd;
87 class TCFit;
88 
89 class TCFitData : public TNamed {
90 friend class Deriv1st;
91 friend class Deriv2nd;
92 public:
93 enum {kMEAS=0,kSLAC=1,kCNSR=2};
94 enum {kMaxId=100};
95 public:
96  TCFitData(const char *name, const char *title="");
97  virtual ~TCFitData();
98  void Reset();
99 
103 int AddPar(int tyPar,int idPar,double *par,int nPars=1,const char *name="",double tiny=0.);
104 
105 int GetId(const char *name) const; //get par id by name
106 
108 
109 int GetId(int jd) const; //get par id by jd
110 int GetJd(int id) const; //get par jd by id
111 const char *GetNam(int idx) const; //get par name by id
112 int GetType(int id) const; //get type(Meas,Slac,Constr) by dd
113 
114 void FixPar (int id,int yes=1); //(dis/en)able parameter/constrain
115 int IsFixed(int id) const ; //
116 
117 
118  virtual int Ready(); //Ready() could be overloaded
119  virtual int Approx()=0; //Approx() must be overloaded
120  virtual double Fcn()=0; //Fcn, must be overloaded
121  virtual double DFcn(int ipar); //dFcn/dPar could be overloaded
122  virtual double DDFcn(int ipar,int jpar); //d2Fcn/dPar1/dPar2 could be overloaded
123 
124  virtual double Con(int icon); //constrain, must be overloaded
125  virtual double DCon(int icon,int ipar); //dCon/dPar could be overloaded
126  //d2Con/dPar1/dPar2 could be overloaded
127  virtual void Update()=0; //Called when status changed,
128  //like parameters were changed after fit
129  //must be overloaded even if user does not need it.
130  //In last case it could be dummy routine
131  virtual void Print(const char *name) const;
132 
133 
134 //-----------------------------------------------------------------------------
135  void SetFitter(const TCFit *fitter) {fFitter = fitter;}
136  void SetFail(int ifail) {fFail = ifail ;}
137  int GetFail() const {return fFail;}
138  int Modified() const {return fModi;}
139  void Modify(int m=1) {fModi = m;}
140  void Evaluate();
141  double &GetPar (int ipar);
142  double GetPar (int ipar) const ;
143  double GetTiny(int ipar) const {return fTiny[ipar];}
144  double GetFcn () const {return fFcn[0];}
145  void SetFcn (double fcn) {fFcn[0] = fcn ;}
146  void SetFcn (double tiny,double big) {fFcn[1] = tiny; fFcn[2]=big;}
147  double GetBigFcn() const {return fFcn[2];}
148  double GetTinyFcn() const {return fFcn[1];}
149 
150  int GetNPars() const {return fNPars[0]+fNPars[1];}
151  int GetNMeas() const {return fNPars[0];}
152  int GetNSlac() const {return fNPars[1];}
153  int GetNCons() const {return fNPars[2];}
154  int GetUPars() const {return fNPars[0]-fNFixs[0]+fNPars[1]-fNFixs[1];}
155  int GetUMeas() const {return fNPars[0]-fNFixs[0];}
156  int GetUSlac() const {return fNPars[1]-fNFixs[1];}
157  int GetUCons() const {return fNPars[2]-fNFixs[2];}
158  int GetNDF() const {return GetUCons()-GetUSlac();}
159  double ErMx(int icol,int irow) const;
160 
161 //-----------------------------------------------------------------------------
162 private:
163 //-----------------------------------------------------------------------------
164 protected:
165 const TCFit *fFitter;
166 
167  char fBeg[1];
168  int fFail;
169  int fModi;
170  int fFlag;
171  int fNPars[3];
172  int fNFixs[3];
173 
174  double *fPars[kMaxId+1];
175  double fTiny[kMaxId+1];
176  short fTyps[kMaxId+1];
177  short fFixs[kMaxId+1];
178  char fIndx[kMaxId+1];
179  char fJndx[kMaxId+1];
180  double fFcn[3];
181  char fEnd[1];
182  TString fNams[kMaxId];
183 private:
184  Deriv1st *fD1st;
185  Deriv2nd *fD2nd;
186 ClassDef(TCFitData,0)
187 };
188 //end TCFitData
189 
190 
191 class TLorentzVector;
192 class THelixTrack;
193 class TVector3;
194 class TkErrs;
195 
196 class TkPars {
197 // Track parameters in the perigee point (near 2D point xvtx,yvtx)
198 public:
199  TkPars() { Reset(); SetHz();}
200 void Reset();
201 void Update(){curv = hz*ptin;}
202 void Print(const char *tit) const;
203  double *Arr() { return &dca;}
204 const double *Arr() const { return &dca;}
205  double P() const { return sqrt(1.+tanl*tanl)/fabs(ptin);}
206  double E() const { return sqrt((1.+tanl*tanl)/(ptin*ptin)+mass*mass);}
207 TLorentzVector P4() const;
208 void P4D(double D[4][5]) const;
209 TVector3 V3() const;
210  void Fill(THelixTrack &hlx);
211  void Set(const TVector3 &v3,const TVector3 &d3 ,double pts );
212  void Get( TVector3 *v3, TVector3 *d3=0,double *pts=0) const;
213  void Rand(const TkErrs &errs);
214  void SetHz(double factor=1.);
215 
216 static const char* Name(int mem);
217 static double Tiny(int mem);
218 
219 
220 TkPars &operator+=(const TkPars &a);
222  double dca,z;
224  double phi;
226  double ptin;
228  double tanl;
230  double curv;
232  double hz;
234  double mass;
235 };
236 
237 class TkErrs {
238 public:
239 
240  TkErrs() {Reset();}
241 void Reset();
242 void Set(int i,int j,double err);
243 double Get(int i,int j) const;
244 void Invert();
245 double Xi2(const TkPars &pars) const;
246 void Mpy(const TkPars &pars,double der[5]) const;
247 public:
248  double emx[15];
249 
250 };
251 
252 class VxPars {
253 public:
254 VxPars() {memset(x,0,sizeof(x));}
255 
256 TVector3 V3() const;
257 
258  double x[3];
259 public:
260 };
261 
262 class VxErrs {
263 public:
264  double emx[6];
265 };
266 
267 class TCFitV0 : public TCFitData {
268 public:
269 enum eTCFitV0 { kDCA_0= 0,kZ_0 = 1,kPHI_0= 2,kPTIN_0= 3,kTANL_0= 4
270  , kDCA_1=10,kZ_1 =11,kPHI_1=12,kPTIN_1=13,kTANL_1=14
271  , kLEN_0=20,kLEN_1=21,kLEN_2=22
272  , kCX_0 =30,kCY_0 =31,kCZ_0 =32
273  , kCX_1 =33,kCY_1 =34,kCZ_1 =35
274  , kCNRJ =36 };
275 public:
276  TCFitV0();
277  ~TCFitV0(){;}
278  virtual double Fcn();
279  virtual double DFcn(int ipar);
280  virtual double DDFcn(int ipar,int jpar);
281 
282  virtual double Con(int icon); //constrain, must be overloaded
283  virtual double DCon(int icon,int ipar); //dCon/dPar could be overloaded
284  int Ready();
285  int Approx();
286  virtual void Update();
287  virtual void Print(const char *name) const;
288 public:
289  void Reset();
290 static void Test(int mode=0);
291 
292 private:
293 public:
294 TkPars mTkBas[2]; //Base parameters of two tracks
295 TkErrs mTEBas[2]; //Errors of track parameters
296 TkPars mTkFit[2]; //Fitted parameters of two tracks
297 TkPars mTkDif[2]; //
298 
299 char mBeg[1];
300 char mReady;
301 
302 VxPars mVx; //Vertex parameters. Not used now
303 VxErrs mVE; //Vertex errorss. Not used now
304 double mLen[3]; //Lengths to V0 vertex along tracks. 2==V0 track
305 
306 double mConr[7]; //Constrains. 6=energy conservation
307 double mDFcn[2][5]; // dFcn/dTkPars
308 double mDConDL[3][3]; // dCon
309 double mP4d[2][4][5];
310 double mMas; //V0 mass
311 char mEnd[1];
312 
313 ClassDef(TCFitV0,0)
314 };
315 
316 #endif //ROOT_TCFit
317 
318 
319 
320 
321 
322 
323 
324 
325 
326 
327 
328 
329 
int fModi
fail flag, fit is impossible
Definition: TCFit.h:169
double tanl
tangent of the track momentum dip angle
Definition: TCFit.h:228
int fNPars[3]
&1=1 fcn calculated from error matrix
Definition: TCFit.h:171
double mass
Mass of track.
Definition: TCFit.h:234
Definition: TCFit.h:196
Definition: TCFit.h:19
short fTyps[kMaxId+1]
Tiny values still modifying fcn.
Definition: TCFit.h:176
Definition: TCFit.h:267
double curv
signed curvature [sign = sign(-qB)]
Definition: TCFit.h:230
double hz
Z component magnetic field in units Pt(Gev) = Hz * RCurv(cm)
Definition: TCFit.h:232
Definition: TCFit.h:237
double dca
point
Definition: TCFit.h:222
double phi
angle between track direction and X-axis in xy plane
Definition: TCFit.h:224
double * fPars[kMaxId+1]
number of "slack" parameters
Definition: TCFit.h:174
char fEnd[1]
Current value of fcn [1].
Definition: TCFit.h:181
double ptin
signed invert pt [sign = sign(-qB)]
Definition: TCFit.h:226
Definition: TCFit.h:252
int fNFixs[3]
number of "measured" ,slack, constrains
Definition: TCFit.h:172
Definition: TCFit.h:262
int AddPar(int tyPar, int idPar, double *par, int nPars=1, const char *name="", double tiny=0.)
Definition: TCFit.cxx:416