Belle II Software development
OPALFitterGSL.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * Forked from https://github.com/iLCSoft/MarlinKinfit *
6 * *
7 * Further information about the fit engine and the user interface *
8 * provided in MarlinKinfit can be found at *
9 * https://www.desy.de/~blist/kinfit/doc/html/ *
10 * and in the LCNotes LC-TOOL-2009-001 and LC-TOOL-2009-004 available *
11 * from http://www-flc.desy.de/lcnotes/ *
12 * *
13 * See git log for contributors and copyright holders. *
14 * This file is licensed under LGPL-3.0, see LICENSE.md. *
15 **************************************************************************/
16
17#ifndef __OPALFITTERGSL_H
18#define __OPALFITTERGSL_H
19
20#include<vector>
21#include "analysis/OrcaKinFit/BaseFitter.h"
22
23#include <gsl/gsl_vector.h>
24#include <gsl/gsl_matrix.h>
25#include <gsl/gsl_permutation.h>
26
27namespace Belle2 {
33 namespace OrcaKinFit {
34
35 class BaseFitObject;
36 class BaseConstraint;
88 class OPALFitterGSL : public BaseFitter {
89 public:
91 virtual ~OPALFitterGSL();
92 virtual double fit() override;
93
95
104 virtual int getError() const override;
105
107 virtual double getProbability() const override;
109 virtual double getChi2() const override;
111 virtual int getDoF() const override;
113 virtual int getIterations() const override;
114
116 virtual int getNcon() const;
117
119 virtual int getNsoft() const;
120
122 virtual int getNpar() const;
123
125 virtual int getNunm() const;
126
127 virtual bool initialize() override;
128
130 virtual void setDebug(int debuglevel);
131
132 protected:
133
134
135 virtual bool updateFitObjects(double etaxi[]);
136 enum {NPARMAX = 50, NCONMAX = 20, NUNMMAX = 20};
137
138 int npar; // total number of parameters
139 int nmea; // total number of measured parameters
140 int nunm; // total number of unmeasured parameters
141 int ncon; // total number of constraints
142 int ierr; // Error status
143 int nit; // Number of iterations
144
145 double fitprob; // fit probability
146 double chi2; // final chi2
147
148 static void ini_gsl_permutation(gsl_permutation*& p, unsigned int size);
149 static void ini_gsl_vector(gsl_vector*& v, int unsigned size);
150 static void ini_gsl_matrix(gsl_matrix*& m, int unsigned size1, unsigned int size2);
151
152 static void debug_print(gsl_matrix* m, const char* name);
153 static void debug_print(gsl_vector* v, const char* name);
154
155 private:
156 gsl_vector* f;
157 gsl_vector* r;
158 gsl_matrix* Fetaxi;
159 gsl_matrix* S;
160 gsl_matrix* Sinv;
161 gsl_matrix* SinvFxi;
162 gsl_matrix* SinvFeta;
163 gsl_matrix* W1;
164 gsl_matrix* G;
165 gsl_matrix* H;
166 gsl_matrix* HU;
167 gsl_matrix* IGV;
168 gsl_matrix* V;
169 gsl_matrix* VLU;
170 gsl_matrix* Vinv;
171 gsl_matrix* Vnew;
172 gsl_matrix* Minv;
173 gsl_matrix* dxdt;
174 gsl_matrix* Vdxdt;
175 gsl_vector* dxi;
176 gsl_vector* Fxidxi;
177 gsl_vector* lambda;
178 gsl_vector* FetaTlambda;
179 gsl_vector* etaxi;
180 gsl_vector* etasv;
181 gsl_vector* y;
182 gsl_vector* y_eta;
183 gsl_vector* Vinvy_eta;
184 // gsl_matrix *Feta; //<-- DJeans: not used/needed
185 gsl_matrix* FetaV;
186
187 gsl_permutation* permS;
188 gsl_permutation* permU;
189 gsl_permutation* permV;
190
191 int debug;
192
193 };
194
195 }// end OrcaKinFit namespace
197} // end Belle2 namespace
198
199#endif // __OPALFITTERGSL_H
Abstract base class for fitting engines of kinematic fits.
Definition: BaseFitter.h:47
Description of the fit algorithm and interface:
Definition: OPALFitterGSL.h:88
gsl_permutation * permU
Helper permutation vector.
virtual int getNsoft() const
Get the number of soft constraints of the last fit.
gsl_permutation * permV
Helper permutation vector.
virtual double fit() override
virtual int getNcon() const
Get the number of hard constraints of the last fit.
virtual int getError() const override
Return error code.
virtual int getNunm() const
Get the number of unmeasured parameters of the last fit.
virtual int getIterations() const override
Get the number of iterations of the last fit.
virtual double getProbability() const override
Get the fit probability of the last fit.
virtual void setDebug(int debuglevel)
Set the Debug Level.
virtual int getDoF() const override
Get the number of degrees of freedom of the last fit.
virtual double getChi2() const override
Get the chi**2 of the last fit.
virtual int getNpar() const
Get the number of all parameters of the last fit.
gsl_permutation * permS
Helper permutation vector.
Abstract base class for different kinds of events.