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

#include <KFPTopoReconstructor.h>

Classes

struct  KFParticleCluster
 

Public Member Functions

 KFPTopoReconstructor ()
 
 ~KFPTopoReconstructor ()
 
void Init (vector< KFPTrack > &tracks, int nParticles)
 
void ReconstructPrimVertex ()
 
KFParticleGetPrimVertex ()
 Accessors. More...
 
KFParticleGetParticle (int i)
 

Private Member Functions

KFPTopoReconstructoroperator= (KFPTopoReconstructor &)
 
 KFPTopoReconstructor (KFPTopoReconstructor &)
 
void FindPrimaryClusters ()
 

Private Attributes

vector< KFParticlefParticles
 
int fNParticles
 
vector< KFParticleClusterfClusters
 
vector< KFVertexfPrimVertices
 

Detailed Description

Definition at line 17 of file KFPTopoReconstructor.h.

Constructor & Destructor Documentation

KFPTopoReconstructor::KFPTopoReconstructor ( )
inline

Definition at line 19 of file KFPTopoReconstructor.h.

19 {};
KFPTopoReconstructor::~KFPTopoReconstructor ( )
inline

Definition at line 20 of file KFPTopoReconstructor.h.

20 {};
KFPTopoReconstructor::KFPTopoReconstructor ( KFPTopoReconstructor )
private

Member Function Documentation

void KFPTopoReconstructor::FindPrimaryClusters ( )
private

Definition at line 53 of file KFPTopoReconstructor.cxx.

References KFPTopoReconstructor::KFParticleCluster::fC, KFPTopoReconstructor::KFParticleCluster::fCluster, fClusters, fNParticles, KFPTopoReconstructor::KFParticleCluster::fP, fParticles, i, and sqrt().

Referenced by ReconstructPrimVertex().

54 {
55  // The function finds a set of clusters of tracks.
56  // Tracks are assumed to be transported to the beam line.
57 
58  vector<unsigned short int> notUsedTracks(fNParticles);
59  for(int iTr=0; iTr<fNParticles; iTr++)
60  notUsedTracks[iTr] = iTr;
61 
62  while(notUsedTracks.size()>0)
63  {
64  short int bestTrack = 0;
65  float bestWeight = -1.f;
66 
67  vector<unsigned short int> notUsedTracksNew;
68 
69  for(unsigned short int iTr = 0; iTr < notUsedTracks.size(); iTr++)
70  {
71  unsigned short int &curTrack = notUsedTracks[iTr];
72  float dX = fParticles[curTrack].CovarianceMatrix()[0];
73  float dY = fParticles[curTrack].CovarianceMatrix()[2];
74  float dZ = fParticles[curTrack].CovarianceMatrix()[5];
75 
76  float weight = 1.f/sqrt(dX + dY + dZ);
77  if (weight > bestWeight)
78  {
79  bestWeight = weight;
80  bestTrack = curTrack;
81  }
82  }
83 // std::cout << "best " << bestWeight << " " << bestTrack << std::endl;
84 
85  KFParticleCluster cluster;
86 
87  const double *rBest = fParticles[bestTrack].Parameters();
88  const double *covBest = fParticles[bestTrack].CovarianceMatrix();
89 
90  float rVertex[3] = {0.f};
91  float covVertex[6] = {0.f};
92  float weightVertex = 0.f;
93 /*std::cout << "tatata " << std::endl;*/
94  for(unsigned short int iTr = 0; iTr < notUsedTracks.size(); iTr++)
95  {
96  unsigned short int &curTrack = notUsedTracks[iTr];
97 
98  double dr[3] = {rBest[0] - fParticles[curTrack].X(),
99  rBest[1] - fParticles[curTrack].Y(),
100  rBest[2] - fParticles[curTrack].Z() };
101 
102  double cov[6] = {covBest[0] + fParticles[curTrack].CovarianceMatrix()[0],
103  covBest[1] + fParticles[curTrack].CovarianceMatrix()[1],
104  covBest[2] + fParticles[curTrack].CovarianceMatrix()[2],
105  covBest[3] + fParticles[curTrack].CovarianceMatrix()[3],
106  covBest[4] + fParticles[curTrack].CovarianceMatrix()[4],
107  covBest[5] + fParticles[curTrack].CovarianceMatrix()[5] };
108 
109  float dr2 = static_cast<float>(dr[0]*dr[0] + dr[1]*dr[1] + dr[2]*dr[2]);
110  float drError2 = static_cast<float>(
111  dr[0]* (cov[0]* dr[0] + cov[1]* dr[1] + cov[3]* dr[2]) +
112  dr[1]* (cov[1]* dr[0] + cov[2]* dr[1] + cov[4]* dr[2]) +
113  dr[2]* (cov[3]* dr[0] + cov[4]* dr[1] + cov[5]* dr[2])
114  );
115 // std::cout << "dr " << sqrt(dr2) << " drError " << sqrt(drError2/dr2) << std::endl;
116 
117  if( dr2*dr2 <= (10*10 * drError2) )
118  {
119  float dX = fParticles[curTrack].CovarianceMatrix()[0];
120  float dY = fParticles[curTrack].CovarianceMatrix()[2];
121  float dZ = fParticles[curTrack].CovarianceMatrix()[5];
122 
123  float weight = 1.f/sqrt(dX + dY + dZ);
124 
125  for(int iP=0; iP<3; iP++)
126  rVertex[iP] += weight * fParticles[curTrack].Parameters()[iP];
127 
128  for(int iC=0; iC<6; iC++)
129  covVertex[iC] += weight * weight *fParticles[curTrack].CovarianceMatrix()[iC];
130 
131  weightVertex += weight;
132  cluster.fCluster.push_back(curTrack);
133  }
134  else
135  {
136 // int iMC = PndFTSCAPerformance::Instance().GetSubPerformance("Global Performance")->GetRecoData()[fParticles[curTrack].DaughterIds()[0]].GetMCTrackId();
137 // int MCPDG = -1;
138 // int PDG = -1;
139 //
140 // if(iMC>0)
141 // {
142 //
143 // int motherID = (*PndFTSCAPerformance::Instance().GetSubPerformance("Global Performance")->fMCTracks)[iMC].MotherId();
144 // PDG = (*PndFTSCAPerformance::Instance().GetSubPerformance("Global Performance")->fMCTracks)[iMC].PDG();
145 //
146 // if(motherID > 0)
147 // MCPDG = (*PndFTSCAPerformance::Instance().GetSubPerformance("Global Performance")->fMCTracks)[motherID].PDG();
148 // }
149 //
150 // std::cout << "Bad PDG " << PDG << " MC PDG "<< MCPDG << std::endl;
151  notUsedTracksNew.push_back(curTrack);
152  }
153  }
154 
155  notUsedTracks = notUsedTracksNew;
156 
157  if(cluster.fCluster.size() < 2) continue;
158  if((cluster.fCluster.size() < 3) && (fNParticles>3)) continue;
159 
160  for(int iP=0; iP<3; iP++)
161  cluster.fP[iP] = rVertex[iP]/weightVertex;
162 
163  for(int iC=0; iC<6; iC++)
164  cluster.fC[iC] = covVertex[iC]/(weightVertex*weightVertex);
165 // std::cout << "clust " << cluster.fP[0] << " " << cluster.fP[1] << " " << cluster.fP[2] << " " << cluster.fCluster.size() << " " << notUsedTracksNew.size() << std::endl;
166 // std::cout << "Cov " << cluster.fC[0] << " "
167 // << cluster.fC[1] << " " << cluster.fC[2] << " "
168 // << cluster.fC[3] << " " << cluster.fC[4] << " " << cluster.fC[5] << std::endl;
169 // int ui;
170 // std::cin >> ui;
171 
172  fClusters.push_back(cluster);
173  }
174 // if(fClusters.size()>1)
175 // {
176 // float dx = fClusters[0].fP[0] - fClusters[1].fP[0];
177 // float dy = fClusters[0].fP[1] - fClusters[1].fP[1];
178 // float dz = fClusters[0].fP[2] - fClusters[1].fP[2];
179 //
180 // float dr[3] = {dx, dy, dz};
181 // float cov[6] = {fClusters[0].fC[0] + fClusters[1].fC[0],
182 // fClusters[0].fC[1] + fClusters[1].fC[1],
183 // fClusters[0].fC[2] + fClusters[1].fC[2],
184 // fClusters[0].fC[3] + fClusters[1].fC[3],
185 // fClusters[0].fC[4] + fClusters[1].fC[4],
186 // fClusters[0].fC[5] + fClusters[1].fC[5] };
187 // float dr2 = dr[0]*dr[0] + dr[1]*dr[1] + dr[2]*dr[2];
188 // float drError2 = dr[0]* (cov[0]* dr[0] + cov[1]* dr[1] + cov[3]* dr[2]) +
189 // dr[1]* (cov[1]* dr[0] + cov[2]* dr[1] + cov[4]* dr[2]) +
190 // dr[2]* (cov[3]* dr[0] + cov[4]* dr[1] + cov[5]* dr[2]);
191 // drError2 /= dr2;
192 // std::cout << " dr2 " << dr2 << " Err " << drError2 << std::endl;
193 // }
194 // int ui;
195 // std::cin >> ui;
196 
197  static int nVert[10]={0};
198  if(fClusters.size()<10)
199  nVert[fClusters.size()]++;
200  std::cout << "N Vert ";
201  for(int i=0; i<10; i++)
202  std::cout << i << ": " << nVert[i] << " ";
203  std::cout << std::endl;
204 
205  static int nPart[10]={0};
206  if(fNParticles<10)
207  nPart[fNParticles]++;
208  std::cout << "N Part ";
209  for(int i=0; i<10; i++)
210  std::cout << i << ": " << nPart[i] << " ";
211  std::cout << std::endl;
212 }
Int_t i
Definition: run_full.C:25
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
vector< KFParticleCluster > fClusters
vector< KFParticle > fParticles
KFParticle& KFPTopoReconstructor::GetParticle ( int  i)
inline

Definition at line 28 of file KFPTopoReconstructor.h.

References fNParticles, fParticles, and i.

28 { assert( i < fNParticles ); return fParticles[i]; };
Int_t i
Definition: run_full.C:25
vector< KFParticle > fParticles
KFParticle& KFPTopoReconstructor::GetPrimVertex ( )
inline

Accessors.

Definition at line 27 of file KFPTopoReconstructor.h.

References fPrimVertices.

Referenced by PndFTSTopoReconstructor::GetPrimVertex().

27 { assert( fPrimVertices.size() > 0 ); return fPrimVertices[0]; };
vector< KFVertex > fPrimVertices
void KFPTopoReconstructor::Init ( vector< KFPTrack > &  tracks,
int  nParticles 
)

Definition at line 30 of file KFPTopoReconstructor.cxx.

References fNParticles, and fParticles.

Referenced by PndFTSTopoReconstructor::Init().

31 {
32  // copy tracks in particles.
33  fNParticles = nParticles;
34  fParticles.resize(fNParticles);
35  for ( int iTr = 0; iTr < fNParticles; iTr++ ) {
36  // { // check cov matrix TODO
37  // bool ok = true;
38  // int k = 0;
39  // for (int i = 0; i < 6; i++) {
40  // for (int j = 0; j <= i; j++, k++) {
41  // ok &= finite( KFPCov[k] );
42  // }
43  // ok &= ( KFPCov[k-1] > 0 );
44  // }
45  // if (!ok) continue;
46  // }
47 
48  fParticles[iTr] = KFParticle( tracks[iTr], 211 ); // pi+ // TODO include PDG in KFPTrack?
49  fParticles[iTr].AddDaughterId(tracks[iTr].Id());
50  }
51 } // void KFPTopoReconstructor::Init
vector< KFParticle > fParticles
KFPTopoReconstructor& KFPTopoReconstructor::operator= ( KFPTopoReconstructor )
private
void KFPTopoReconstructor::ReconstructPrimVertex ( )

Definition at line 214 of file KFPTopoReconstructor.cxx.

References KFVertex::ConstructPrimaryVertex(), f, FindPrimaryClusters(), fNParticles, fParticles, fPrimVertices, KFParticle::GetNDF(), p, KFPVertex::SetChi2(), KFPVertex::SetCovarianceMatrix(), KFPVertex::SetNContributors(), and KFPVertex::SetXYZ().

Referenced by PndFTSTopoReconstructor::ReconstructPrimVertex().

215 {
216 
218 
219  // select particles for Primary vertex reco
220  const KFParticle **pParticles = new const KFParticle*[fNParticles]; // tmp array
221 
222  int nPrimCand = 0;
223  for ( int iP = 0; iP < fNParticles; iP++ ) {
224  KFParticle &p = fParticles[iP];
225  // if ( p->GetPt() < 1.f ) continue;
226  // if ( p->GetNDF() < 40-5 ) continue; // at least 20 hits in track // dbg
227  // if ( p->GetChi2() > 2.f * p->GetNDF() ) continue; // dbg
228 
229  pParticles[nPrimCand++] = &p;
230  }
231 
232  // find prim vertex
233  KFVertex primVtx;
234 
235  if(fNParticles>1)
236  {
237  bool *vFlags = new bool[fNParticles]; // flags returned by the vertex finder
238  primVtx.ConstructPrimaryVertex( pParticles, 2, vFlags, 10.f );
239  if( primVtx.GetNDF() >= 1 )
240  fPrimVertices.push_back(primVtx);
241 
242  if(vFlags) delete [] vFlags;
243  }
244 
245 
246  if ( fPrimVertices.size() == 0 ) { // fill prim vertex by dummy values
247  KFPVertex primVtx_tmp;
248  primVtx_tmp.SetXYZ(0, 0, 0);
249  primVtx_tmp.SetCovarianceMatrix( 0, 0, 0, 0, 0, 0 );
250  primVtx_tmp.SetNContributors(0);
251  primVtx_tmp.SetChi2(-100);
252 
253  fPrimVertices.push_back(primVtx_tmp);
254  }
255 
256  if(pParticles) delete [] pParticles;
257  // if ( fPrimVertices.size() > 0 ) { // dbg
258  // std::cout << std::endl << " PrimVertex " << std::endl;
259  // std::cout << fPrimVertices[0] << std::endl << std::endl;
260  // }
261 } // void KFPTopoReconstructor::Run()
Int_t GetNDF() const
Definition: KFParticle.h:497
void SetChi2(float chi)
Definition: KFPVertex.h:44
void SetNContributors(int nc)
Definition: KFPVertex.h:46
Double_t p
Definition: anasim.C:58
void SetXYZ(float *position)
Definition: KFPVertex.h:39
TFile * f
Definition: bump_analys.C:12
void ConstructPrimaryVertex(const KFParticle *vDaughters[], int NDaughters, Bool_t vtxFlag[], float ChiCut=3.5)
Definition: KFVertex.cxx:67
vector< KFVertex > fPrimVertices
void SetCovarianceMatrix(float *C)
Definition: KFPVertex.h:48
vector< KFParticle > fParticles

Member Data Documentation

vector< KFParticleCluster > KFPTopoReconstructor::fClusters
private

Definition at line 45 of file KFPTopoReconstructor.h.

Referenced by FindPrimaryClusters().

int KFPTopoReconstructor::fNParticles
private
vector<KFParticle> KFPTopoReconstructor::fParticles
private
vector<KFVertex> KFPTopoReconstructor::fPrimVertices
private

Definition at line 46 of file KFPTopoReconstructor.h.

Referenced by GetPrimVertex(), and ReconstructPrimVertex().


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