ATLAS Offline Software
Loading...
Searching...
No Matches
ConstraintFit.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
6//
7// Z mass constraint fit
8//
9// original version by:
10// Konstantinos Nikolopoulos
11// Univ. of Athens
12//
14#ifndef ZMASSCONSTRAINT_CONSTRAINTFIT_H
15#define ZMASSCONSTRAINT_CONSTRAINTFIT_H
16
17// Framework include(s):
18#include "AsgTools/AsgTool.h"
19
21
23#include "xAODEgamma/Photon.h"
24#include "AsgTools/ToolHandle.h"
27
28namespace ZMassConstraint
29{
30 class ConstraintFit : public virtual IConstraintFit,
31 public asg::AsgTool
32 {
33
36
37 public:
38
40 ConstraintFit( const std::string& name );
41
42 ~ConstraintFit(void); // default destructor
43
45 StatusCode initialize();
46
48 StatusCode doMassFit (const ConstraintFitInput& input,
49 ConstraintFitOutput& output);
50
52 double getMassError (const ConstraintFitInput& firstInput,
53 const ConstraintFitInput& secondInput = ConstraintFitInput());
54
57 double getMassError (const ConstraintFitOutput& fitOutput,
58 const ConstraintFitInput& extraInput = ConstraintFitInput());
59
61 double getMassError (const ConstraintFitOutput& fitOutput,
62 const ConstraintFitOutput& secondFitOutput);
63
65 void addParticle (const xAOD::Muon& part,
66 ConstraintFitInput& input,
68
70 void addParticle (const xAOD::Electron& part,
71 float elEnergyRes,
72 ConstraintFitInput& input);
73
75 void addFSRParticle(const xAOD::IParticle & part,
76 const TLorentzVector& fsr4Vec,
77 ConstraintFitInput& input);
78
79 private:
80
81 // Enum list for tracking covariance
83 d0 = 0,
84 z0 = 1,
85 phi0 = 2,
86 theta = 3,
87 qOverP = 4,
88 qP = 4,
89 x = 0,
90 y = 1,
91 z = 2
92 };
93
94 int massFitInterface(const ConstraintFitInput& theInput);
95 void massFitRun(ConstraintFitOutput& output, double zresol = -1.);
96 double getMassError(const std::vector<TLorentzVector>& particles,
97 const std::vector<AmgMatrix(3,3)>& covariances);
98 double likelihoodMass2(void);
99 double likelihoodMass(double);
100 void setIgnore(bool val){ m_ignoreInputChecks=val; }
101 bool doSanityChecksOnCovariance(const TLorentzVector& vector, const AmgMatrix(5,5)& covar) const;
102 void convertCovd0z0PhiThetaPToXYZ(const TLorentzVector& fourVec,
103 const AmgMatrix(5,5)& covard0z0PhiThetaP,
104 AmgMatrix(5,5)& covarXYZ) const;
105 void convertCovXYZTod0z0PhiThetaP(const std::vector<TLorentzVector>& particleList,
106 const Amg::MatrixX& covarXYZ,
107 Amg::MatrixX& covard0z0PhiThetaP) const;
108 private:
110 double m_conMass;
114
115 // Calibration tools
116 ToolHandle<CP::IEgammaCalibrationAndSmearingTool> m_energyRescaler; // electron resolution
117 ToolHandle<CP::IMuonCalibrationAndSmearingTool> m_mu_resolSFTool; // Muon error rescaling
118
120 unsigned int m_parameters;
121 unsigned int m_nobj;
122 std::vector<double> m_objmass;
127
128 // double calculateChi2(const Amg::MatrixX* p0, Amg::MatrixX* var, const double Mass);
129 double massFitCalculation(const Amg::MatrixX& var, double mass, Amg::MatrixX& p0);
130 double massFit(const Amg::MatrixX& /*p0*/, const Amg::MatrixX& var, double mass, Amg::MatrixX& pOut,
131 Amg::MatrixX& /*varout*/);
132 float retrieve_eta_calo(const xAOD::IParticle& part) const;
133 };
134}
135#endif // ZMASSCONSTRAINT_CONSTRAINTFIT_H
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
#define AmgMatrix(rows, cols)
StatusCode doMassFit(const ConstraintFitInput &input, ConstraintFitOutput &output)
Perform the constrained mass fit.
void addFSRParticle(const xAOD::IParticle &part, const TLorentzVector &fsr4Vec, ConstraintFitInput &input)
Add in FSR photon to input, (energy resolution is obtain in method)
double getMassError(const ConstraintFitInput &firstInput, const ConstraintFitInput &secondInput=ConstraintFitInput())
Calculate the mass error without fit - use just the inputs.
int massFitInterface(const ConstraintFitInput &theInput)
void massFitRun(ConstraintFitOutput &output, double zresol=-1.)
void convertCovd0z0PhiThetaPToXYZ(const TLorentzVector &fourVec, const AmgMatrix(5, 5)&covard0z0PhiThetaP, AmgMatrix(5, 5)&covarXYZ) const
float retrieve_eta_calo(const xAOD::IParticle &part) const
ConstraintFit(const std::string &name)
Create a proper constructor for Athena.
StatusCode initialize()
Initialize constraint fit.
bool doSanityChecksOnCovariance(const TLorentzVector &vector, const AmgMatrix(5, 5)&covar) const
void addParticle(const xAOD::Muon &part, ConstraintFitInput &input, MassConstraintMuonType muonType=isCombMCMT)
Add muon to input, must provide the resolution Scale Factor.
ToolHandle< CP::IEgammaCalibrationAndSmearingTool > m_energyRescaler
double massFitCalculation(const Amg::MatrixX &var, double mass, Amg::MatrixX &p0)
void convertCovXYZTod0z0PhiThetaP(const std::vector< TLorentzVector > &particleList, const Amg::MatrixX &covarXYZ, Amg::MatrixX &covard0z0PhiThetaP) const
double massFit(const Amg::MatrixX &, const Amg::MatrixX &var, double mass, Amg::MatrixX &pOut, Amg::MatrixX &)
std::vector< double > m_objmass
ToolHandle< CP::IMuonCalibrationAndSmearingTool > m_mu_resolSFTool
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
Class providing the definition of the 4-vector interface.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Muon_v1 Muon
Reference the current persistent version:
Electron_v1 Electron
Definition of the current "egamma version".