Belle II Software development
BaseFitObject Class Referenceabstract
Inheritance diagram for BaseFitObject:
ParticleFitObject ISRPhotonFitObject JetFitObject NeutrinoFitObject PxPyPzMFitObject

Public Member Functions

 BaseFitObject ()
 Default constructor.
 
 BaseFitObject (const BaseFitObject &rhs)
 Copy constructor.
 
BaseFitObjectoperator= (const BaseFitObject &rhs)
 Assignment.
 
virtual ~BaseFitObject ()
 Virtual destructor.
 
virtual BaseFitObjectcopy () const =0
 Return a new copy of itself.
 
virtual BaseFitObjectassign (const BaseFitObject &source)
 Assign from anther object, if of same type.
 
virtual bool setParam (int ilocal, double par_, bool measured_, bool fixed_=false)
 Set value and measured flag of parameter i; return: significant change.
 
virtual bool setParam (int ilocal, double par_)
 Set value of parameter ilocal; return: significant change.
 
virtual bool updateParams (double p[], int idim)
 Read values from global vector, readjust vector; return: significant change.
 
virtual bool setMParam (int ilocal, double mpar_)
 Set measured value of parameter ilocal; return: success.
 
virtual bool setError (int ilocal, double err_)
 Set error of parameter ilocal; return: success.
 
virtual bool setCov (int ilocal, int jlocal, double cov_)
 Set covariance of parameters ilocal and jlocal; return: success.
 
virtual bool setGlobalParNum (int ilocal, int iglobal)
 Set number of parameter ilocal in global list return true signals OK.
 
virtual bool fixParam (int ilocal, bool fix=true)
 Fix a parameter (fix=true), or release it (fix=false)
 
virtual bool releaseParam (int ilocal)
 Release a parameter.
 
virtual bool isParamFixed (int ilocal) const
 Returns whether parameter is fixed.
 
virtual double getParam (int ilocal) const
 Get current value of parameter ilocal.
 
virtual const char * getParamName (int ilocal) const =0
 Get name of parameter ilocal.
 
virtual const char * getName () const
 Get object's name.
 
void setName (const char *name_)
 Set object's name.
 
virtual double getMParam (int ilocal) const
 Get measured value of parameter ilocal.
 
virtual double getError (int ilocal) const
 Get error of parameter ilocal.
 
virtual double getCov (int ilocal, int jlocal) const
 Get covariance between parameters ilocal and jlocal.
 
virtual double getRho (int ilocal, int jlocal) const
 Get correlation coefficient between parameters ilocal and jlocal.
 
virtual bool isParamMeasured (int ilocal) const
 Get measured flag for parameter ilocal.
 
virtual int getGlobalParNum (int ilocal) const
 Get global parameter number of parameter ilocal.
 
virtual int getNPar () const =0
 Get total number of parameters of this FitObject.
 
virtual int getNMeasured () const
 Get number of measured parameters of this FitObject.
 
virtual int getNUnmeasured () const
 Get number of unmeasured parameters of this FitObject.
 
virtual int getNFree () const
 Get number of free parameters of this FitObject.
 
virtual int getNFixed () const
 Get number of fixed parameters of this FitObject.
 
virtual double getChi2 () const
 Get chi squared from measured and fitted parameters.
 
virtual double getDChi2DParam (int ilocal) const
 Get derivative of chi squared w.r.t. parameter ilocal.
 
virtual double getD2Chi2DParam2 (int ilocal, int jlocal) const
 Get second derivative of chi squared w.r.t. parameters ilocal1 and ilocal2.
 
virtual std::ostream & printParams (std::ostream &os) const
 print the parameters and errors
 
virtual std::ostream & printRhoValues (std::ostream &os) const
 print the correlation coefficients
 
virtual std::ostream & print1stDerivatives (std::ostream &os) const
 print the 1st derivatives wrt metaSet 0 (E, px, py, pz)
 
virtual std::ostream & print2ndDerivatives (std::ostream &os) const
 print the 2nd derivatives wrt metaSet 0 (E, px, py, pz)
 
virtual std::ostream & print (std::ostream &os) const =0
 print object to ostream
 
void invalidateCache () const
 invalidate any cached quantities
 
virtual void updateCache () const =0
 
virtual void addToGlobCov (double *glcov, int idim) const
 Add covariance matrix elements to global covariance matrix of size idim x idim.
 
virtual int addToGlobalChi2DerVector (double *y, int idim) const
 Add derivatives of chi squared to global derivative vector.
 
virtual void addToGlobalChi2DerMatrix (double *M, int idim) const
 Add 2nd derivatives of chi squared to global derivative matrix.
 
virtual void addToGlobalChi2DerVector (double *y, int idim, double lambda, double der[], int metaSet) const
 Add derivatives of momentum vector to global derivative vector.
 
virtual void addTo1stDerivatives (double M[], int idim, double der[], int kglobal, int metaSet) const
 
virtual void addTo2ndDerivatives (double der2[], int idim, double factor[], int metaSet) const
 
virtual void addTo2ndDerivatives (double M[], int idim, double lambda, double der[], int metaSet) const
 
virtual double getFirstDerivative_Meta_Local (int iMeta, int ilocal, int metaSet) const =0
 
virtual double getSecondDerivative_Meta_Local (int iMeta, int ilocal, int jlocal, int metaSet) const =0
 
virtual void initCov ()
 
virtual double getError2 (double der[], int metaset) const
 
virtual void getDerivatives (double der[], int idim) const =0
 

Protected Member Functions

virtual bool calculateCovInv () const
 Calculate the inverse of the covariance matrix.
 

Protected Attributes

char * name
 
double par [BaseDefs::MAXPAR]
 fit parameters
 
double mpar [BaseDefs::MAXPAR]
 measured parameters
 
bool measured [BaseDefs::MAXPAR]
 measured flag
 
bool fixed [BaseDefs::MAXPAR]
 fixed flag
 
int globalParNum [BaseDefs::MAXPAR]
 global parameter number for each parameter
 
double cov [BaseDefs::MAXPAR][BaseDefs::MAXPAR]
 local covariance matrix
 
double covinv [BaseDefs::MAXPAR][BaseDefs::MAXPAR]
 inverse pf local covariance matrix
 
bool covinvvalid
 flag for valid inverse covariance matrix
 
bool cachevalid
 flag for valid cache
 

Static Protected Attributes

static const double eps2 = 0.0001
 

Related Functions

(Note that these are not member functions.)

std::ostream & operator<< (std::ostream &os, const BaseFitObject &bfo)
 Prints out a BaseFitObject, using its print method.
 

Detailed Description

Definition at line 91 of file BaseFitObject.h.

Constructor & Destructor Documentation

◆ BaseFitObject() [1/2]

Default constructor.

Definition at line 38 of file BaseFitObject.cc.

38 : name(nullptr), par{}, mpar{}, measured{}, fixed{}, globalParNum{}, cov{}, covinv{}, covinvvalid(
39 false),
40 cachevalid(false)
41 {
42 setName("???");
44
45 for (int ilocal = 0; ilocal < BaseDefs::MAXPAR; ++ilocal) {
46 globalParNum[ilocal] = -1;
47 fixed[ilocal] = false;
48 for (int jlocal = 0; jlocal < BaseDefs::MAXPAR; ++jlocal)
49 cov[ilocal][jlocal] = 0;
50 }
51
52 }
bool covinvvalid
flag for valid inverse covariance matrix
bool fixed[BaseDefs::MAXPAR]
fixed flag
void setName(const char *name_)
Set object's name.
double covinv[BaseDefs::MAXPAR][BaseDefs::MAXPAR]
inverse pf local covariance matrix
bool cachevalid
flag for valid cache
int globalParNum[BaseDefs::MAXPAR]
global parameter number for each parameter
double par[BaseDefs::MAXPAR]
fit parameters
void invalidateCache() const
invalidate any cached quantities
double mpar[BaseDefs::MAXPAR]
measured parameters
double cov[BaseDefs::MAXPAR][BaseDefs::MAXPAR]
local covariance matrix
bool measured[BaseDefs::MAXPAR]
measured flag

◆ BaseFitObject() [2/2]

BaseFitObject ( const BaseFitObject rhs)

Copy constructor.

Parameters
rhsright hand side

Definition at line 54 of file BaseFitObject.cc.

55 : name(nullptr), par{}, mpar{}, measured{}, fixed{}, globalParNum{}, cov{}, covinv{}, covinvvalid(false), cachevalid(false)
56 {
58 }
virtual BaseFitObject & assign(const BaseFitObject &source)
Assign from anther object, if of same type.

◆ ~BaseFitObject()

~BaseFitObject ( )
virtual

Virtual destructor.

Definition at line 91 of file BaseFitObject.cc.

92 {
93 delete[] name;
94 }

Member Function Documentation

◆ addTo1stDerivatives()

void addTo1stDerivatives ( double  M[],
int  idim,
double  der[],
int  kglobal,
int  metaSet 
) const
virtual

Definition at line 515 of file BaseFitObject.cc.

517 {
518 if (!cachevalid) updateCache();
519 for (int ilocal = 0; ilocal < getNPar(); ilocal++) {
520 int iglobal = globalParNum[ilocal];
521 if (iglobal >= 0) {
522 for (int j = 0; j < BaseDefs::nMetaVars[metaSet]; j++) {
523 double x = der[j] * getFirstDerivative_Meta_Local(j, ilocal, metaSet);
524 M[idim * kglobal + iglobal] += x;
525 M[idim * iglobal + kglobal] += x;
526 }
527 }
528 }
529 return;
530 }
virtual int getNPar() const =0
Get total number of parameters of this FitObject.

◆ addTo2ndDerivatives() [1/2]

void addTo2ndDerivatives ( double  der2[],
int  idim,
double  factor[],
int  metaSet 
) const
virtual

Definition at line 534 of file BaseFitObject.cc.

539 {
540 if (!cachevalid) updateCache();
541 for (int ilocal = 0; ilocal < getNPar(); ilocal++) {
542 int iglobal = getGlobalParNum(ilocal);
543 if (iglobal < 0) continue;
544 for (int jlocal = ilocal; jlocal < getNPar(); jlocal++) {
545 int jglobal = getGlobalParNum(jlocal);
546 if (jglobal < 0) continue;
547 double sum(0);
548 for (int imeta = 0; imeta < BaseDefs::nMetaVars[metaSet]; imeta++) {
549 sum += factor[imeta] * getSecondDerivative_Meta_Local(imeta, ilocal, jlocal, metaSet);
550 }
551 der2[idim * iglobal + jglobal] += sum;
552 if (iglobal != jglobal) der2[idim * jglobal + iglobal] += sum;
553 }
554 }
555 return;
556 }
virtual int getGlobalParNum(int ilocal) const
Get global parameter number of parameter ilocal.

◆ addTo2ndDerivatives() [2/2]

void addTo2ndDerivatives ( double  M[],
int  idim,
double  lambda,
double  der[],
int  metaSet 
) const
virtual

Definition at line 559 of file BaseFitObject.cc.

560 {
561 double factor[BaseDefs::MAXINTERVARS];
562 for (int i = 0; i < BaseDefs::nMetaVars[metaSet]; i++) factor[i] = lambda * der[i];
563 addTo2ndDerivatives(M, idim, factor, metaSet);
564 return;
565 }

◆ addToGlobalChi2DerMatrix()

void addToGlobalChi2DerMatrix ( double *  M,
int  idim 
) const
virtual

Add 2nd derivatives of chi squared to global derivative matrix.

Parameters
MGlobal derivative matrix
idimFirst dimension of global derivative matrix

Definition at line 462 of file BaseFitObject.cc.

463 {
465 for (int ilocal = 0; ilocal < getNPar(); ++ilocal) {
466 if (!isParamFixed(ilocal) && isParamMeasured(ilocal)) {
467 int iglobal = getGlobalParNum(ilocal);
468 assert(iglobal >= 0 && iglobal < idim);
469 int ioffs = idim * iglobal;
470 for (int jlocal = 0; jlocal < getNPar(); ++jlocal) {
471 if (!isParamFixed(jlocal) && isParamMeasured(jlocal)) {
472 int jglobal = getGlobalParNum(jlocal);
473 assert(jglobal >= 0 && jglobal < idim);
474 M[ioffs + jglobal] += getD2Chi2DParam2(ilocal, jlocal);
475 }
476 }
477 }
478 }
479 }
virtual double getD2Chi2DParam2(int ilocal, int jlocal) const
Get second derivative of chi squared w.r.t. parameters ilocal1 and ilocal2.
virtual bool isParamFixed(int ilocal) const
Returns whether parameter is fixed.
virtual bool calculateCovInv() const
Calculate the inverse of the covariance matrix.
virtual bool isParamMeasured(int ilocal) const
Get measured flag for parameter ilocal.

◆ addToGlobalChi2DerVector() [1/2]

int addToGlobalChi2DerVector ( double *  y,
int  idim 
) const
virtual

Add derivatives of chi squared to global derivative vector.

Parameters
yVector of chi2 derivatives
idimVector size

Definition at line 482 of file BaseFitObject.cc.

483 {
484 // This adds the dChi2/dpar piece
485 assert(getNPar() <= BaseDefs::MAXPAR);
486 if (not covinvvalid and not calculateCovInv()) return 0;
487 for (int ilocal = 0; ilocal < getNPar(); ++ilocal) {
488 if (!isParamFixed(ilocal) && isParamMeasured(ilocal)) {
489 int iglobal = getGlobalParNum(ilocal);
490 assert(iglobal >= 0 && iglobal < idim);
491 y[iglobal] += getDChi2DParam(ilocal);
492 }
493 }
494 return 0;
495 }
virtual double getDChi2DParam(int ilocal) const
Get derivative of chi squared w.r.t. parameter ilocal.

◆ addToGlobalChi2DerVector() [2/2]

void addToGlobalChi2DerVector ( double *  y,
int  idim,
double  lambda,
double  der[],
int  metaSet 
) const
virtual

Add derivatives of momentum vector to global derivative vector.

Parameters
yVector of chi2 derivatives
idimVector size
lambdaThe lambda value
derderivatives of constraint wrt intermediate variables (e.g. 4-vector with dg/dE, dg/dpx, dg/dpy, dg/dpz)
metaSetwhich set of intermediate variables

Definition at line 498 of file BaseFitObject.cc.

500 {
501 // This adds the lambda * dConst/dpar piece
502 if (!cachevalid) updateCache();
503 for (int ilocal = 0; ilocal < getNPar(); ilocal++) {
504 int iglobal = globalParNum[ilocal];
505 if (iglobal >= 0) {
506 assert(iglobal < idim);
507 for (int j = 0; j < BaseDefs::nMetaVars[metaSet]; j++) {
508 y[iglobal] += lambda * der[j] * getFirstDerivative_Meta_Local(j, ilocal, metaSet);
509 }
510 }
511 }
512 }

◆ addToGlobCov()

void addToGlobCov ( double *  glcov,
int  idim 
) const
virtual

Add covariance matrix elements to global covariance matrix of size idim x idim.

Parameters
glcovGlobal covariance matrix
idimFirst dimension of global derivative matrix

Definition at line 217 of file BaseFitObject.cc.

218 {
219 for (int ilocal = 0; ilocal < getNPar(); ++ilocal) {
220 if (!isParamFixed(ilocal) && isParamMeasured(ilocal)) {
221 int iglobal = getGlobalParNum(ilocal);
222 assert(iglobal >= 0 && iglobal < idim);
223 int ioffs = idim * iglobal;
224 for (int jlocal = 0; jlocal < getNPar(); ++jlocal) {
225 if (!isParamFixed(jlocal) && isParamMeasured(jlocal)) {
226 int jglobal = getGlobalParNum(jlocal);
227 assert(jglobal >= 0 && jglobal < idim);
228 globCov[ioffs + jglobal] += getCov(ilocal, jlocal);
229 }
230 }
231 }
232 }
233 }
virtual double getCov(int ilocal, int jlocal) const
Get covariance between parameters ilocal and jlocal.

◆ assign()

BaseFitObject & assign ( const BaseFitObject source)
virtual

Assign from anther object, if of same type.

Parameters
sourceThe source object

Reimplemented in ISRPhotonFitObject, JetFitObject, NeutrinoFitObject, ParticleFitObject, and PxPyPzMFitObject.

Definition at line 68 of file BaseFitObject.cc.

69 {
70 if (&source != this) {
71 name = nullptr;
72 setName(source.name);
73 for (int i = 0; i < BaseDefs::MAXPAR; ++i) {
74 par[i] = source.par[i];
75 mpar[i] = source.mpar[i];
76 measured[i] = source.measured[i];
77 fixed[i] = source.fixed[i];
78 globalParNum[i] = source.globalParNum[i];
79 for (int j = 0; j < BaseDefs::MAXPAR; ++j) {
80 cov[i][j] = source.cov[i][j];
81 covinv[i][j] = source.covinv[i][j];
82 }
83 }
84 covinvvalid = false;
85 cachevalid = false;
86 }
87 return *this;
88 }

◆ calculateCovInv()

bool calculateCovInv ( ) const
protectedvirtual

Calculate the inverse of the covariance matrix.

Definition at line 235 of file BaseFitObject.cc.

236 {
237
238 // DANIEL added
239
240 int n = getNPar();
241
242 gsl_matrix* covm = gsl_matrix_alloc(n, n);
243 gsl_matrix_set_identity(covm);
244
245 for (int i = 0; i < n; ++i) {
246 if (isParamMeasured(i)) {
247 for (int j = 0; j < n; ++j) {
248 if (isParamMeasured(j)) {
249 gsl_matrix_set(covm, i, j, cov[i][j]);
250 }
251 }
252 }
253 }
254 gsl_error_handler_t* e = gsl_set_error_handler_off();
255 int result = gsl_linalg_cholesky_decomp(covm);
256 if (result == 0) result = gsl_linalg_cholesky_invert(covm);
257 gsl_set_error_handler(e);
258
259 for (int i = 0; i < n; ++i) {
260 for (int j = 0; j < i; ++j) {
261 covinv[i][j] = covinv[j][i] = gsl_matrix_get(covm, i, j);
262 }
263 covinv[i][i] = gsl_matrix_get(covm, i, i);
264 }
265
266 gsl_matrix_free(covm);
267 covinvvalid = (result == 0);
268
269 if (!covinvvalid) {
270 B2WARNING("ERROR, COULD NOT INVERT COV MATR!");
271
272 B2INFO("COV ");
273 for (int i = 0; i < n; ++i) {
274 if (isParamMeasured(i)) {
275 for (int j = 0; j < n; ++j) {
276 if (isParamMeasured(j)) {
277 B2INFO(cov[i][j] << " ");
278 }
279 }
280 }
281 }
282
283 B2INFO("CORREL ");
284 for (int i = 0; i < n; ++i) {
285 if (isParamMeasured(i)) {
286 for (int j = 0; j < n; ++j) {
287 if (isParamMeasured(j)) {
288 B2INFO(cov[i][j] / sqrt(cov[i][i]*cov[j][j]) << " ");
289 }
290 }
291 }
292 }
293
294 }
295
296 return covinvvalid;
297 }
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ copy()

virtual BaseFitObject * copy ( ) const
pure virtual

Return a new copy of itself.

Implemented in ISRPhotonFitObject, JetFitObject, NeutrinoFitObject, and PxPyPzMFitObject.

◆ fixParam()

bool fixParam ( int  ilocal,
bool  fix = true 
)
virtual

Fix a parameter (fix=true), or release it (fix=false)

Parameters
ilocalLocal parameter number
fixfix if true, release if false

Definition at line 352 of file BaseFitObject.cc.

353 {
354 assert(ilocal >= 0 && ilocal < getNPar());
355 return fixed [ilocal] = fix;
356 }

◆ getChi2()

double getChi2 ( ) const
virtual

Get chi squared from measured and fitted parameters.

Reimplemented in ParticleFitObject.

Definition at line 414 of file BaseFitObject.cc.

415 {
416 if (not covinvvalid and not calculateCovInv()) return -1;
417 double chi2 = 0;
418 static double resid[BaseDefs::MAXPAR];
419 static bool chi2contr[BaseDefs::MAXPAR];
420 for (int i = 0; i < getNPar(); ++i) {
421 resid[i] = par[i] - mpar[i];
422
423 B2INFO(" BaseFitObject::getChi2() " << i << " " << resid[i] << " " << covinv[i][i]);
424
425 if ((chi2contr[i] = (isParamMeasured(i) && !isParamFixed(i)))) {
426 chi2 += resid[i] * covinv[i][i] * resid[i];
427 for (int j = 0; j < getNPar(); ++j) {
428 if (j < i && chi2contr[j]) chi2 += 2 * resid[i] * covinv[i][j] * resid[j];
429 }
430 }
431 }
432 return chi2;
433 }

◆ getCov()

double getCov ( int  ilocal,
int  jlocal 
) const
virtual

Get covariance between parameters ilocal and jlocal.

Parameters
ilocalLocal parameter number i
jlocalLocal parameter number j

Reimplemented in JetFitObject.

Definition at line 389 of file BaseFitObject.cc.

390 {
391 assert(ilocal >= 0 && ilocal < getNPar());
392 assert(jlocal >= 0 && jlocal < getNPar());
393 return cov[ilocal][jlocal];
394 }

◆ getD2Chi2DParam2()

double getD2Chi2DParam2 ( int  ilocal,
int  jlocal 
) const
virtual

Get second derivative of chi squared w.r.t. parameters ilocal1 and ilocal2.

Parameters
ilocalLocal parameter number i
jlocalLocal parameter number j

Definition at line 449 of file BaseFitObject.cc.

450 {
451 assert(ilocal >= 0 && ilocal < getNPar());
452 assert(jlocal >= 0 && jlocal < getNPar());
453 if (isParamFixed(ilocal) || !isParamMeasured(ilocal) ||
454 isParamFixed(jlocal) || !isParamMeasured(jlocal))
455 return 0;
456 if (not covinvvalid and not calculateCovInv()) return 0;
457 return 2 * covinv[ilocal][jlocal]; // JL: ok in absence of soft constraints
458 }

◆ getDChi2DParam()

double getDChi2DParam ( int  ilocal) const
virtual

Get derivative of chi squared w.r.t. parameter ilocal.

Parameters
ilocalLocal parameter number

Definition at line 435 of file BaseFitObject.cc.

436 {
437 assert(ilocal >= 0 && ilocal < getNPar());
438 if (isParamFixed(ilocal) || !isParamMeasured(ilocal)) return 0;
439
440 if (not covinvvalid and not calculateCovInv()) return 0;
441
442 double result = 0;
443 for (int jlocal = 0; jlocal < getNPar(); jlocal++)
444 if (!isParamFixed(jlocal) && isParamMeasured(jlocal))
445 result += covinv[ilocal][jlocal] * (par[jlocal] - mpar[jlocal]);
446 return 2 * result;
447 }

◆ getError()

double getError ( int  ilocal) const
virtual

Get error of parameter ilocal.

Parameters
ilocalLocal parameter number

Reimplemented in JetFitObject.

Definition at line 383 of file BaseFitObject.cc.

384 {
385 assert(ilocal >= 0 && ilocal < getNPar());
386 return std::sqrt(cov[ilocal][ilocal]);
387 }

◆ getError2()

double getError2 ( double  der[],
int  metaset 
) const
virtual

Definition at line 597 of file BaseFitObject.cc.

598 {
599 if (!cachevalid) updateCache();
600 double totError(0);
601 for (int i = 0; i < BaseDefs::nMetaVars[metaSet]; i++) {
602 for (int j = 0; j < BaseDefs::nMetaVars[metaSet]; j++) {
603 double cov_i_j = 0; // this will become the covariance of intermediate variables i and j
604 for (int k = 0; k < getNPar(); k++) {
605 for (int l = 0; l < getNPar(); l++) {
606 cov_i_j += getFirstDerivative_Meta_Local(i, k, metaSet) * cov[k][l] * getFirstDerivative_Meta_Local(j, l, metaSet);
607 }
608 }
609 totError += der[i] * der[j] * cov_i_j;
610 }
611 }
612 return totError;
613 }

◆ getFirstDerivative_Meta_Local()

virtual double getFirstDerivative_Meta_Local ( int  iMeta,
int  ilocal,
int  metaSet 
) const
pure virtual

Implemented in JetFitObject.

◆ getGlobalParNum()

int getGlobalParNum ( int  ilocal) const
virtual

Get global parameter number of parameter ilocal.

Parameters
ilocalLocal parameter number

Definition at line 365 of file BaseFitObject.cc.

366 {
367 if (ilocal < 0 || ilocal >= getNPar()) return -1;
368 return globalParNum[ilocal];
369 }

◆ getMParam()

double getMParam ( int  ilocal) const
virtual

Get measured value of parameter ilocal.

Parameters
ilocalLocal parameter number

Definition at line 377 of file BaseFitObject.cc.

378 {
379 assert(ilocal >= 0 && ilocal < getNPar());
380 return mpar[ilocal];
381 }

◆ getName()

const char * getName ( ) const
virtual

Get object's name.

Definition at line 107 of file BaseFitObject.cc.

108 {
109 return name ? name : "???";
110 }

◆ getNFixed()

int getNFixed ( ) const
virtual

Get number of fixed parameters of this FitObject.

Definition at line 130 of file BaseFitObject.cc.

131 {
132 int nfixed = 0;
133 for (int i = 0; i < getNPar(); ++i) if (isParamFixed(i)) ++nfixed;
134 return nfixed;
135 }

◆ getNFree()

int getNFree ( ) const
virtual

Get number of free parameters of this FitObject.

Definition at line 124 of file BaseFitObject.cc.

125 {
126 int nfree = 0;
127 for (int i = 0; i < getNPar(); ++i) if (!isParamFixed(i)) ++nfree;
128 return nfree;
129 }

◆ getNMeasured()

int getNMeasured ( ) const
virtual

Get number of measured parameters of this FitObject.

Definition at line 112 of file BaseFitObject.cc.

113 {
114 int nmeasured = 0;
115 for (int i = 0; i < getNPar(); ++i) if (isParamMeasured(i) && !isParamFixed(i)) ++nmeasured;
116 return nmeasured;
117 }

◆ getNPar()

virtual int getNPar ( ) const
pure virtual

Get total number of parameters of this FitObject.

Implemented in ISRPhotonFitObject, JetFitObject, NeutrinoFitObject, and PxPyPzMFitObject.

◆ getNUnmeasured()

int getNUnmeasured ( ) const
virtual

Get number of unmeasured parameters of this FitObject.

Definition at line 118 of file BaseFitObject.cc.

119 {
120 int nunmeasrd = 0;
121 for (int i = 0; i < getNPar(); ++i) if (!isParamMeasured(i) && !isParamFixed(i)) ++nunmeasrd;
122 return nunmeasrd;
123 }

◆ getParam()

double getParam ( int  ilocal) const
virtual

Get current value of parameter ilocal.

Parameters
ilocalLocal parameter number

Definition at line 371 of file BaseFitObject.cc.

372 {
373 assert(ilocal >= 0 && ilocal < getNPar());
374 return par[ilocal];
375 }

◆ getParamName()

virtual const char * getParamName ( int  ilocal) const
pure virtual

Get name of parameter ilocal.

Parameters
ilocalLocal parameter number

Implemented in ISRPhotonFitObject, JetFitObject, NeutrinoFitObject, and PxPyPzMFitObject.

◆ getRho()

double getRho ( int  ilocal,
int  jlocal 
) const
virtual

Get correlation coefficient between parameters ilocal and jlocal.

Parameters
ilocalLocal parameter number i
jlocalLocal parameter number j

Definition at line 396 of file BaseFitObject.cc.

397 {
398 assert(ilocal >= 0 && ilocal < getNPar());
399 assert(jlocal >= 0 && jlocal < getNPar());
400 return cov[ilocal][jlocal] / std::sqrt(cov[ilocal][ilocal] * cov[jlocal][jlocal]);
401 }

◆ initCov()

void initCov ( )
virtual

Definition at line 587 of file BaseFitObject.cc.

588 {
589 for (int i = 0; i < getNPar(); ++i) {
590 for (int j = 0; j < getNPar(); ++j) {
591 cov[i][j] = static_cast<double>(i == j);
592 }
593 }
594 }

◆ invalidateCache()

void invalidateCache ( ) const
inline

invalidate any cached quantities

Definition at line 241 of file BaseFitObject.h.

241{cachevalid = false;};

◆ isParamFixed()

bool isParamFixed ( int  ilocal) const
virtual

Returns whether parameter is fixed.

Parameters
ilocalLocal parameter number

Definition at line 408 of file BaseFitObject.cc.

409 {
410 assert(ilocal >= 0 && ilocal < getNPar());
411 return fixed[ilocal];
412 }

◆ isParamMeasured()

bool isParamMeasured ( int  ilocal) const
virtual

Get measured flag for parameter ilocal.

Parameters
ilocalLocal parameter number

Definition at line 402 of file BaseFitObject.cc.

403 {
404 assert(ilocal >= 0 && ilocal < getNPar());
405 return measured[ilocal];
406 }

◆ operator=()

BaseFitObject & operator= ( const BaseFitObject rhs)

Assignment.

Parameters
rhsright hand side

Definition at line 60 of file BaseFitObject.cc.

61 {
62 if (this != &rhs) {
63 assign(rhs); // calls virtual function assign of derived class
64 }
65 return *this;
66 }

◆ print()

virtual std::ostream & print ( std::ostream &  os) const
pure virtual

print object to ostream

Parameters
osThe output stream

Implemented in ParticleFitObject.

◆ print1stDerivatives()

std::ostream & print1stDerivatives ( std::ostream &  os) const
virtual

print the 1st derivatives wrt metaSet 0 (E, px, py, pz)

Parameters
osThe output stream

Definition at line 164 of file BaseFitObject.cc.

165 {
166 int metaSet = 0;
167 os << "#";
168 for (int i = 0; i < BaseDefs::nMetaVars[metaSet]; ++i) {
169 for (int j = 0; j < getNPar(); ++j) {
170 if (i > 0 && j == 0) os << ",";
171 os << " " << getFirstDerivative_Meta_Local(i, j, metaSet);
172 }
173 }
174 os << "#" << std::endl;
175 return os;
176 }

◆ print2ndDerivatives()

std::ostream & print2ndDerivatives ( std::ostream &  os) const
virtual

print the 2nd derivatives wrt metaSet 0 (E, px, py, pz)

Parameters
osThe output stream

Definition at line 178 of file BaseFitObject.cc.

179 {
180 int metaSet = 0;
181 os << "##";
182 for (int i = 0; i < BaseDefs::nMetaVars[metaSet]; ++i) {
183 if (i > 0) os << std::endl << " ";
184 for (int j = 0; j < getNPar(); ++j) {
185 for (int k = 0; k < getNPar(); ++k) {
186// if (i>0 && k==0) os << ",";
187 os << " " << getSecondDerivative_Meta_Local(i, j, k, metaSet);
188 if (k == getNPar() - 1) os << ",";
189 }
190 }
191 }
192 os << "##" << std::endl;
193 return os;
194 }

◆ printParams()

std::ostream & printParams ( std::ostream &  os) const
virtual

print the parameters and errors

Parameters
osThe output stream

Definition at line 137 of file BaseFitObject.cc.

138 {
139 os << "(";
140 for (int i = 0; i < getNPar(); ++i) {
141 if (i > 0) os << ", ";
142 os << " " << getParam(i);
143 if (isParamFixed(i)) os << " fix";
144// else if (getError(i)>0) os << " \261 " << getError(i);
145 else if (getError(i) > 0) os << " +- " << getError(i); // " \261 " appeared as <B1> in ASCII ... Graham
146 }
147 os << ")";
148 return os;
149 }
virtual double getError(int ilocal) const
Get error of parameter ilocal.
virtual double getParam(int ilocal) const
Get current value of parameter ilocal.

◆ printRhoValues()

std::ostream & printRhoValues ( std::ostream &  os) const
virtual

print the correlation coefficients

Parameters
osThe output stream

Definition at line 151 of file BaseFitObject.cc.

152 {
153 os << "{";
154 for (int i = 0; i < getNPar(); ++i) {
155 for (int j = 0; j < getNPar(); ++j) {
156 if (i > 0 && j == 0) os << ",";
157 os << " " << getRho(i, j);
158 }
159 }
160 os << "}" << std::endl;
161 return os;
162 }
virtual double getRho(int ilocal, int jlocal) const
Get correlation coefficient between parameters ilocal and jlocal.

◆ releaseParam()

virtual bool releaseParam ( int  ilocal)
inlinevirtual

Release a parameter.

Parameters
ilocalLocal parameter number

Definition at line 157 of file BaseFitObject.h.

159 { return fixParam(ilocal, false); }
virtual bool fixParam(int ilocal, bool fix=true)
Fix a parameter (fix=true), or release it (fix=false)

◆ setCov()

bool setCov ( int  ilocal,
int  jlocal,
double  cov_ 
)
virtual

Set covariance of parameters ilocal and jlocal; return: success.

Parameters
ilocalLocal parameter number
jlocalLocal parameter number
cov_New error value

Definition at line 340 of file BaseFitObject.cc.

341 {
342 if (!isfinite(cov_)) return false;
343 assert(ilocal >= 0 && ilocal < getNPar());
344 assert(jlocal >= 0 && jlocal < getNPar());
346 covinvvalid = false;
347 cov[ilocal][jlocal] = cov[jlocal][ilocal] = cov_;
348 return true;
349 }

◆ setError()

bool setError ( int  ilocal,
double  err_ 
)
virtual

Set error of parameter ilocal; return: success.

Parameters
ilocalLocal parameter number
err_New error value

Definition at line 330 of file BaseFitObject.cc.

331 {
332 if (!isfinite(err_)) return false;
333 assert(ilocal >= 0 && ilocal < getNPar());
335 covinvvalid = false;
336 cov[ilocal][ilocal] = err_ * err_;
337 return true;
338 }

◆ setGlobalParNum()

bool setGlobalParNum ( int  ilocal,
int  iglobal 
)
virtual

Set number of parameter ilocal in global list return true signals OK.

Parameters
ilocalLocal parameter number
iglobalNew global parameter number

Definition at line 358 of file BaseFitObject.cc.

359 {
360 if (ilocal < 0 || ilocal >= getNPar()) return false;
361 globalParNum[ilocal] = iglobal;
362 return true;
363 }

◆ setMParam()

bool setMParam ( int  ilocal,
double  mpar_ 
)
virtual

Set measured value of parameter ilocal; return: success.

Parameters
ilocalLocal parameter number
mpar_New measured parameter value

Definition at line 320 of file BaseFitObject.cc.

321 {
322 if (!isfinite(mpar_)) return false;
323 assert(ilocal >= 0 && ilocal < getNPar());
324 if (mpar[ilocal] == mpar_) return true;
326 mpar[ilocal] = mpar_;
327 return true;
328 }

◆ setName()

void setName ( const char *  name_)

Set object's name.

Definition at line 98 of file BaseFitObject.cc.

99 {
100 if (name_ == nullptr) return;
101 size_t l = strlen(name_);
102 if (name) delete[] name;
103 name = new char[l + 1];
104 strcpy(name, name_);
105 }

◆ setParam() [1/2]

bool setParam ( int  ilocal,
double  par_ 
)
virtual

Set value of parameter ilocal; return: significant change.

Parameters
ilocalLocal parameter number
par_New parameter value

Definition at line 309 of file BaseFitObject.cc.

310 {
311 if (!isfinite(par_)) return true;
312 assert(ilocal >= 0 && ilocal < getNPar());
313 if (par[ilocal] == par_) return false;
315 bool result = pow((par_ - par[ilocal]), 2) > eps2 * cov[ilocal][ilocal];
316 par[ilocal] = par_;
317 return result;
318 }

◆ setParam() [2/2]

bool setParam ( int  ilocal,
double  par_,
bool  measured_,
bool  fixed_ = false 
)
virtual

Set value and measured flag of parameter i; return: significant change.

Parameters
ilocalLocal parameter number
par_New parameter value
measured_New "measured" flag
fixed_New "fixed" flag

Definition at line 299 of file BaseFitObject.cc.

301 {
302 assert(ilocal >= 0 && ilocal < getNPar());
303 if (measured[ilocal] != measured_ || fixed[ilocal] != fixed_) invalidateCache();
304 measured[ilocal] = measured_;
305 fixed[ilocal] = fixed_;
306 return setParam(ilocal, par_);
307 }
virtual bool setParam(int ilocal, double par_, bool measured_, bool fixed_=false)
Set value and measured flag of parameter i; return: significant change.

◆ updateParams()

bool updateParams ( double  p[],
int  idim 
)
virtual

Read values from global vector, readjust vector; return: significant change.

Parameters
pThe parameter vector
idimLength of the vector

Reimplemented in ISRPhotonFitObject, JetFitObject, NeutrinoFitObject, and PxPyPzMFitObject.

Definition at line 196 of file BaseFitObject.cc.

197 {
198 bool result = false;
200 for (int ilocal = 0; ilocal < getNPar(); ++ilocal) {
201 if (!isParamFixed(ilocal)) { // daniel added this
202 int iglobal = getGlobalParNum(ilocal);
203 assert(iglobal >= 0 && iglobal < idim);
204 // result = result || setParam (ilocal, p[iglobal]); // daniel thinks this is a BUG !!! if first param is successfully updated, no further ones are
205 // result = setParam (ilocal, p[iglobal]) || result; // daniel thinks this would be OK
206 // but to makes it clearer, does the following:
207 bool thisresult = setParam(ilocal, p[iglobal]);
208 result = result || thisresult;
209 // if illegal value: read back legal value
210 if (!thisresult) // daniel added
211 p[iglobal] = getParam(ilocal);
212 }
213 }
214 return result;
215 }

Friends And Related Function Documentation

◆ operator<<()

std::ostream & operator<< ( std::ostream &  os,
const BaseFitObject bfo 
)
related

Prints out a BaseFitObject, using its print method.

Parameters
osThe output stream
bfoThe object to print

Definition at line 341 of file BaseFitObject.h.

344 {
345 return bfo.print(os);
346 }

Member Data Documentation

◆ cachevalid

bool cachevalid
mutableprotected

flag for valid cache

Definition at line 333 of file BaseFitObject.h.

◆ cov

double cov[BaseDefs::MAXPAR][BaseDefs::MAXPAR]
protected

local covariance matrix

Definition at line 327 of file BaseFitObject.h.

◆ covinv

double covinv[BaseDefs::MAXPAR][BaseDefs::MAXPAR]
mutableprotected

inverse pf local covariance matrix

Definition at line 329 of file BaseFitObject.h.

◆ covinvvalid

bool covinvvalid
mutableprotected

flag for valid inverse covariance matrix

Definition at line 331 of file BaseFitObject.h.

◆ eps2

const double eps2 = 0.0001
staticprotected

Definition at line 307 of file BaseFitObject.h.

◆ fixed

bool fixed[BaseDefs::MAXPAR]
protected

fixed flag

Definition at line 323 of file BaseFitObject.h.

◆ globalParNum

int globalParNum[BaseDefs::MAXPAR]
protected

global parameter number for each parameter

Definition at line 325 of file BaseFitObject.h.

◆ measured

bool measured[BaseDefs::MAXPAR]
protected

measured flag

Definition at line 321 of file BaseFitObject.h.

◆ mpar

double mpar[BaseDefs::MAXPAR]
protected

measured parameters

Definition at line 319 of file BaseFitObject.h.

◆ name

char* name
protected

Definition at line 306 of file BaseFitObject.h.

◆ par

double par[BaseDefs::MAXPAR]
protected

fit parameters

Definition at line 317 of file BaseFitObject.h.


The documentation for this class was generated from the following files: