#include <EventShape.h>
Definition at line 6 of file EventShape.h.
Definition at line 14 of file EventShape.cc.
References eigenvalues, i, j, convertSQLitetoXML_cfg::output, p, and python.multivaluedict::sort().
17 for(reco::TrackCollection::const_iterator itTrack =
tracks.begin(); itTrack<
tracks.end(); ++itTrack) {
18 p.push_back(TVector3(itTrack->px(),itTrack->py(),itTrack->pz()));
22 TMatrixDSym MomentumTensor(3);
23 for(std::vector<TVector3>::const_iterator momentum =
p.begin();momentum<
p.end();++momentum) {
24 for(
unsigned int i=0;
i<3;
i++)
25 for(
unsigned int j=0;
j<=
i;
j++) {
26 MomentumTensor[
i][
j] += momentum[
i]*momentum[
j];
29 MomentumTensor*=1/(MomentumTensor[0][0]+MomentumTensor[1][1]+MomentumTensor[2][2]);
31 TMatrixDSymEigen eigen(MomentumTensor);
32 TVectorD eigenvals = eigen.GetEigenValues();
std::vector< float > eigenvalues
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
std::vector< TVector3 > p
Definition at line 244 of file EventShape.cc.
References i, j, and python.multivaluedict::sort().
Referenced by TrackerDpgAnalysis::analyze().
247 if (
tracks.size()==0)
return 0;
249 TMatrixDSym MomentumTensor(3);
250 for(reco::TrackCollection::const_iterator itTrack =
tracks.begin(); itTrack<
tracks.end(); ++itTrack) {
251 std::vector<double> momentum(3);
252 momentum[0] = itTrack->px();
253 momentum[1] = itTrack->py();
254 momentum[2] = itTrack->pz();
255 for(
unsigned int i=0;
i<3;
i++)
256 for(
unsigned int j=0;
j<=
i;
j++) {
257 MomentumTensor[
i][
j] += momentum[
i]*momentum[
j];
260 MomentumTensor*=1/(MomentumTensor[0][0]+MomentumTensor[1][1]+MomentumTensor[2][2]);
262 TMatrixDSymEigen eigen(MomentumTensor);
263 TVectorD eigenvals = eigen.GetEigenValues();
264 vector<float> eigenvaluess(3);
265 eigenvaluess[0] = eigenvals[0];
266 eigenvaluess[1] = eigenvals[1];
267 eigenvaluess[2] = eigenvals[2];
268 sort(eigenvaluess.begin(),eigenvaluess.end());
270 return ( 1.5*eigenvaluess[0]);
float EventShape::aplanarity |
( |
| ) |
const |
Definition at line 273 of file EventShape.cc.
References i, j, and python.multivaluedict::sort().
Referenced by TrackerDpgAnalysis::analyze().
276 if (
tracks.size()==0)
return 0;
278 TMatrixDSym MomentumTensor(3);
279 for(reco::TrackCollection::const_iterator itTrack =
tracks.begin(); itTrack<
tracks.end(); ++itTrack) {
280 std::vector<double> momentum(3);
281 momentum[0] = itTrack->px();
282 momentum[1] = itTrack->py();
283 momentum[2] = itTrack->pz();
284 for(
unsigned int i=0;
i<3;
i++)
285 for(
unsigned int j=0;
j<=
i;
j++) {
286 MomentumTensor[
i][
j] += momentum[
i]*momentum[
j];
289 MomentumTensor*=1/(MomentumTensor[0][0]+MomentumTensor[1][1]+MomentumTensor[2][2]);
291 TMatrixDSymEigen eigen(MomentumTensor);
292 TVectorD eigenvals = eigen.GetEigenValues();
293 vector<float> eigenvaluess(3);
294 eigenvaluess[0] = eigenvals[0];
295 eigenvaluess[1] = eigenvals[1];
296 eigenvaluess[2] = eigenvals[2];
297 sort(eigenvaluess.begin(),eigenvaluess.end());
299 return (eigenvaluess[0]/eigenvaluess[1]);
float EventShape::planarity |
( |
| ) |
const |
Definition at line 213 of file EventShape.cc.
References i, j, and python.multivaluedict::sort().
Referenced by TrackerDpgAnalysis::analyze().
216 if(
tracks.size()==0)
return 0;
219 TMatrixDSym MomentumTensor(3);
220 for(reco::TrackCollection::const_iterator itTrack =
tracks.begin(); itTrack<
tracks.end(); ++itTrack) {
221 std::vector<double> momentum(3);
222 momentum[0] = itTrack->px();
223 momentum[1] = itTrack->py();
224 momentum[2] = itTrack->pz();
225 for(
unsigned int i=0;
i<3;
i++)
226 for(
unsigned int j=0;
j<=
i;
j++) {
227 MomentumTensor[
i][
j] += momentum[
i]*momentum[
j];
230 MomentumTensor*=1/(MomentumTensor[0][0]+MomentumTensor[1][1]+MomentumTensor[2][2]);
232 TMatrixDSymEigen eigen(MomentumTensor);
233 TVectorD eigenvals = eigen.GetEigenValues();
234 vector<float> eigenvaluess(3);
235 eigenvaluess[0] = eigenvals[0];
236 eigenvaluess[1] = eigenvals[1];
237 eigenvaluess[2] = eigenvals[2];
238 sort(eigenvaluess.begin(),eigenvaluess.end());
240 float sph = ( 1.5*(1-eigenvaluess[2]));
float EventShape::sphericity |
( |
| ) |
const |
Definition at line 124 of file EventShape.cc.
References i, j, gen::k, convertSQLitetoXML_cfg::output, createTree::pp, v, x, detailsBasic3DVector::y, detailsBasic3DVector::z, and zero.
Referenced by TrackerDpgAnalysis::analyze().
126 std::vector<TVector3>
pp;
127 uint32_t Np =
tracks.size();
129 for(reco::TrackCollection::const_iterator itTrack =
tracks.begin(); itTrack<
tracks.end(); ++itTrack) {
130 pp.push_back(TVector3(itTrack->px(),itTrack->py(),itTrack->pz()));
133 TVector3
zero(0.,0.,0.);
139 TVector3 vn, vm, vc, vl;
140 for(
unsigned int i=0;
i< Np-1;
i++)
141 for(
unsigned int j=
i+1;
j < Np;
j++) {
142 vc = pp[
i].Cross(pp[
j]);
144 for(
unsigned int k=0;
k<Np;
k++)
145 if ((
k !=
i) && (
k != j)) {
146 if (pp[
k].Dot(vc) >= 0.) vl = vl + pp[
k];
147 else vl = vl - pp[
k];
150 vn = vl + pp[
j] + pp[
i];
156 vn = vl + pp[
j] - pp[
i];
162 vn = vl - pp[
j] + pp[
i];
168 vn = vl - pp[
j] - pp[
i];
176 for(
int iter=1; iter<=4; iter++) {
178 for(
unsigned int i=0;
i< Np;
i++)
179 if (vm.Dot(pp[
i]) >= 0.)
184 if (vnew == vmax)
break;
191 if (pp[0].Dot(pp[1]) >= 0.)
192 qtbo = pp[0] + pp[1];
194 qtbo = pp[0] - pp[1];
203 for(
unsigned int i=0; i < Np; i++) vsum = vsum + pp[i].Mag();
206 float x = qtbo.X()/vnew;
207 float y = qtbo.Y()/vnew;
208 float z = qtbo.Z()/vnew;
209 output.SetPxPyPzE(x, y, z, v);
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition at line 39 of file EventShape.cc.
References i, j, gen::k, convertSQLitetoXML_cfg::output, p, v, x, detailsBasic3DVector::y, detailsBasic3DVector::z, and zero.
43 TVector3
zero(0.,0.,0.);
45 uint32_t Np =
p.size();
50 TVector3 vn, vm, vc, vl;
51 for(
unsigned int i=0;
i< Np-1;
i++)
52 for(
unsigned int j=
i+1;
j < Np;
j++) {
53 vc =
p[
i].Cross(
p[
j]);
55 for(
unsigned int k=0;
k<Np;
k++)
56 if ((
k !=
i) && (
k != j)) {
57 if (
p[
k].Dot(vc) >= 0.) vl = vl +
p[
k];
61 vn = vl +
p[
j] +
p[
i];
67 vn = vl + p[
j] - p[
i];
73 vn = vl - p[
j] + p[
i];
79 vn = vl - p[
j] - p[
i];
87 for(
int iter=1; iter<=4; iter++) {
89 for(
unsigned int i=0;
i< Np;
i++)
90 if (vm.Dot(p[
i]) >= 0.)
95 if (vnew == vmax)
break;
102 if (p[0].Dot(p[1]) >= 0.)
114 for(
unsigned int i=0; i < Np; i++) vsum = vsum + p[i].Mag();
117 float x = qtbo.X()/vnew;
118 float y = qtbo.Y()/vnew;
119 float z = qtbo.Z()/vnew;
120 output.SetPxPyPzE(x, y, z, v);
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
std::vector< TVector3 > p
std::vector<float> EventShape::eigenvalues |
|
private |
std::vector<TVector3> EventShape::p |
|
private |