Belle II Software development
ISRPhotonFitObject Class Reference
Inheritance diagram for ISRPhotonFitObject:
ParticleFitObject BaseFitObject

Public Member Functions

 ISRPhotonFitObject (double px, double py, double pz, double b_, double PzMaxB_, double PzMinB_=0.)
 
 ISRPhotonFitObject (const ISRPhotonFitObject &rhs)
 Copy constructor.
 
ISRPhotonFitObjectoperator= (const ISRPhotonFitObject &rhs)
 Assignment.
 
virtual ISRPhotonFitObjectcopy () const override
 Return a new copy of itself.
 
virtual ISRPhotonFitObjectassign (const BaseFitObject &source) override
 Assign from anther object, if of same type.
 
virtual const char * getParamName (int ilocal) const override
 Get name of parameter ilocal.
 
virtual bool updateParams (double p[], int idim) override
 Read values from global vector, readjust vector; return: significant change.
 
virtual double getDPx (int ilocal) const override
 Return d p_x / d par_ilocal (derivative of px w.r.t. local parameter ilocal)
 
virtual double getDPy (int ilocal) const override
 Return d p_y / d par_ilocal (derivative of py w.r.t. local parameter ilocal)
 
virtual double getDPz (int ilocal) const override
 Return d p_z / d par_ilocal (derivative of pz w.r.t. local parameter ilocal)
 
virtual double getDE (int ilocal) const override
 Return d E / d par_ilocal (derivative of E w.r.t. local parameter ilocal)
 
virtual double getFirstDerivative_Meta_Local (int iMeta, int ilocal, int metaSet) const override
 
virtual double getSecondDerivative_Meta_Local (int iMeta, int ilocal, int jlocal, int metaSet) const override
 
virtual int getNPar () const override
 Get total number of parameters of this FitObject.
 
virtual bool setMass (double mass_)
 Set mass of particle; return=success.
 
virtual double getMass () const
 Get mass of particle.
 
virtual std::ostream & print4Vector (std::ostream &os) const
 print the four-momentum (E, px, py, pz)
 
virtual FourVector getFourMomentum () const
 
virtual double getE () const
 Return E.
 
virtual double getPx () const
 Return px.
 
virtual double getPy () const
 Return py.
 
virtual double getPz () const
 Return pz.
 
virtual double getP () const
 Return p (momentum)
 
virtual double getP2 () const
 Return p (momentum) squared.
 
virtual double getPt () const
 Return pt (transverse momentum)
 
virtual double getPt2 () const
 Return pt (transverse momentum) squared.
 
virtual void getDerivatives (double der[], int idim) const override
 
virtual void addToGlobalChi2DerMatrixNum (double *M, int idim, double eps)
 Add numerically determined derivatives of chi squared to global covariance matrix.
 
virtual void addToGlobalChi2DerVectorNum (double *y, int idim, double eps)
 Add numerically determined derivatives of chi squared to global derivative vector.
 
virtual std::ostream & print (std::ostream &os) const override
 print object to ostream
 
void test1stDerivatives ()
 
void test2ndDerivatives ()
 
double num1stDerivative (int ilocal, double eps)
 Evaluates numerically the 1st derivative of chi2 w.r.t. a parameter.
 
double num2ndDerivative (int ilocal1, double eps1, int ilocal2, double eps2)
 Evaluates numerically the 2nd derivative of chi2 w.r.t. 2 parameters.
 
virtual double getChi2 () const override
 Get chi squared from measured and fitted parameters.
 
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 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 * 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 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 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)
 
void invalidateCache () const
 invalidate any cached quantities
 
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 addToGlobalChi2DerVector (double *y, int idim, double lambda, double der[], int metaSet) const
 Add derivatives of momentum vector to global derivative vector.
 
virtual void addToGlobalChi2DerMatrix (double *M, int idim) const
 Add 2nd derivatives of chi squared to global derivative matrix.
 
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 void initCov ()
 
virtual double getError2 (double der[], int metaset) const
 

Protected Types

enum  { NPAR = 3 }
 

Protected Member Functions

double PgFromPz (double pz)
 
void updateCache () const override
 
virtual bool calculateCovInv () const
 Calculate the inverse of the covariance matrix.
 

Protected Attributes

bool cachevalid
 
double pt2
 
double p2
 
double p
 
double pz
 
double dpx0
 
double dpy0
 
double dpz0
 
double dE0
 
double dpx1
 
double dpy1
 
double dpz1
 
double dE1
 
double dpx2
 
double dpy2
 
double dpz2
 
double dE2
 
double d2pz22
 
double d2E22
 
double chi2
 
double b
 
double PzMinB
 
double PzMaxB
 
double dp2zFact
 
double mass
 mass of particle
 
FourVector fourMomentum
 
double paramCycl [BaseDefs::MAXPAR]
 
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
 

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 29 of file ISRPhotonFitObject.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
protected

Definition at line 77 of file ISRPhotonFitObject.h.

77{NPAR = 3}; // well, it's actually 1...Daniel should update

Constructor & Destructor Documentation

◆ ISRPhotonFitObject() [1/2]

ISRPhotonFitObject ( double  px,
double  py,
double  pz,
double  b_,
double  PzMaxB_,
double  PzMinB_ = 0. 
)
Parameters
pxinitial px value for photon (fixed)
pyinitial py value for photon (fixed)
pzinitial pz value for photon
b_photon spectrum parametrization
PzMaxB_photon spectrum parametrization
PzMinB_photon spectrum parametrization

Definition at line 49 of file ISRPhotonFitObject.cc.

51 : cachevalid(false),
52 pt2(0), p2(0), p(0), pz(0),
53 dpx0(0), dpy0(0), dpz0(0), dE0(0), dpx1(0), dpy1(0), dpz1(0), dE1(0),
54 dpx2(0), dpy2(0), dpz2(0), dE2(0), d2pz22(0), d2E22(0),
55 chi2(0), b(0), PzMinB(0), PzMaxB(0), dp2zFact(0)
56 {
57
58 assert(int(NPAR) <= int(BaseDefs::MAXPAR));
59
60 initCov();
61 b = b_;
62 PzMinB = PzMinB_;
63 PzMaxB = PzMaxB_;
64#ifdef DEBUG
65 B2INFO("ISRPhotonFitObject: b: " << b << " PzMinB: " << PzMinB << " PzMaxB: " << PzMaxB);
66#endif
67
68 if (b <= 0. || b >= 1.) {
69 B2INFO("ISRPhotonFitObject: b must be from ]0,1[ ");
70 }
71 assert(b > 0. && b < 1.);
72 if (PzMinB < 0. || PzMaxB <= PzMinB) {
73 B2INFO("ISRPhotonFitObject: PzMinB and PzMaxB must be chosen such that 0 <= PzMinB < PzMaxB");
74 }
75 assert(PzMinB >= 0.);
76 assert(PzMaxB > PzMinB);
77 dp2zFact = (PzMaxB - PzMinB) / b * sqrt(2. / pi_);
78 double pg = PgFromPz(ppz); // using internally Gauss-distributed parameter p_g instead of p_z
79 setParam(0, px, true, true);
80 setParam(1, py, true, true);
81 setParam(2, pg, true);
82 setMParam(0, 0.); // all measured parameters
83 setMParam(1, 0.); // are assumed to be zero
84 setMParam(2, 0.); // in this photon parametrization
85#ifdef DEBUG
86 B2INFO("ISRPhotonFitObject: Initial pg: " << pg);
87#endif
88 setError(2, 1.);
89 setMass(0.);
91 }
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 setMParam(int ilocal, double mpar_)
Set measured value of parameter ilocal; return: success.
void invalidateCache() const
invalidate any cached quantities
virtual bool setError(int ilocal, double err_)
Set error of parameter ilocal; return: success.
virtual bool setMass(double mass_)
Set mass of particle; return=success.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ ISRPhotonFitObject() [2/2]

Copy constructor.

Parameters
rhsright hand side

Definition at line 97 of file ISRPhotonFitObject.cc.

98 : ParticleFitObject(rhs), cachevalid(false),
99 pt2(0), p2(0), p(0), pz(0),
100 dpx0(0), dpy0(0), dpz0(0), dE0(0), dpx1(0), dpy1(0), dpz1(0), dE1(0),
101 dpx2(0), dpy2(0), dpz2(0), dE2(0), d2pz22(0), d2E22(0),
102 chi2(0), b(0), PzMinB(0), PzMaxB(0), dp2zFact(0)
103 {
105 }
virtual ISRPhotonFitObject & assign(const BaseFitObject &source) override
Assign from anther object, if of same type.

Member Function Documentation

◆ addTo1stDerivatives()

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

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.
bool cachevalid
flag for valid cache
int globalParNum[BaseDefs::MAXPAR]
global parameter number for each parameter

◆ addTo2ndDerivatives() [1/2]

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

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
virtualinherited

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
virtualinherited

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 }
bool covinvvalid
flag for valid inverse covariance matrix
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.

◆ addToGlobalChi2DerMatrixNum()

void addToGlobalChi2DerMatrixNum ( double *  M,
int  idim,
double  eps 
)
virtualinherited

Add numerically determined derivatives of chi squared to global covariance matrix.

Parameters
MGlobal covariance matrix
idimFirst dimension of global covariance matrix
epsParameter variation

Definition at line 162 of file ParticleFitObject.cc.

163 {
164 for (int ilocal1 = 0; ilocal1 < getNPar(); ++ilocal1) {
165 int iglobal1 = getGlobalParNum(ilocal1);
166 for (int ilocal2 = ilocal1; ilocal2 < getNPar(); ++ilocal2) {
167 int iglobal2 = getGlobalParNum(ilocal2);
168 M[idim * iglobal1 + iglobal2] += num2ndDerivative(ilocal1, eps, ilocal2, eps);
169 }
170 }
171 }
double num2ndDerivative(int ilocal1, double eps1, int ilocal2, double eps2)
Evaluates numerically the 2nd derivative of chi2 w.r.t. 2 parameters.

◆ addToGlobalChi2DerVector() [1/2]

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

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
virtualinherited

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 }

◆ addToGlobalChi2DerVectorNum()

void addToGlobalChi2DerVectorNum ( double *  y,
int  idim,
double  eps 
)
virtualinherited

Add numerically determined derivatives of chi squared to global derivative vector.

Parameters
yVector of chi2 derivatives
idimVector size
epsParameter variation

Definition at line 152 of file ParticleFitObject.cc.

153 {
154 for (int ilocal = 0; ilocal < getNPar(); ++ilocal) {
155 int iglobal = getGlobalParNum(ilocal);
156 assert(iglobal >= 0 && iglobal < idim);
157 y[iglobal] += num1stDerivative(ilocal, eps);
158 }
159 }
double num1stDerivative(int ilocal, double eps)
Evaluates numerically the 1st derivative of chi2 w.r.t. a parameter.

◆ addToGlobCov()

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

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()

ISRPhotonFitObject & assign ( const BaseFitObject source)
overridevirtual

Assign from anther object, if of same type.

Parameters
sourceThe source object

Reimplemented from ParticleFitObject.

Definition at line 120 of file ISRPhotonFitObject.cc.

121 {
122 if (const auto* psource = dynamic_cast<const ISRPhotonFitObject*>(&source)) {
123 if (psource != this) {
125 // only mutable data members, need not to be copied, if cache is invalid
126 }
127 } else {
128 assert(0);
129 }
130 return *this;
131 }
ISRPhotonFitObject(double px, double py, double pz, double b_, double PzMaxB_, double PzMinB_=0.)
virtual ParticleFitObject & assign(const BaseFitObject &source) override
Assign from anther object, if of same type.

◆ calculateCovInv()

bool calculateCovInv ( ) const
protectedvirtualinherited

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 covinv[BaseDefs::MAXPAR][BaseDefs::MAXPAR]
inverse pf local covariance matrix
double cov[BaseDefs::MAXPAR][BaseDefs::MAXPAR]
local covariance matrix

◆ copy()

ISRPhotonFitObject * copy ( ) const
overridevirtual

Return a new copy of itself.

Implements BaseFitObject.

Definition at line 115 of file ISRPhotonFitObject.cc.

116 {
117 return new ISRPhotonFitObject(*this);
118 }

◆ fixParam()

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

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 }
bool fixed[BaseDefs::MAXPAR]
fixed flag

◆ getChi2()

double getChi2 ( ) const
overridevirtualinherited

Get chi squared from measured and fitted parameters.

Reimplemented from BaseFitObject.

Definition at line 276 of file ParticleFitObject.cc.

277 {
278 // reimplemented here to take account of cyclical variables e.g azimuthal angle phi - DJeans
279
280 if (not covinvvalid and not calculateCovInv()) return 0;
281
282 double resid[BaseDefs::MAXPAR] = {0};
283 for (int i = 0; i < getNPar(); i++) {
284 if (isParamMeasured(i) && !isParamFixed(i)) {
285 resid[i] = par[i] - mpar[i];
286 if (paramCycl[i] > 0) {
287 resid[i] = fmod(resid[i], paramCycl[i]);
288 if (resid[i] > paramCycl[i] / 2) resid[i] -= paramCycl[i];
289 if (resid[i] < -paramCycl[i] / 2) resid[i] += paramCycl[i];
290 }
291 }
292 }
293
294
295 double chi2 = 0;
296 for (int i = 0; i < getNPar(); i++) {
297 if (isParamMeasured(i) && !isParamFixed(i)) {
298 for (int j = 0; j < getNPar(); j++) {
299 if (isParamMeasured(j) && !isParamFixed(j)) {
300 chi2 += resid[i] * covinv[i][j] * resid[j];
301 // B2INFO( getName () << " === " << i << " " << j << " : " <<
302 // resid[i] << "*" << covinv[i][j] << "*" << resid[j] << " = " << resid[i]*covinv[i][j]*resid[j] << " , sum " << chi2);
303 }
304 }
305 }
306 }
307
308 // B2INFO("ParticleFitObject::getChi2 () chi2 = " << chi2 );
309 return chi2;
310 }
double par[BaseDefs::MAXPAR]
fit parameters
double mpar[BaseDefs::MAXPAR]
measured parameters

◆ getCov()

double getCov ( int  ilocal,
int  jlocal 
) const
virtualinherited

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
virtualinherited

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
virtualinherited

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 }

◆ getDE()

double getDE ( int  ilocal) const
overridevirtual

Return d E / d par_ilocal (derivative of E w.r.t. local parameter ilocal)

Parameters
ilocalLocal parameter number

Implements ParticleFitObject.

Definition at line 194 of file ISRPhotonFitObject.cc.

195 {
196 assert(ilocal >= 0 && ilocal < NPAR);
197 if (!cachevalid) updateCache();
198 switch (ilocal) {
199 case 0: return dE0;
200 case 1: return dE1;
201 case 2: return dE2;
202 }
203 return 0;
204 }

◆ getDerivatives()

void getDerivatives ( double  der[],
int  idim 
) const
overridevirtualinherited

Implements BaseFitObject.

Definition at line 173 of file ParticleFitObject.cc.

174 {
175 for (int ilocal = 0; ilocal < getNPar(); ++ilocal) {
176 assert(ilocal < idim);
177 der [4 * ilocal] = getDE(ilocal);
178 der [4 * ilocal + 1] = getDPx(ilocal);
179 der [4 * ilocal + 2] = getDPy(ilocal);
180 der [4 * ilocal + 3] = getDPz(ilocal);
181 }
182 }
virtual double getDE(int ilocal) const =0
Return d E / d par_ilocal (derivative of E w.r.t. local parameter ilocal)
virtual double getDPx(int ilocal) const =0
Return d p_x / d par_ilocal (derivative of px w.r.t. local parameter ilocal)
virtual double getDPz(int ilocal) const =0
Return d p_z / d par_ilocal (derivative of pz w.r.t. local parameter ilocal)
virtual double getDPy(int ilocal) const =0
Return d p_y / d par_ilocal (derivative of py w.r.t. local parameter ilocal)

◆ getDPx()

double getDPx ( int  ilocal) const
overridevirtual

Return d p_x / d par_ilocal (derivative of px w.r.t. local parameter ilocal)

Parameters
ilocalLocal parameter number

Implements ParticleFitObject.

Definition at line 158 of file ISRPhotonFitObject.cc.

159 {
160 assert(ilocal >= 0 && ilocal < NPAR);
161 if (!cachevalid) updateCache();
162 switch (ilocal) {
163 case 0: return dpx0;
164 case 1: return dpx1;
165 case 2: return dpx2;
166 }
167 return 0;
168 }

◆ getDPy()

double getDPy ( int  ilocal) const
overridevirtual

Return d p_y / d par_ilocal (derivative of py w.r.t. local parameter ilocal)

Parameters
ilocalLocal parameter number

Implements ParticleFitObject.

Definition at line 170 of file ISRPhotonFitObject.cc.

171 {
172 assert(ilocal >= 0 && ilocal < NPAR);
173 if (!cachevalid) updateCache();
174 switch (ilocal) {
175 case 0: return dpy0;
176 case 1: return dpy1;
177 case 2: return dpy2;
178 }
179 return 0;
180 }

◆ getDPz()

double getDPz ( int  ilocal) const
overridevirtual

Return d p_z / d par_ilocal (derivative of pz w.r.t. local parameter ilocal)

Parameters
ilocalLocal parameter number

Implements ParticleFitObject.

Definition at line 182 of file ISRPhotonFitObject.cc.

183 {
184 assert(ilocal >= 0 && ilocal < NPAR);
185 if (!cachevalid) updateCache();
186 switch (ilocal) {
187 case 0: return dpz0;
188 case 1: return dpz1;
189 case 2: return dpz2;
190 }
191 return 0;
192 }

◆ getE()

double getE ( ) const
virtualinherited

Return E.

Definition at line 106 of file ParticleFitObject.cc.

107 {
108 return getFourMomentum().getE();
109 }
double getE() const
Returns the energy / 0 component.
Definition: FourVector.h:127

◆ getError()

double getError ( int  ilocal) const
virtualinherited

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
virtualinherited

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()

double getFirstDerivative_Meta_Local ( int  iMeta,
int  ilocal,
int  metaSet 
) const
overridevirtual

Implements BaseFitObject.

Definition at line 207 of file ISRPhotonFitObject.cc.

208 {
209 assert(metaSet == 0);
210 switch (iMeta) {
211 case 0:
212 return getDE(ilocal);
213 break;
214 case 1:
215 return getDPx(ilocal);
216 break;
217 case 2:
218 return getDPy(ilocal);
219 break;
220 case 3:
221 return getDPz(ilocal);
222 break;
223 default:
224 assert(0);
225 }
226 return -999;
227 }
virtual double getDPx(int ilocal) const override
Return d p_x / d par_ilocal (derivative of px w.r.t. local parameter ilocal)
virtual double getDPy(int ilocal) const override
Return d p_y / d par_ilocal (derivative of py w.r.t. local parameter ilocal)
virtual double getDE(int ilocal) const override
Return d E / d par_ilocal (derivative of E w.r.t. local parameter ilocal)
virtual double getDPz(int ilocal) const override
Return d p_z / d par_ilocal (derivative of pz w.r.t. local parameter ilocal)

◆ getFourMomentum()

FourVector getFourMomentum ( ) const
virtualinherited

Definition at line 101 of file ParticleFitObject.cc.

102 {
103 if (!cachevalid) updateCache();
104 return fourMomentum;
105 }

◆ getGlobalParNum()

int getGlobalParNum ( int  ilocal) const
virtualinherited

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 }

◆ getMass()

double getMass ( ) const
virtualinherited

Get mass of particle.

Definition at line 89 of file ParticleFitObject.cc.

90 {
91 return mass;
92 }

◆ getMParam()

double getMParam ( int  ilocal) const
virtualinherited

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
virtualinherited

Get object's name.

Definition at line 107 of file BaseFitObject.cc.

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

◆ getNFixed()

int getNFixed ( ) const
virtualinherited

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
virtualinherited

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
virtualinherited

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
inlineoverridevirtual

Get total number of parameters of this FitObject.

Implements BaseFitObject.

Definition at line 73 of file ISRPhotonFitObject.h.

73{return NPAR;}

◆ getNUnmeasured()

int getNUnmeasured ( ) const
virtualinherited

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 }

◆ getP()

double getP ( ) const
virtualinherited

Return p (momentum)

Definition at line 122 of file ParticleFitObject.cc.

123 {
124 return getFourMomentum().getP();
125 }
double getP() const
Returns the momentum / magnitude of the three vector.
Definition: FourVector.h:136

◆ getP2()

double getP2 ( ) const
virtualinherited

Return p (momentum) squared.

Definition at line 126 of file ParticleFitObject.cc.

127 {
128 return getFourMomentum().getP2();
129 }
double getP2() const
Returns the momentum squared / magnitude of the three vector squared.
Definition: FourVector.h:135

◆ getParam()

double getParam ( int  ilocal) const
virtualinherited

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()

const char * getParamName ( int  ilocal) const
overridevirtual

Get name of parameter ilocal.

Parameters
ilocalLocal parameter number

Implements BaseFitObject.

Definition at line 133 of file ISRPhotonFitObject.cc.

134 {
135 switch (ilocal) {
136 case 0: return "P_x";
137 case 1: return "P_y";
138 case 2: return "P_g";
139 }
140 return "undefined";
141 }

◆ getPt()

double getPt ( ) const
virtualinherited

Return pt (transverse momentum)

Definition at line 130 of file ParticleFitObject.cc.

131 {
132 return getFourMomentum().getPt();
133 }
double getPt() const
Returns the transverse momentum / magnitude of the 1 and 2 component vector.
Definition: FourVector.h:133

◆ getPt2()

double getPt2 ( ) const
virtualinherited

Return pt (transverse momentum) squared.

Definition at line 134 of file ParticleFitObject.cc.

135 {
136 return getFourMomentum().getPt2();
137 }
double getPt2() const
Returns the transverse momentum squared / magnitude of the 1 and 2 component vector squared.
Definition: FourVector.h:132

◆ getPx()

double getPx ( ) const
virtualinherited

Return px.

Definition at line 110 of file ParticleFitObject.cc.

111 {
112 return getFourMomentum().getPx();
113 }
double getPx() const
Returns the x momentum / 1 component.
Definition: FourVector.h:128

◆ getPy()

double getPy ( ) const
virtualinherited

Return py.

Definition at line 114 of file ParticleFitObject.cc.

115 {
116 return getFourMomentum().getPy();
117 }
double getPy() const
Returns the y momentum / 2 component.
Definition: FourVector.h:129

◆ getPz()

double getPz ( ) const
virtualinherited

Return pz.

Definition at line 118 of file ParticleFitObject.cc.

119 {
120 return getFourMomentum().getPz();
121 }
double getPz() const
Returns the z momentum / 3 component.
Definition: FourVector.h:130

◆ getRho()

double getRho ( int  ilocal,
int  jlocal 
) const
virtualinherited

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 }

◆ getSecondDerivative_Meta_Local()

double getSecondDerivative_Meta_Local ( int  iMeta,
int  ilocal,
int  jlocal,
int  metaSet 
) const
overridevirtual

Implements BaseFitObject.

Definition at line 229 of file ISRPhotonFitObject.cc.

230 {
231 assert(metaSet == 0);
232 if (!cachevalid) updateCache();
233
234 if (jlocal < ilocal) {
235 int temp = jlocal;
236 jlocal = ilocal;
237 ilocal = temp;
238 }
239
240 // daniel hasn't checked these, copied from orig code
241 switch (iMeta) {
242
243 case 0:
244 if (ilocal == 2 && jlocal == 2) return d2E22;
245 else return 0;
246 break;
247 case 1:
248 case 2:
249 return 0;
250 break;
251 case 3:
252 if (ilocal == 2 && jlocal == 2) return d2pz22;
253 else return 0;
254 break;
255 default:
256 return 0;
257 }
258
259 }

◆ initCov()

void initCov ( )
virtualinherited

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
inlineinherited

invalidate any cached quantities

Definition at line 241 of file BaseFitObject.h.

241{cachevalid = false;};

◆ isParamFixed()

bool isParamFixed ( int  ilocal) const
virtualinherited

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
virtualinherited

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 }
bool measured[BaseDefs::MAXPAR]
measured flag

◆ num1stDerivative()

double num1stDerivative ( int  ilocal,
double  eps 
)
inherited

Evaluates numerically the 1st derivative of chi2 w.r.t. a parameter.

Parameters
ilocalLocal parameter number
epsvariation of local parameter

Definition at line 230 of file ParticleFitObject.cc.

231 {
232 double save = getParam(ilocal);
233 setParam(ilocal, save + eps);
234 double v1 = getChi2();
235 setParam(ilocal, save - eps);
236 double v2 = getChi2();
237 double result = (v1 - v2) / (2 * eps);
238 setParam(ilocal, save);
239 return result;
240 }
virtual double getParam(int ilocal) const
Get current value of parameter ilocal.
virtual double getChi2() const override
Get chi squared from measured and fitted parameters.
const std::vector< double > v2
MATLAB generated random vector.
const std::vector< double > v1
MATLAB generated random vector.

◆ num2ndDerivative()

double num2ndDerivative ( int  ilocal1,
double  eps1,
int  ilocal2,
double  eps2 
)
inherited

Evaluates numerically the 2nd derivative of chi2 w.r.t. 2 parameters.

Parameters
ilocal11st local parameter number
eps1variation of 1st local parameter
ilocal21st local parameter number
eps2variation of 2nd local parameter

Definition at line 242 of file ParticleFitObject.cc.

244 {
245 double result;
246
247 if (ilocal1 == ilocal2) {
248 double save = getParam(ilocal1);
249 double v0 = getChi2();
250 setParam(ilocal1, save + eeps1);
251 double v1 = getChi2();
252 setParam(ilocal1, save - eeps1);
253 double v2 = getChi2();
254 result = (v1 + v2 - 2 * v0) / (eeps1 * eeps1);
255 setParam(ilocal1, save);
256 } else {
257 double save1 = getParam(ilocal1);
258 double save2 = getParam(ilocal2);
259 setParam(ilocal1, save1 + eeps1);
260 setParam(ilocal2, save2 + eeps2);
261 double v11 = getChi2();
262 setParam(ilocal2, save2 - eeps2);
263 double v12 = getChi2();
264 setParam(ilocal1, save1 - eeps1);
265 double v22 = getChi2();
266 setParam(ilocal2, save2 + eeps2);
267 double v21 = getChi2();
268 result = (v11 + v22 - v12 - v21) / (4 * eeps1 * eeps2);
269 setParam(ilocal1, save1);
270 setParam(ilocal2, save2);
271 }
272 return result;
273 }

◆ operator=()

ISRPhotonFitObject & operator= ( const ISRPhotonFitObject rhs)

Assignment.

right hand side

Definition at line 107 of file ISRPhotonFitObject.cc.

108 {
109 if (this != &rhs) {
110 assign(rhs); // calls virtual function assign of derived class
111 }
112 return *this;
113 }

◆ PgFromPz()

double PgFromPz ( double  pz)
protected

Definition at line 263 of file ISRPhotonFitObject.cc.

264 {
265
266 int sign = (ppz > 0.) - (ppz < 0.);
267 double u = (pow(fabs(ppz), b) - PzMinB) / (PzMaxB - PzMinB);
268
269 if (u < 0.) {
270 B2INFO("ISRPhotonFitObject: Initial pz with abs(pz) < pzMin adjusted to zero.");
271 u = 0.;
272 }
273
274 if (u >= 1.) {
275 B2INFO("ISRPhotonFitObject: Initial pz with abs(pz) >= pzMax adjusted.");
276 u = 0.99999999;
277 }
278
279 double g = std::log(1. - u * u);
280 double g4pa = g + 4. / pi_ / a;
281 return sign * sqrt(-g4pa + sqrt(g4pa * g4pa - 4. / a * g)) ;
282 }

◆ print()

std::ostream & print ( std::ostream &  os) const
overridevirtualinherited

print object to ostream

Parameters
osThe output stream

Implements BaseFitObject.

Definition at line 140 of file ParticleFitObject.cc.

141 {
142
143 if (!cachevalid) updateCache();
144
145 printParams(os);
146 os << " => ";
147 print4Vector(os);
148 os << std::endl;
149 return os;
150 }
virtual std::ostream & printParams(std::ostream &os) const
print the parameters and errors
virtual std::ostream & print4Vector(std::ostream &os) const
print the four-momentum (E, px, py, pz)

◆ print1stDerivatives()

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

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
virtualinherited

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 }

◆ print4Vector()

std::ostream & print4Vector ( std::ostream &  os) const
virtualinherited

print the four-momentum (E, px, py, pz)

Parameters
osThe output stream

Definition at line 94 of file ParticleFitObject.cc.

95 {
96 os << "[" << getE() << ", " << getPx() << ", "
97 << getPy() << ", " << getPz() << "]";
98 return os;
99 }
virtual double getPx() const
Return px.
virtual double getPz() const
Return pz.
virtual double getPy() const
Return py.
virtual double getE() const
Return E.

◆ printParams()

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

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.

◆ printRhoValues()

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

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)
inlinevirtualinherited

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_ 
)
virtualinherited

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_ 
)
virtualinherited

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 
)
virtualinherited

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 }

◆ setMass()

bool setMass ( double  mass_)
virtualinherited

Set mass of particle; return=success.

Definition at line 63 of file ParticleFitObject.cc.

64 {
65 if (!isfinite(mass_)) return false;
66 if (mass == mass_) return true;
68 mass = std::abs(mass_);
69 return true;
70 }

◆ setMParam()

bool setMParam ( int  ilocal,
double  mpar_ 
)
virtualinherited

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_)
inherited

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_ 
)
virtualinherited

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 
)
virtualinherited

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 }

◆ test1stDerivatives()

void test1stDerivatives ( )
inherited

Definition at line 184 of file ParticleFitObject.cc.

185 {
186 B2INFO("ParticleFitObject::test1stDerivatives, object " << getName() << "\n");
187 double ycalc[100], ynum[100];
188 for (int i = 0; i < 100; ++i) ycalc[i] = ynum[i] = 0;
189 addToGlobalChi2DerVector(ycalc, 100);
190 double eps = 0.00001;
191 addToGlobalChi2DerVectorNum(ynum, 100, eps);
192 for (int ilocal = 0; ilocal < getNPar(); ++ilocal) {
193 int iglobal = getGlobalParNum(ilocal);
194 double calc = ycalc[iglobal];
195 double num = ynum[iglobal];
196 B2INFO("fo: " << getName() << " par " << ilocal << "/"
197 << iglobal << " (" << getParamName(ilocal)
198 << ") calc: " << calc << " - num: " << num << " = " << calc - num);
199 }
200 }
virtual const char * getParamName(int ilocal) const =0
Get name of parameter ilocal.
virtual const char * getName() const
Get object's name.
virtual int addToGlobalChi2DerVector(double *y, int idim) const
Add derivatives of chi squared to global derivative vector.
virtual void addToGlobalChi2DerVectorNum(double *y, int idim, double eps)
Add numerically determined derivatives of chi squared to global derivative vector.

◆ test2ndDerivatives()

void test2ndDerivatives ( )
inherited

Definition at line 202 of file ParticleFitObject.cc.

203 {
204 B2INFO("ParticleFitObject::test2ndDerivatives, object " << getName() << "\n");
205 const int idim = 100;
206 auto* Mnum = new double[idim * idim];
207 auto* Mcalc = new double[idim * idim];
208 for (int i = 0; i < idim * idim; ++i) Mnum[i] = Mcalc[i] = 0;
209 addToGlobalChi2DerMatrix(Mcalc, idim);
210 double eps = 0.0001;
211 B2INFO("eps=" << eps);
212 addToGlobalChi2DerMatrixNum(Mnum, idim, eps);
213 for (int ilocal1 = 0; ilocal1 < getNPar(); ++ilocal1) {
214 int iglobal1 = getGlobalParNum(ilocal1);
215 for (int ilocal2 = ilocal1; ilocal2 < getNPar(); ++ilocal2) {
216 int iglobal2 = getGlobalParNum(ilocal2);
217 double calc = Mcalc[idim * iglobal1 + iglobal2];
218 double num = Mnum[idim * iglobal1 + iglobal2];
219 B2INFO("fo: " << getName() << " par " << ilocal1 << "/"
220 << iglobal1 << " (" << getParamName(ilocal1)
221 << "), par " << ilocal2 << "/"
222 << iglobal2 << " (" << getParamName(ilocal2)
223 << ") calc: " << calc << " - num: " << num << " = " << calc - num);
224 }
225 }
226 delete[] Mnum;
227 delete[] Mcalc;
228 }
virtual void addToGlobalChi2DerMatrix(double *M, int idim) const
Add 2nd derivatives of chi squared to global derivative matrix.
virtual void addToGlobalChi2DerMatrixNum(double *M, int idim, double eps)
Add numerically determined derivatives of chi squared to global covariance matrix.

◆ updateCache()

void updateCache ( ) const
overrideprotectedvirtual

Implements BaseFitObject.

Definition at line 285 of file ISRPhotonFitObject.cc.

286 {
287 double px = par[0];
288 double py = par[1];
289 double pg = par[2];
290
291 int sign = (pg > 0.) - (pg < 0.);
292 double pg2h = pg * pg / 2.;
293 double exponent = -pg2h * (4. / pi_ + a * pg2h) / (1. + a * pg2h);
294 double u = sqrt((exponent < -1.e-14) ? 1. - exp(exponent) : -exponent); // approximation to avoid numerical problem
295 pz = sign * pow((PzMinB + (PzMaxB - PzMinB) * u), (1. / b));
296
297 pt2 = px * px + py * py;
298 p2 = pt2 + pz * pz;
299 p = std::sqrt(p2);
300
301 fourMomentum.setValues(p, px, py, pz);
302
303 dpx0 = 1.;
304 dpx1 = 0.;
305 dpx2 = 0.;
306 dpy0 = 0.;
307 dpy1 = 1.;
308 dpy2 = 0.;
309 dpz0 = 0.;
310 dpz1 = 0.;
311 dpz2 = dp2zFact * pow(fabs(pz), (1. - b)) * exp(-pg * pg / 2.);
312
313 // if p,pz==0, derivatives are zero (catch up 1/0)
314 if (pz) {
315 d2pz22 = dpz2 * ((1. - b) * dpz2 / pz - par[2]);
316 } else {
317 d2pz22 = 0.;
318 }
319 if (p) {
320 dE0 = px / p;
321 dE1 = py / p;
322 dE2 = pz / p * dpz2;
323 d2E22 = pz / p * d2pz22;
324 if (pt2) {
325 d2E22 += pt2 / p / p / p * dpz2 * dpz2; // NOT using /p2/p to avoid numerical problems
326 }
327 } else {
328 dE0 = 0.;
329 dE1 = 0.;
330 dE2 = 0.;
331 d2E22 = 0.;
332 }
333
334#ifdef DEBUG
335 B2INFO("ISRPhotonFitObject::updateCache: pg: " << pg << " pz: " << pz << " p: " << p << " p^2: " << p2 << "\n"
336 << " Dpz/Dpg: " << dpz2 << " DE/Dpg: " << dE2 << " D^2pz/Dpg^2: " << d2pz22
337 << " D^2E/Dpg^2: " << d2E22);
338#endif
339
340 cachevalid = true;
341 }

◆ updateParams()

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

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

Parameters
pThe parameter vector
idimLength of the vector

Reimplemented from BaseFitObject.

Definition at line 143 of file ISRPhotonFitObject.cc.

144 {
146 int i2 = getGlobalParNum(2);
147 assert(i2 >= 0 && i2 < idim);
148 double pp2 = pp[i2];
149#ifdef DEBUG
150 B2INFO("ISRPhotonFitObject::updateParams: p2(new) = " << pp[i2] << " par[2](old) = " << par[2]);
151#endif
152 bool result = ((pp2 - par[2]) * (pp2 - par[2]) > eps2 * cov[2][2]);
153 par[2] = pp2;
154 pp[i2] = par[2];
155 return result;
156 }

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

◆ b

double b
protected

Definition at line 89 of file ISRPhotonFitObject.h.

◆ cachevalid

bool cachevalid
mutableprotected

Definition at line 83 of file ISRPhotonFitObject.h.

◆ chi2

double chi2
protected

Definition at line 88 of file ISRPhotonFitObject.h.

◆ cov

double cov[BaseDefs::MAXPAR][BaseDefs::MAXPAR]
protectedinherited

local covariance matrix

Definition at line 327 of file BaseFitObject.h.

◆ covinv

double covinv[BaseDefs::MAXPAR][BaseDefs::MAXPAR]
mutableprotectedinherited

inverse pf local covariance matrix

Definition at line 329 of file BaseFitObject.h.

◆ covinvvalid

bool covinvvalid
mutableprotectedinherited

flag for valid inverse covariance matrix

Definition at line 331 of file BaseFitObject.h.

◆ d2E22

double d2E22
protected

Definition at line 87 of file ISRPhotonFitObject.h.

◆ d2pz22

double d2pz22
protected

Definition at line 87 of file ISRPhotonFitObject.h.

◆ dE0

double dE0
protected

Definition at line 86 of file ISRPhotonFitObject.h.

◆ dE1

double dE1
protected

Definition at line 86 of file ISRPhotonFitObject.h.

◆ dE2

double dE2
protected

Definition at line 87 of file ISRPhotonFitObject.h.

◆ dp2zFact

double dp2zFact
protected

Definition at line 89 of file ISRPhotonFitObject.h.

◆ dpx0

double dpx0
protected

Definition at line 86 of file ISRPhotonFitObject.h.

◆ dpx1

double dpx1
protected

Definition at line 86 of file ISRPhotonFitObject.h.

◆ dpx2

double dpx2
protected

Definition at line 87 of file ISRPhotonFitObject.h.

◆ dpy0

double dpy0
protected

Definition at line 86 of file ISRPhotonFitObject.h.

◆ dpy1

double dpy1
protected

Definition at line 86 of file ISRPhotonFitObject.h.

◆ dpy2

double dpy2
protected

Definition at line 87 of file ISRPhotonFitObject.h.

◆ dpz0

double dpz0
protected

Definition at line 86 of file ISRPhotonFitObject.h.

◆ dpz1

double dpz1
protected

Definition at line 86 of file ISRPhotonFitObject.h.

◆ dpz2

double dpz2
protected

Definition at line 87 of file ISRPhotonFitObject.h.

◆ eps2

const double eps2 = 0.0001
staticprotectedinherited

Definition at line 307 of file BaseFitObject.h.

◆ fixed

bool fixed[BaseDefs::MAXPAR]
protectedinherited

fixed flag

Definition at line 323 of file BaseFitObject.h.

◆ fourMomentum

FourVector fourMomentum
mutableprotectedinherited

Definition at line 181 of file ParticleFitObject.h.

◆ globalParNum

int globalParNum[BaseDefs::MAXPAR]
protectedinherited

global parameter number for each parameter

Definition at line 325 of file BaseFitObject.h.

◆ mass

double mass
protectedinherited

mass of particle

Definition at line 179 of file ParticleFitObject.h.

◆ measured

bool measured[BaseDefs::MAXPAR]
protectedinherited

measured flag

Definition at line 321 of file BaseFitObject.h.

◆ mpar

double mpar[BaseDefs::MAXPAR]
protectedinherited

measured parameters

Definition at line 319 of file BaseFitObject.h.

◆ name

char* name
protectedinherited

Definition at line 306 of file BaseFitObject.h.

◆ p

double p
protected

Definition at line 85 of file ISRPhotonFitObject.h.

◆ p2

double p2
protected

Definition at line 85 of file ISRPhotonFitObject.h.

◆ par

double par[BaseDefs::MAXPAR]
protectedinherited

fit parameters

Definition at line 317 of file BaseFitObject.h.

◆ paramCycl

double paramCycl[BaseDefs::MAXPAR]
protectedinherited

Definition at line 184 of file ParticleFitObject.h.

◆ pt2

double pt2
mutableprotected

Definition at line 85 of file ISRPhotonFitObject.h.

◆ pz

double pz
protected

Definition at line 85 of file ISRPhotonFitObject.h.

◆ PzMaxB

double PzMaxB
protected

Definition at line 89 of file ISRPhotonFitObject.h.

◆ PzMinB

double PzMinB
protected

Definition at line 89 of file ISRPhotonFitObject.h.


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