ATLAS Offline Software
Loading...
Searching...
No Matches
NeutralParticle_v1.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Misc includes
6#include <bitset>
7#include <vector>
8
9// EDM include(s):
11
12// Local include(s):
14// #include "xAODTracking/VertexContainer.h" FIXME - need to get ELs working to vertices for neutrals - currently causes compilation failure EJWM
16
17
18namespace xAOD {
19
23
28
30 if(this == &tp) return *this;
31
32 if(!hasStore() ) makePrivateStore();
33 this->IParticle::operator=( tp );
34#ifndef XAOD_ANALYSIS
35 m_perigeeParameters.reset();
36#endif // not XAOD_ANALYSIS
37 return *this;
38 }
39
42
43 double NeutralParticle_v1::pt() const {
44 return genvecP4().Pt();
45 }
46
47 double NeutralParticle_v1::eta() const {
48 return genvecP4().Eta();
49 }
50
52
53
57
59
60
62 return genvecP4().Rapidity();
63 }
64
66 using namespace std;
67 float p = 1/fabs(oneOverP());
68 float thetaT = theta();
69 float phiT = phi();
70 float sinTheta= sin(thetaT);
71 float px = p*sinTheta*cos(phiT);
72 float py = p*sinTheta*sin(phiT);
73 float pz = p*cos(thetaT);
74 return GenVecFourMom_t(px, py, pz, m());
75 }
76
78 using namespace std;
80 float p = 1/fabs(oneOverP());
81 float thetaT = theta();
82 float phiT = phi();
83 float sinTheta= sin(thetaT);
84 float px = p*sinTheta*cos(phiT);
85 float py = p*sinTheta*sin(phiT);
86 float pz = p*cos(thetaT);
87 float e = pow (m(),2) +
88 pow( px,2) + pow( py,2) + pow( pz,2);
89 p4.SetPxPyPzE( px, py, pz, sqrt(e) );
90 return p4;
91 }
92
93 Type::ObjectType NeutralParticle_v1::type() const {
94 return Type::NeutralParticle;
95 }
96
97
103
104 const DefiningParameters_t NeutralParticle_v1::definingParameters() const{
105 DefiningParameters_t tmp;
106 tmp << d0(),z0(),phi0(),theta(),oneOverP();
107 return tmp;
108 }
109
110 void NeutralParticle_v1::setDefiningParameters(float d0, float z0, float phi0, float theta, float oneOverP) {
111#ifndef XAOD_ANALYSIS
112 if(m_perigeeParameters.isValid()) {
113 m_perigeeParameters.reset();
114 }
115#endif // not XAOD_ANALYSIS
116 static const Accessor< float > acc1( "d0" );
117 acc1( *this ) = d0;
118
119 static const Accessor< float > acc2( "z0" );
120 acc2( *this ) = z0;
121
122 static const Accessor< float > acc3( "phi" );
123 acc3( *this ) = phi0;
124
125 static const Accessor< float > acc4( "theta" );
126 acc4( *this ) = theta;
127
128 static const Accessor< float > acc5( "oneOverP" );
129 acc5( *this ) = oneOverP;
130
131 return;
132 }
133
134 void NeutralParticle_v1::setDefiningParametersCovMatrix(const xAOD::ParametersCovMatrix_t& cov){
135#ifndef XAOD_ANALYSIS
136 if(m_perigeeParameters.isValid()) {
137 m_perigeeParameters.reset();
138 }
139#endif // not XAOD_ANALYSIS
140
141 static const Accessor< std::vector<float> > acc( "definingParametersCovMatrix" );
142 std::vector<float>& v = acc(*this);
143 v.reserve(15);
144 for (size_t irow = 0; irow<5; ++irow)
145 for (size_t icol =0; icol<=irow; ++icol)
146 v.push_back(cov(icol,irow));
147 }
148
149 const xAOD::ParametersCovMatrix_t NeutralParticle_v1::definingParametersCovMatrix() const {
150 static const Accessor< std::vector<float> > acc( "definingParametersCovMatrix" );
151 std::vector<float> v = acc(*this);
152 std::vector<float>::const_iterator it = v.begin();
153 xAOD::ParametersCovMatrix_t cov;
154 cov.setZero();
155 for (size_t irow = 0; irow<5; ++irow){
156 for (size_t icol =0; icol<=irow; ++icol){
157 cov.fillSymmetric(icol,irow, *it++);
158 }
159 }
160 return cov;
161 }
162
163 const std::vector<float>& NeutralParticle_v1::definingParametersCovMatrixVec() const {
164 static const Accessor< std::vector<float> > acc( "definingParametersCovMatrix" );
165 return acc(*this);
166 }
167
168 void NeutralParticle_v1::setDefiningParametersCovMatrixVec(const std::vector<float>& cov){
169 static const Accessor< std::vector<float> > acc( "definingParametersCovMatrix" );
170 acc(*this)=cov;
171 }
172
176
177 void NeutralParticle_v1::setParametersOrigin(float x, float y, float z){
178 static const Accessor< float > acc1( "vx" );
179 acc1( *this ) = x;
180
181 static const Accessor< float > acc2( "vy" );
182 acc2( *this ) = y;
183
184 static const Accessor< float > acc3( "vz" );
185 acc3( *this ) = z;
186 }
187
188#ifndef XAOD_ANALYSIS
190
191 // Require the cache to be valid and check if the cached pointer has been set
192 if(m_perigeeParameters.isValid()){
193 return *(m_perigeeParameters.ptr());
194 }
195 static const Accessor< float > acc1( "d0" );
196 static const Accessor< float > acc2( "z0" );
197 static const Accessor< float > acc3( "phi" );
198 static const Accessor< float > acc4( "theta" );
199 static const Accessor< float > acc5( "oneOverP" );
200 static const Accessor< std::vector<float> > acc6( "definingParametersCovMatrix" );
201 ParametersCovMatrix_t cov;
202 cov.setZero();
203 auto it= acc6(*this).begin();
204 for (size_t irow = 0; irow<5; ++irow){
205 for (size_t icol =0; icol<=irow; ++icol){
206 cov.fillSymmetric(irow,icol,*it++) ;
207 }
208 }
209 Trk::NeutralPerigee tmpPerigeeParameters(acc1(*this),acc2(*this),acc3(*this),acc4(*this),acc5(*this),
210 Trk::PerigeeSurface(Amg::Vector3D(vx(),vy(),vz())),std::move(cov));
211 m_perigeeParameters.set(tmpPerigeeParameters);
212 return *(m_perigeeParameters.ptr());
213 }
214#endif // not XAOD_ANALYSIS
215
217#ifndef XAOD_ANALYSIS
218 m_perigeeParameters.reset();
219#endif // not XAOD_ANALYSIS
220 }
221
222} // namespace xAOD
Scalar theta() const
theta method
#define AUXSTORE_PRIMITIVE_GETTER_WITH_CAST(CL, PERSTYPE, TRANSTYPE, NAME)
Macro creating a getter function with a type conversion.
#define AUXSTORE_PRIMITIVE_GETTER(CL, TYPE, NAME)
Macro creating the reader function for a primitive auxiliary property.
A number of constexpr particle constants to avoid hardcoding them directly in various places.
#define y
#define x
#define z
constexpr int pow(int base, int exp) noexcept
void makePrivateStore()
Create a new (empty) private store for this object.
bool hasStore() const
Return true if this object has an associated store.
Class describing the Line to which the Perigee refers to.
IParticle & operator=(const IParticle &)=default
IParticle()=default
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
Class describing a NeutralParticle.
const Trk::NeutralPerigee & perigeeParameters() const
Returns the Trk::NeutralPerigee track parameters.
float d0() const
Returns the parameter.
const DefiningParameters_t definingParameters() const
Returns a SVector of the Perigee track parameters.
float oneOverP() const
Returns the parameter.
void setDefiningParametersCovMatrixVec(const std::vector< float > &cov)
Set the defining parameters covariance matrix using a length 15 vector.
const ParametersCovMatrix_t definingParametersCovMatrix() const
Returns the 5x5 symmetric matrix containing the defining parameters covariance matrix.
float theta() const
Returns the parameter, which has range 0 to .
float vz() const
The z origin for the parameters.
virtual double m() const
The invariant mass of the particle..
void resetCache()
Reset the internal cache of the object.
virtual double pt() const
The transverse momentum ( ) of the particle.
NeutralParticle_v1 & operator=(const NeutralParticle_v1 &tp)
Assignment operator. This can involve creating and copying an Auxilary store, and so should be used s...
void setDefiningParameters(float d0, float z0, float phi0, float theta, float qOverP)
Set the defining parameters.
IParticle::FourMom_t FourMom_t
Definition of the 4-momentum type.
virtual Type::ObjectType type() const
The type of the object as a simple enumeration.
void setDefiningParametersCovMatrix(const ParametersCovMatrix_t &cov)
Set the defining parameters covariance matrix.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
float phi0() const
Returns the parameter, which has range to .
CxxUtils::CachedValue< Trk::NeutralPerigee > m_perigeeParameters
Cached NeutralPerigee, built from this object.
float vx() const
The x origin for the parameters.
NeutralParticle_v1()
Default constructor.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double e() const
The total energy of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
void setParametersOrigin(float x, float y, float z)
Set the origin for the parameters.
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzM4D< double > > GenVecFourMom_t
Base 4 Momentum type for TrackParticle.
GenVecFourMom_t genvecP4() const
The full 4-momentum of the particle : GenVector form.
const std::vector< float > & definingParametersCovMatrixVec() const
Returns the vector of the covariance values - 15 elements.
float vy() const
The y origin for the parameters.
float z0() const
Returns the parameter.
virtual double rapidity() const
The true rapidity (y) of the particle.
Eigen::Matrix< double, 3, 1 > Vector3D
constexpr double chargedPionMassInMeV
the mass of the charged pion (in MeV)
ParametersT< NeutralParametersDim, Neutral, PerigeeSurface > NeutralPerigee
STL namespace.
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
static const SG::AuxElement::Accessor< ElementLink< IParticleContainer > > acc("originalObjectLink")
Object used for setting/getting the dynamic decoration in question.