FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
EventShape Class Reference

#include <EventShape.h>

Public Member Functions

 EventShape (const CandList &l, TLorentzVector cms)
 
 ~EventShape ()
 
int NParticles () const
 
int NCharged () const
 
int NNeutral () const
 
double PmaxLab () const
 
double PmaxCms () const
 
double PminLab () const
 
double PminCms () const
 
double Ptmax () const
 
double Ptmin () const
 
double PRapmax () const
 
double PtSumLab () const
 
double NeutEtSumLab () const
 
double NeutESumLab () const
 
double ChrgPtSumLab () const
 
double ChrgPSumLab () const
 
double PtSumCms () const
 
double NeutEtSumCms () const
 
double NeutESumCms () const
 
double ChrgPtSumCms () const
 
double ChrgPSumCms () const
 
int MultPminLab (double pmin)
 
int MultPmaxLab (double pmax)
 
int MultPminCms (double pmin)
 
int MultPmaxCms (double pmax)
 
int MultPtminLab (double ptmin)
 
int MultPtmaxLab (double ptmax)
 
int MultPtminCms (double ptmin)
 
int MultPtmaxCms (double ptmax)
 
int MultNeutEminLab (double emin)
 
int MultNeutEmaxLab (double emax)
 
int MultNeutEminCms (double emin)
 
int MultNeutEmaxCms (double emax)
 
int MultChrgPminLab (double pmin)
 
int MultChrgPmaxLab (double pmax)
 
int MultChrgPminCms (double pmin)
 
int MultChrgPmaxCms (double pmax)
 
double SumPminLab (double pmin)
 
double SumPmaxLab (double pmax)
 
double SumPminCms (double pmin)
 
double SumPmaxCms (double pmax)
 
double SumPtminLab (double ptmin)
 
double SumPtmaxLab (double ptmax)
 
double SumPtminCms (double ptmin)
 
double SumPtmaxCms (double ptmax)
 
double SumNeutEminLab (double emin)
 
double SumNeutEmaxLab (double emax)
 
double SumNeutEminCms (double emin)
 
double SumNeutEmaxCms (double emax)
 
double SumChrgPminLab (double pmin)
 
double SumChrgPmaxLab (double pmax)
 
double SumChrgPminCms (double pmin)
 
double SumChrgPmaxCms (double pmax)
 
double Sphericity ()
 
double Aplanarity ()
 
double Planarity ()
 
double Circularity ()
 
double FoxWolfMomH (int order)
 
double FoxWolfMomR (int order)
 
double Thrust ()
 
TVector3 ThrustVector ()
 

Private Member Functions

void ComputeSphericity ()
 
double Eps (const TVector3 v1, const TVector3 v2)
 
double Legendre (int l, double x)
 

Private Attributes

std::vector< TLorentzVector > fLabList
 
std::vector< TLorentzVector > fCmsList
 
std::vector< int > fCharge
 
int fnChrg
 
int fnNeut
 
int fN
 
double fpmaxlab
 
double fpmaxcms
 
double fpminlab
 
double fpmincms
 
double fptmax
 
double fptmin
 
double fprapmax
 
double fptsumlab
 
double fneutetsumlab
 
double fneutesumlab
 
double fchrgptsumlab
 
double fchrgpsumlab
 
double fptsumcms
 
double fneutetsumcms
 
double fneutesumcms
 
double fchrgptsumcms
 
double fchrgpsumcms
 
double fsph
 
double fapl
 
double fpla
 
double fcir
 
double fFWmom [FWMAX]
 
bool fFWready
 
double fthr
 
TVector3 fThrVect
 
TVector3 fBoost
 

Detailed Description

Definition at line 9 of file EventShape.h.

Constructor & Destructor Documentation

EventShape::EventShape ( const CandList l,
TLorentzVector  cms 
)

Definition at line 152 of file EventShape.h.

References fBoost, fCharge, fchrgpsumcms, fchrgpsumlab, fchrgptsumcms, fchrgptsumlab, fCmsList, fFWmom, fLabList, fN, fnChrg, fneutesumcms, fneutesumlab, fneutetsumcms, fneutetsumlab, fnNeut, fpmaxcms, fpmaxlab, fpmincms, fpminlab, fprapmax, fptmax, fptmin, fptsumcms, fptsumlab, fThrVect, FWMAX, and i.

152  :
153  fnChrg(0), fnNeut(0), fN(0), fpmaxlab(0.), fpmaxcms(0.),fpminlab(0.), fpmincms(0.), fptmax(0.), fptmin(0.),
156  fsph(-1.), fapl(-1.), fpla(-1.), fcir(-1.), fFWready(false), fthr(-1.)
157 {
158  int i;
159 
160  // initialize more complex things
161  fThrVect.SetXYZ(0.,0.,0.);
162  fBoost=cms.BoostVector();
163  for (i=0;i<=FWMAX;++i) fFWmom[i]=0.;
164 
165 
166  double pmax=0., ptmax=0., pmaxcms=0.;
167  double pmin=1000., ptmin=1000., pmincms=1000.;
168  double prapmax=0.;
169 
170  fLabList.clear();
171  fCmsList.clear();
172 
173  for (i=0;i<l.size();++i)
174  {
175  TLorentzVector lv(l[i].P4());
176  int chrg(l[i].Charge());
177 
178  fN++;
179 
180  // cache multiplicities
181  if (chrg==0) fnNeut++;
182  else fnChrg++;
183 
184  // cache unboosted 4-vectors
185  fLabList.push_back(lv);
186  // cache charges
187  fCharge.push_back(chrg);
188 
189  // sum momentum variables (lab)
190  fptsumlab += lv.Pt();
191  if (chrg==0)
192  {
193  fneutetsumlab += lv.Pt();
194  fneutesumlab += lv.E();
195  }
196  else
197  {
198  fchrgptsumlab += lv.Pt();
199  fchrgpsumlab += lv.P();
200  }
201 
202  // cache maximum momenta in lab
203  if (lv.P()>pmax) pmax=lv.P();
204  if (lv.Pt()>ptmax) ptmax=lv.Pt();
205 
206  // cache minimum momenta in lab
207  if (lv.P()<pmin) pmin=lv.P();
208  if (lv.Pt()<ptmin) ptmin=lv.Pt();
209 
210  // cache pseudorapidity
211  if (lv.PseudoRapidity()>prapmax) prapmax = lv.PseudoRapidity();
212 
213  // cache boosted vectors
214  lv.Boost(-fBoost);
215  fCmsList.push_back(lv);
216 
217  // sum momentum variables (cms)
218  fptsumcms += lv.Pt();
219  if (chrg==0)
220  {
221  fneutetsumcms += lv.Pt();
222  fneutesumcms += lv.E();
223  }
224  else
225  {
226  fchrgptsumcms += lv.Pt();
227  fchrgpsumcms += lv.P();
228  }
229 
230  // cache maximum momenta in cms
231  if (lv.P()>pmaxcms) pmaxcms=lv.P();
232  // cache minimum momenta in cms
233  if (lv.P()<pmincms) pmincms=lv.P();
234  }
235 
236  fpmaxlab = pmax;
237  fptmax = ptmax;
238  fpmaxcms = pmaxcms;
239  fpminlab = pmin;
240  fptmin = ptmin;
241  fpmincms = pmincms;
242  fprapmax = prapmax;
243 }
double fprapmax
Definition: EventShape.h:122
double fsph
Definition: EventShape.h:136
double fthr
Definition: EventShape.h:144
Int_t i
Definition: run_full.C:25
double fFWmom[FWMAX]
Definition: EventShape.h:141
double fpmaxlab
Definition: EventShape.h:116
double fptsumlab
Definition: EventShape.h:124
TVector3 fThrVect
Definition: EventShape.h:145
double fneutesumlab
Definition: EventShape.h:126
double fpmaxcms
Definition: EventShape.h:117
double fneutetsumlab
Definition: EventShape.h:125
double fcir
Definition: EventShape.h:139
double fneutesumcms
Definition: EventShape.h:132
double fpla
Definition: EventShape.h:138
TVector3 fBoost
Definition: EventShape.h:146
double fpminlab
Definition: EventShape.h:118
double fapl
Definition: EventShape.h:137
bool fFWready
Definition: EventShape.h:142
double fchrgptsumlab
Definition: EventShape.h:127
double fchrgpsumcms
Definition: EventShape.h:134
double fptmax
Definition: EventShape.h:120
std::vector< int > fCharge
Definition: EventShape.h:110
#define FWMAX
Definition: EventShape.h:7
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
double fptmin
Definition: EventShape.h:121
double fptsumcms
Definition: EventShape.h:130
double fpmincms
Definition: EventShape.h:119
double fneutetsumcms
Definition: EventShape.h:131
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double fchrgpsumlab
Definition: EventShape.h:128
double fchrgptsumcms
Definition: EventShape.h:133
EventShape::~EventShape ( )
inline

Definition at line 12 of file EventShape.h.

12 {};

Member Function Documentation

double EventShape::Aplanarity ( )

Definition at line 299 of file EventShape.h.

References ComputeSphericity(), and fapl.

Referenced by softtrigger_kin5(), and writeTuple().

300 {
301  if (fapl<0.) ComputeSphericity();
302 
303  return fapl;
304 }
void ComputeSphericity()
Definition: EventShape.h:247
double fapl
Definition: EventShape.h:137
double EventShape::ChrgPSumCms ( ) const
inline

Definition at line 40 of file EventShape.h.

References fchrgpsumcms.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

40 {return fchrgpsumcms;} // sum of momenta of charged (cms)
double fchrgpsumcms
Definition: EventShape.h:134
double EventShape::ChrgPSumLab ( ) const
inline

Definition at line 33 of file EventShape.h.

References fchrgpsumlab.

Referenced by softtrigger_kin5().

33 {return fchrgpsumlab;} // sum of momenta of charged (lab)
double fchrgpsumlab
Definition: EventShape.h:128
double EventShape::ChrgPtSumCms ( ) const
inline

Definition at line 39 of file EventShape.h.

References fchrgptsumcms.

Referenced by softtrigger_kin5().

39 {return fchrgptsumcms;} // sum of pt of charged (cms)
double fchrgptsumcms
Definition: EventShape.h:133
double EventShape::ChrgPtSumLab ( ) const
inline

Definition at line 32 of file EventShape.h.

References fchrgptsumlab.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

32 {return fchrgptsumlab;} // sum of pt of charged (lab)
double fchrgptsumlab
Definition: EventShape.h:127
double EventShape::Circularity ( )

Definition at line 308 of file EventShape.h.

References ComputeSphericity(), and fcir.

Referenced by softtrigger_kin5(), and writeTuple().

309 {
310  if (fcir<0.) ComputeSphericity();
311 
312  return fcir;
313 }
void ComputeSphericity()
Definition: EventShape.h:247
double fcir
Definition: EventShape.h:139
void EventShape::ComputeSphericity ( )
private

Definition at line 247 of file EventShape.h.

References fapl, fcir, fCmsList, fN, fpla, fsph, i, and v.

Referenced by Aplanarity(), Circularity(), Planarity(), and Sphericity().

248 {
249  if( fN==0 ) return;
250 
251  double stot=0, sxx=0, sxy=0, sxz=0, syy=0, syz=0, szz=0;
252  int i;
253 
254  for (i=0;i<fN;++i)
255  {
256  TVector3 v(fCmsList[i].Vect());
257 
258  sxx += v.X()*v.X(); sxy += v.X()*v.Y(); sxz += v.X()*v.Z();
259  syy += v.Y()*v.Y(); syz += v.Y()*v.Z();
260  szz += v.Z()*v.Z();
261 
262  stot += v.Mag2();
263  }
264 
265  TMatrixD sm(3,3);
266  sm(0,0) = sxx/stot; sm(0,1) = sxy/stot; sm(0,2) = sxz/stot;
267  sm(1,0) = sxy/stot; sm(1,1) = syy/stot; sm(1,2) = syz/stot;
268  sm(2,0) = sxz/stot; sm(2,1) = syz/stot; sm(2,2) = szz/stot;
269 
270  TMatrixDEigen ei(sm);
271  TMatrixD eiv=ei.GetEigenValues();
272 
273  fsph = 1.5 * (eiv(1,1) + eiv(2,2)); // 3/2 (lam_2+lam_3)
274  fapl = 1.5 * eiv(2,2); // 3/2 lam_3
275  fpla = eiv(1,1) - eiv(2,2); // lam_2-lam_3
276  fcir = 2.*eiv(1,1)/(eiv(1,1)+eiv(0,0)); // 2* min(lam_1, lam_2)/(lam_1+lam_2)
277 }
double fsph
Definition: EventShape.h:136
Int_t i
Definition: run_full.C:25
double fcir
Definition: EventShape.h:139
double fpla
Definition: EventShape.h:138
__m128 v
Definition: P4_F32vec4.h:4
double fapl
Definition: EventShape.h:137
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
double EventShape::Eps ( const TVector3  v1,
const TVector3  v2 
)
inlineprivate

Definition at line 105 of file EventShape.h.

Referenced by Thrust().

105 {return (v1*v2)>0. ? 1. : -1.;} // aux for Thrust
TVector3 v1
Definition: bump_analys.C:40
TVector3 v2
Definition: bump_analys.C:40
double EventShape::FoxWolfMomH ( int  order)

Definition at line 404 of file EventShape.h.

References cos(), fCmsList, fFWmom, fFWready, fN, FWMAX, i, Legendre(), p1, p2, and s.

Referenced by FoxWolfMomR().

405 {
406  if (order>FWMAX) return -1.;
407 
408  // already computed FW moments
409  if (fFWready) return fFWmom[order];
410 
411  if( fN==0 ) return -1.;
412 
413  double s = 0.;
414  int i,j,l;
415 
416  for (i=0; i<fN-1; ++i)
417  {
418  // this candidate's 3-momentum
419  TVector3 p1(fCmsList[i].Vect());
420  double pmag1 = p1.Mag();
421 
422  // loop over other candidates, starting at the next one in the list
423 
424  for (j=i+1; j<fN; ++j)
425  {
426  // this candidate's 3-momentum
427  TVector3 p2(fCmsList[j].Vect());
428  double pmag2 = p2.Mag();
429 
430  // the cosine of the angle between the two candidates
431  double cosPhi = cos ( p1.Angle(p2) );
432 
433  // the contribution of this pair of track
434  // (note the factor 2 : the pair enters the sum twice)
435  for( l=0; l<=FWMAX; l++ )
436  fFWmom[l] += 2 * pmag1 * pmag2 * Legendre( l, cosPhi );
437 
438  }
439  // contribution of this track
440  for( l=0; l<=FWMAX; l++ )
441  fFWmom[l] += pmag1 * pmag1 * Legendre( l, 1. );
442 
443  // total energy
444  s += fCmsList[i].Energy();
445  }
446 
447  // well ...
448  if( s<=0. ) return -1.;
449  double s2=s*s;
450 
451  // normalize Fox Wolfram Moments
452  for( i=0; i<=FWMAX; i++) fFWmom[i]/=s2 ;
453 
454  fFWready = true;
455 
456  return fFWmom[order];
457 }
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
Int_t i
Definition: run_full.C:25
double fFWmom[FWMAX]
Definition: EventShape.h:141
TLorentzVector s
Definition: Pnd2DStar.C:50
double Legendre(int l, double x)
Definition: EventShape.h:381
bool fFWready
Definition: EventShape.h:142
TPad * p2
Definition: hist-t7.C:117
#define FWMAX
Definition: EventShape.h:7
TPad * p1
Definition: hist-t7.C:116
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::FoxWolfMomR ( int  order)

Definition at line 462 of file EventShape.h.

References FoxWolfMomH().

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

463 {
464  return FoxWolfMomH(order)/FoxWolfMomH(0);
465 }
double FoxWolfMomH(int order)
Definition: EventShape.h:404
double EventShape::Legendre ( int  l,
double  x 
)
private

Definition at line 381 of file EventShape.h.

References fabs(), and x.

Referenced by FoxWolfMomH().

382 {
383  if (fabs(x)>1.) return -999.;
384 
385  if (l==0) return 1.;
386 
387  double pmm = 1.;
388  double pmmp1 = x;
389 
390  if (l>1)
391  {
392  for(int ll=2; ll<=l; ll++)
393  {
394  double pll = (x * (2 * ll - 1) * pmmp1 - (ll - 1) * pmm) / (ll);
395  pmm = pmmp1;
396  pmmp1 = pll;
397  }
398  }
399  return pmmp1;
400 }
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
Double_t x
int EventShape::MultChrgPmaxCms ( double  pmax)

Definition at line 641 of file EventShape.h.

References cnt, fCharge, fCmsList, fN, i, and P.

642 {
643  int cnt=0;
644  for (int i=0;i<fN;++i)
645  if (fCharge[i]!=0 && fCmsList[i].P()<pmax) cnt++;
646 
647  return cnt;
648 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
GeV c P
Int_t cnt
Definition: hist-t7.C:106
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
int EventShape::MultChrgPmaxLab ( double  pmax)

Definition at line 619 of file EventShape.h.

References cnt, fCharge, fLabList, fN, i, and P.

620 {
621  int cnt=0;
622  for (int i=0;i<fN;++i)
623  if (fCharge[i]!=0 && fLabList[i].P()<pmax) cnt++;
624 
625  return cnt;
626 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
Int_t cnt
Definition: hist-t7.C:106
int EventShape::MultChrgPminCms ( double  pmin)

Definition at line 630 of file EventShape.h.

References cnt, fCharge, fCmsList, fN, i, and P.

631 {
632  int cnt=0;
633  for (int i=0;i<fN;++i)
634  if (fCharge[i]!=0 && fCmsList[i].P()>pmin) cnt++;
635 
636  return cnt;
637 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
GeV c P
Int_t cnt
Definition: hist-t7.C:106
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
int EventShape::MultChrgPminLab ( double  pmin)

Definition at line 608 of file EventShape.h.

References cnt, fCharge, fLabList, fN, i, and P.

609 {
610  int cnt=0;
611  for (int i=0;i<fN;++i)
612  if (fCharge[i]!=0 && fLabList[i].P()>pmin) cnt++;
613 
614  return cnt;
615 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
Int_t cnt
Definition: hist-t7.C:106
int EventShape::MultNeutEmaxCms ( double  emax)

Definition at line 595 of file EventShape.h.

References cnt, fCharge, fCmsList, fN, and i.

596 {
597  int cnt=0;
598  for (int i=0;i<fN;++i)
599  if (fCharge[i]==0 && fCmsList[i].E()<emax) cnt++;
600 
601  return cnt;
602 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
Int_t cnt
Definition: hist-t7.C:106
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
int EventShape::MultNeutEmaxLab ( double  emax)

Definition at line 573 of file EventShape.h.

References cnt, fCharge, fLabList, fN, and i.

574 {
575  int cnt=0;
576  for (int i=0;i<fN;++i)
577  if (fCharge[i]==0 && fLabList[i].E()<emax) cnt++;
578 
579  return cnt;
580 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
Int_t cnt
Definition: hist-t7.C:106
int EventShape::MultNeutEminCms ( double  emin)

Definition at line 584 of file EventShape.h.

References cnt, fCharge, fCmsList, fN, and i.

585 {
586  int cnt=0;
587  for (int i=0;i<fN;++i)
588  if (fCharge[i]==0 && fCmsList[i].E()>emin) cnt++;
589 
590  return cnt;
591 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
Int_t cnt
Definition: hist-t7.C:106
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
int EventShape::MultNeutEminLab ( double  emin)

Definition at line 562 of file EventShape.h.

References cnt, fCharge, fLabList, fN, and i.

563 {
564  int cnt=0;
565  for (int i=0;i<fN;++i)
566  if (fCharge[i]==0 && fLabList[i].E()>emin) cnt++;
567 
568  return cnt;
569 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
Int_t cnt
Definition: hist-t7.C:106
int EventShape::MultPmaxCms ( double  pmax)

Definition at line 504 of file EventShape.h.

References cnt, fCmsList, fN, i, and P.

505 {
506  int cnt=0;
507  for (int i=0;i<fN;++i)
508  if (fCmsList[i].P()<pmax) cnt++;
509 
510  return cnt;
511 }
Int_t i
Definition: run_full.C:25
GeV c P
Int_t cnt
Definition: hist-t7.C:106
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
int EventShape::MultPmaxLab ( double  pmax)

Definition at line 482 of file EventShape.h.

References cnt, fLabList, fN, i, and P.

483 {
484  int cnt=0;
485  for (int i=0;i<fN;++i)
486  if (fLabList[i].P()<pmax) cnt++;
487 
488  return cnt;
489 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
Int_t cnt
Definition: hist-t7.C:106
int EventShape::MultPminCms ( double  pmin)

Definition at line 493 of file EventShape.h.

References cnt, fCmsList, fN, i, and P.

Referenced by softtrigger_kin5(), softtrigger_toy12(), and toy_core().

494 {
495  int cnt=0;
496  for (int i=0;i<fN;++i)
497  if (fCmsList[i].P()>pmin) cnt++;
498 
499  return cnt;
500 }
Int_t i
Definition: run_full.C:25
GeV c P
Int_t cnt
Definition: hist-t7.C:106
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
int EventShape::MultPminLab ( double  pmin)

Definition at line 471 of file EventShape.h.

References cnt, fLabList, fN, i, and P.

Referenced by softtrigger_kin5().

472 {
473  int cnt=0;
474  for (int i=0;i<fN;++i)
475  if (fLabList[i].P()>pmin) cnt++;
476 
477  return cnt;
478 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
Int_t cnt
Definition: hist-t7.C:106
int EventShape::MultPtmaxCms ( double  ptmax)

Definition at line 549 of file EventShape.h.

References cnt, fCmsList, fN, and i.

550 {
551  int cnt=0;
552  for (int i=0;i<fN;++i)
553  if (fCmsList[i].Pt()<ptmax) cnt++;
554 
555  return cnt;
556 }
Int_t i
Definition: run_full.C:25
Int_t cnt
Definition: hist-t7.C:106
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
int EventShape::MultPtmaxLab ( double  ptmax)

Definition at line 527 of file EventShape.h.

References cnt, fLabList, fN, and i.

528 {
529  int cnt=0;
530  for (int i=0;i<fN;++i)
531  if (fLabList[i].Pt()<ptmax) cnt++;
532 
533  return cnt;
534 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
Int_t cnt
Definition: hist-t7.C:106
int EventShape::MultPtminCms ( double  ptmin)

Definition at line 538 of file EventShape.h.

References cnt, fCmsList, fN, and i.

539 {
540  int cnt=0;
541  for (int i=0;i<fN;++i)
542  if (fCmsList[i].Pt()>ptmin) cnt++;
543 
544  return cnt;
545 }
Int_t i
Definition: run_full.C:25
Int_t cnt
Definition: hist-t7.C:106
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
int EventShape::MultPtminLab ( double  ptmin)

Definition at line 516 of file EventShape.h.

References cnt, fLabList, fN, and i.

517 {
518  int cnt=0;
519  for (int i=0;i<fN;++i)
520  if (fLabList[i].Pt()>ptmin) cnt++;
521 
522  return cnt;
523 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
Int_t cnt
Definition: hist-t7.C:106
int EventShape::NCharged ( ) const
inline

Definition at line 16 of file EventShape.h.

References fnChrg.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

16 {return fnChrg;} // number of charged candidates
double EventShape::NeutESumCms ( ) const
inline

Definition at line 38 of file EventShape.h.

References fneutesumcms.

Referenced by softtrigger_kin5(), softtrigger_toy12(), and toy_core().

38 {return fneutesumcms;} // sum of energys of neutrals (cms)
double fneutesumcms
Definition: EventShape.h:132
double EventShape::NeutESumLab ( ) const
inline

Definition at line 31 of file EventShape.h.

References fneutesumlab.

Referenced by softtrigger_kin5().

31 {return fneutesumlab;} // sum of energys of neutrals (lab)
double fneutesumlab
Definition: EventShape.h:126
double EventShape::NeutEtSumCms ( ) const
inline

Definition at line 37 of file EventShape.h.

References fneutetsumcms.

Referenced by softtrigger_kin5().

37 {return fneutetsumcms;} // sum of transvers energys of neutrals (cms)
double fneutetsumcms
Definition: EventShape.h:131
double EventShape::NeutEtSumLab ( ) const
inline

Definition at line 30 of file EventShape.h.

References fneutetsumlab.

Referenced by softtrigger_kin5().

30 {return fneutetsumlab;} // sum of transvers energys of neutrals (lab)
double fneutetsumlab
Definition: EventShape.h:125
int EventShape::NNeutral ( ) const
inline

Definition at line 17 of file EventShape.h.

References fnNeut.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

17 {return fnNeut;} // number of neutral candidates
int EventShape::NParticles ( ) const
inline

Definition at line 15 of file EventShape.h.

References fN.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

15 {return fN;} // number of particle candidates
double EventShape::Planarity ( )

Definition at line 290 of file EventShape.h.

References ComputeSphericity(), and fpla.

Referenced by softtrigger_kin5(), and writeTuple().

291 {
292  if (fpla<0.) ComputeSphericity();
293 
294  return fpla;
295 }
void ComputeSphericity()
Definition: EventShape.h:247
double fpla
Definition: EventShape.h:138
double EventShape::PmaxCms ( ) const
inline

Definition at line 21 of file EventShape.h.

References fpmaxcms.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

21 {return fpmaxcms;} // max momentum im cms system
double fpmaxcms
Definition: EventShape.h:117
double EventShape::PmaxLab ( ) const
inline

Definition at line 20 of file EventShape.h.

References fpmaxlab.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

20 {return fpmaxlab;} // max momentum in lab system
double fpmaxlab
Definition: EventShape.h:116
double EventShape::PminCms ( ) const
inline

Definition at line 23 of file EventShape.h.

References fpmincms.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

23 {return fpmincms;} // min momentum im cms system
double fpmincms
Definition: EventShape.h:119
double EventShape::PminLab ( ) const
inline

Definition at line 22 of file EventShape.h.

References fpminlab.

Referenced by softtrigger_kin5().

22 {return fpminlab;} // min momentum in lab system
double fpminlab
Definition: EventShape.h:118
double EventShape::PRapmax ( ) const
inline

Definition at line 26 of file EventShape.h.

References fprapmax.

Referenced by softtrigger_kin5().

26 {return fprapmax;} // max pseudorapidity (lab)
double fprapmax
Definition: EventShape.h:122
double EventShape::Ptmax ( ) const
inline

Definition at line 24 of file EventShape.h.

References fptmax.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

24 {return fptmax;} // max pt (same for lab and cms)
double fptmax
Definition: EventShape.h:120
double EventShape::Ptmin ( ) const
inline

Definition at line 25 of file EventShape.h.

References fptmin.

Referenced by softtrigger_kin5().

25 {return fptmin;} // min pt (same for lab and cms)
double fptmin
Definition: EventShape.h:121
double EventShape::PtSumCms ( ) const
inline

Definition at line 36 of file EventShape.h.

References fptsumcms.

Referenced by softtrigger_kin5().

36 {return fptsumcms;} // sum of pt in (cms)
double fptsumcms
Definition: EventShape.h:130
double EventShape::PtSumLab ( ) const
inline

Definition at line 29 of file EventShape.h.

References fptsumlab.

Referenced by softtrigger_kin5(), softtrigger_toy12(), toy_core(), and writeTuple().

29 {return fptsumlab;} // sum of pt in (lab)
double fptsumlab
Definition: EventShape.h:124
double EventShape::Sphericity ( )

Definition at line 281 of file EventShape.h.

References ComputeSphericity(), and fsph.

Referenced by softtrigger_kin5(), and writeTuple().

282 {
283  if (fsph<0.) ComputeSphericity();
284 
285  return fsph;
286 }
double fsph
Definition: EventShape.h:136
void ComputeSphericity()
Definition: EventShape.h:247
double EventShape::SumChrgPmaxCms ( double  pmax)

Definition at line 825 of file EventShape.h.

References fCharge, fCmsList, fLabList, fN, i, and P.

826 {
827  double sum=0;
828  for (int i=0;i<fN;++i)
829  if (fCharge[i]!=0 && fCmsList[i].P()<pmax) sum+=fLabList[i].P();
830 
831  return sum;
832 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::SumChrgPmaxLab ( double  pmax)

Definition at line 803 of file EventShape.h.

References fCharge, fLabList, fN, i, and P.

804 {
805  double sum=0;
806  for (int i=0;i<fN;++i)
807  if (fCharge[i]!=0 && fLabList[i].P()<pmax) sum+=fLabList[i].P();
808 
809  return sum;
810 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
double EventShape::SumChrgPminCms ( double  pmin)

Definition at line 814 of file EventShape.h.

References fCharge, fCmsList, fLabList, fN, i, and P.

Referenced by softtrigger_kin5(), softtrigger_toy12(), and toy_core().

815 {
816  double sum=0;
817  for (int i=0;i<fN;++i)
818  if (fCharge[i]!=0 && fCmsList[i].P()>pmin) sum+=fLabList[i].P();
819 
820  return sum;
821 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::SumChrgPminLab ( double  pmin)

Definition at line 792 of file EventShape.h.

References fCharge, fLabList, fN, i, and P.

Referenced by softtrigger_kin5().

793 {
794  double sum=0;
795  for (int i=0;i<fN;++i)
796  if (fCharge[i]!=0 && fLabList[i].P()>pmin) sum+=fLabList[i].P();
797 
798  return sum;
799 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
double EventShape::SumNeutEmaxCms ( double  emax)

Definition at line 779 of file EventShape.h.

References fCharge, fCmsList, fLabList, fN, and i.

780 {
781  double sum=0;
782  for (int i=0;i<fN;++i)
783  if (fCharge[i]==0 && fCmsList[i].E()<emax) sum+=fLabList[i].E();
784 
785  return sum;
786 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::SumNeutEmaxLab ( double  emax)

Definition at line 757 of file EventShape.h.

References fCharge, fLabList, fN, and i.

758 {
759  double sum=0;
760  for (int i=0;i<fN;++i)
761  if (fCharge[i]==0 && fLabList[i].E()<emax) sum+=fLabList[i].E();
762 
763  return sum;
764 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
double EventShape::SumNeutEminCms ( double  emin)

Definition at line 768 of file EventShape.h.

References fCharge, fCmsList, fLabList, fN, and i.

Referenced by softtrigger_kin5().

769 {
770  double sum=0;
771  for (int i=0;i<fN;++i)
772  if (fCharge[i]==0 && fCmsList[i].E()>emin) sum+=fLabList[i].E();
773 
774  return sum;
775 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::SumNeutEminLab ( double  emin)

Definition at line 746 of file EventShape.h.

References fCharge, fLabList, fN, and i.

Referenced by softtrigger_kin5().

747 {
748  double sum=0;
749  for (int i=0;i<fN;++i)
750  if (fCharge[i]==0 && fLabList[i].E()>emin) sum+=fLabList[i].E();
751 
752  return sum;
753 }
Int_t i
Definition: run_full.C:25
std::vector< int > fCharge
Definition: EventShape.h:110
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
double EventShape::SumPmaxCms ( double  pmax)

Definition at line 687 of file EventShape.h.

References fCmsList, fLabList, fN, i, and P.

688 {
689  double sum=0;
690  for (int i=0;i<fN;++i)
691  if (fCmsList[i].P()<pmax) sum+=fLabList[i].P();
692 
693  return sum;
694 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::SumPmaxLab ( double  pmax)

Definition at line 665 of file EventShape.h.

References fLabList, fN, i, and P.

666 {
667  double sum=0;
668  for (int i=0;i<fN;++i)
669  if (fLabList[i].P()<pmax) sum+=fLabList[i].P();
670 
671  return sum;
672 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
double EventShape::SumPminCms ( double  pmin)

Definition at line 676 of file EventShape.h.

References fCmsList, fLabList, fN, i, and P.

677 {
678  double sum=0;
679  for (int i=0;i<fN;++i)
680  if (fCmsList[i].P()>pmin) sum+=fLabList[i].P();
681 
682  return sum;
683 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::SumPminLab ( double  pmin)

Definition at line 654 of file EventShape.h.

References fLabList, fN, i, and P.

655 {
656  double sum=0;
657  for (int i=0;i<fN;++i)
658  if (fLabList[i].P()>pmin) sum+=fLabList[i].P();
659 
660  return sum;
661 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
GeV c P
double EventShape::SumPtmaxCms ( double  ptmax)

Definition at line 733 of file EventShape.h.

References fCmsList, fLabList, fN, and i.

734 {
735  double sum=0;
736  for (int i=0;i<fN;++i)
737  if (fCmsList[i].Pt()<ptmax) sum+=fLabList[i].Pt();
738 
739  return sum;
740 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::SumPtmaxLab ( double  ptmax)

Definition at line 711 of file EventShape.h.

References fLabList, fN, and i.

712 {
713  double sum=0;
714  for (int i=0;i<fN;++i)
715  if (fLabList[i].Pt()<ptmax) sum+=fLabList[i].Pt();
716 
717  return sum;
718 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
double EventShape::SumPtminCms ( double  ptmin)

Definition at line 722 of file EventShape.h.

References fCmsList, fLabList, fN, and i.

723 {
724  double sum=0;
725  for (int i=0;i<fN;++i)
726  if (fCmsList[i].Pt()>ptmin) sum+=fLabList[i].Pt();
727 
728  return sum;
729 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
double EventShape::SumPtminLab ( double  ptmin)

Definition at line 700 of file EventShape.h.

References fLabList, fN, and i.

Referenced by softtrigger_kin5().

701 {
702  double sum=0;
703  for (int i=0;i<fN;++i)
704  if (fLabList[i].Pt()>ptmin) sum+=fLabList[i].Pt();
705 
706  return sum;
707 }
Int_t i
Definition: run_full.C:25
std::vector< TLorentzVector > fLabList
Definition: EventShape.h:108
double EventShape::Thrust ( )

Definition at line 324 of file EventShape.h.

References Eps(), fabs(), fCmsList, fN, fthr, fThrVect, and i.

Referenced by softtrigger_kin5(), ThrustVector(), and writeTuple().

325 {
326  // did we already compute?
327  if (fthr>-1.) return fthr;
328 
329  TVector3 n0(0,0,0);
330  if( fN==0 ) return -1.;
331 
332  int i,j;
333  double pmax=0;
334 
335  //find starting vector as maximum momentum vector
336  for (i=0;i<fN;++i)
337  {
338  if (fCmsList[i].Vect().Mag()>pmax)
339  {
340  n0=fCmsList[i].Vect();
341  pmax=fCmsList[i].Vect().Mag();
342  }
343  }
344 
345  TVector3 nnew(0,0,0);
346 
347  // find thrust axis (10 iterations)
348  for (i=0;i<10;++i)
349  {
350  for (j=0;j<fN;++j)
351  nnew += Eps(n0,fCmsList[j].Vect()) * fCmsList[j].Vect();
352 
353  n0=nnew.Unit();
354  }
355 
356  fThrVect = n0;
357 
358  double thr=0, sum=0;
359  for (i=0;i<fN;++i)
360  {
361  thr += fabs(fThrVect.Dot(fCmsList[i].Vect()));
362  sum += fCmsList[i].Vect().Mag();
363  }
364 
365  fthr = thr/sum;
366 
367  return fthr;
368 }
double fthr
Definition: EventShape.h:144
Int_t i
Definition: run_full.C:25
TVector3 fThrVect
Definition: EventShape.h:145
double Eps(const TVector3 v1, const TVector3 v2)
Definition: EventShape.h:105
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
std::vector< TLorentzVector > fCmsList
Definition: EventShape.h:109
TVector3 EventShape::ThrustVector ( )

Definition at line 372 of file EventShape.h.

References fthr, fThrVect, and Thrust().

373 {
374  if (fthr<0.) Thrust();
375 
376  return fThrVect;
377 }
double fthr
Definition: EventShape.h:144
TVector3 fThrVect
Definition: EventShape.h:145
double Thrust()
Definition: EventShape.h:324

Member Data Documentation

double EventShape::fapl
private

Definition at line 137 of file EventShape.h.

Referenced by Aplanarity(), and ComputeSphericity().

TVector3 EventShape::fBoost
private

Definition at line 146 of file EventShape.h.

Referenced by EventShape().

std::vector<int> EventShape::fCharge
private
double EventShape::fchrgpsumcms
private

Definition at line 134 of file EventShape.h.

Referenced by ChrgPSumCms(), and EventShape().

double EventShape::fchrgpsumlab
private

Definition at line 128 of file EventShape.h.

Referenced by ChrgPSumLab(), and EventShape().

double EventShape::fchrgptsumcms
private

Definition at line 133 of file EventShape.h.

Referenced by ChrgPtSumCms(), and EventShape().

double EventShape::fchrgptsumlab
private

Definition at line 127 of file EventShape.h.

Referenced by ChrgPtSumLab(), and EventShape().

double EventShape::fcir
private

Definition at line 139 of file EventShape.h.

Referenced by Circularity(), and ComputeSphericity().

std::vector<TLorentzVector> EventShape::fCmsList
private
double EventShape::fFWmom[FWMAX]
private

Definition at line 141 of file EventShape.h.

Referenced by EventShape(), and FoxWolfMomH().

bool EventShape::fFWready
private

Definition at line 142 of file EventShape.h.

Referenced by FoxWolfMomH().

std::vector<TLorentzVector> EventShape::fLabList
private
int EventShape::fN
private
int EventShape::fnChrg
private

Definition at line 112 of file EventShape.h.

Referenced by EventShape(), and NCharged().

double EventShape::fneutesumcms
private

Definition at line 132 of file EventShape.h.

Referenced by EventShape(), and NeutESumCms().

double EventShape::fneutesumlab
private

Definition at line 126 of file EventShape.h.

Referenced by EventShape(), and NeutESumLab().

double EventShape::fneutetsumcms
private

Definition at line 131 of file EventShape.h.

Referenced by EventShape(), and NeutEtSumCms().

double EventShape::fneutetsumlab
private

Definition at line 125 of file EventShape.h.

Referenced by EventShape(), and NeutEtSumLab().

int EventShape::fnNeut
private

Definition at line 113 of file EventShape.h.

Referenced by EventShape(), and NNeutral().

double EventShape::fpla
private

Definition at line 138 of file EventShape.h.

Referenced by ComputeSphericity(), and Planarity().

double EventShape::fpmaxcms
private

Definition at line 117 of file EventShape.h.

Referenced by EventShape(), and PmaxCms().

double EventShape::fpmaxlab
private

Definition at line 116 of file EventShape.h.

Referenced by EventShape(), and PmaxLab().

double EventShape::fpmincms
private

Definition at line 119 of file EventShape.h.

Referenced by EventShape(), and PminCms().

double EventShape::fpminlab
private

Definition at line 118 of file EventShape.h.

Referenced by EventShape(), and PminLab().

double EventShape::fprapmax
private

Definition at line 122 of file EventShape.h.

Referenced by EventShape(), and PRapmax().

double EventShape::fptmax
private

Definition at line 120 of file EventShape.h.

Referenced by EventShape(), and Ptmax().

double EventShape::fptmin
private

Definition at line 121 of file EventShape.h.

Referenced by EventShape(), and Ptmin().

double EventShape::fptsumcms
private

Definition at line 130 of file EventShape.h.

Referenced by EventShape(), and PtSumCms().

double EventShape::fptsumlab
private

Definition at line 124 of file EventShape.h.

Referenced by EventShape(), and PtSumLab().

double EventShape::fsph
private

Definition at line 136 of file EventShape.h.

Referenced by ComputeSphericity(), and Sphericity().

double EventShape::fthr
private

Definition at line 144 of file EventShape.h.

Referenced by Thrust(), and ThrustVector().

TVector3 EventShape::fThrVect
private

Definition at line 145 of file EventShape.h.

Referenced by EventShape(), Thrust(), and ThrustVector().


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