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

#include <PndMvdRiemannVertexFinderTask.h>

Inheritance diagram for PndMvdRiemannVertexFinderTask:

Public Member Functions

 PndMvdRiemannVertexFinderTask ()
 
virtual ~PndMvdRiemannVertexFinderTask ()
 
 PndMvdRiemannVertexFinderTask (const PndMvdRiemannVertexFinderTask &)=delete
 
PndMvdRiemannVertexFinderTaskoperator= (const PndMvdRiemannVertexFinderTask &)=delete
 
virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
virtual void Exec (Option_t *opt)
 
virtual void FinishEvent ()
 
void SetVerbose (Int_t verbose)
 
void SetVertexCut (double cut)
 

Public Attributes

TH1F * delta
 
TH1F * wrongV
 
std::pair< double, double > eff
 
std::pair< double, double > ghosts
 

Private Member Functions

bool CheckRecoTrack (PndTrackCand *cand, PndMCTrack *myTrack)
 
bool CheckVertex (std::vector< int > Combination, std::vector< std::pair< int, int > > PairCand)
 
bool CheckTwoCands (int first, int second)
 
int FoundCandInMCCands (int candN)
 
void refit (std::vector< int > &CheckedCand)
 
void FindVertex (std::vector< int > CheckedCand, std::vector< std::pair< int, int > > &PairCand, std::vector< std::pair< int, int > > &TrueMCCand, std::vector< std::pair< int, int > > &FalseMCCand, std::vector< std::pair< int, int > > &MCCand, int &MaxIndex)
 
void CalcEfficiency (std::vector< std::pair< int, int > > TrueMCCand, std::vector< std::pair< int, int > > FalseMCCand, std::vector< std::pair< int, int > > MCCand)
 
void Register ()
 
void Reset ()
 
void ProduceHits ()
 
 ClassDef (PndMvdRiemannVertexFinderTask, 1)
 

Private Attributes

TString fHitBranch
 
TString fHitBranch2
 
TString fTrackBranch
 
TString fIdealTrackCandBranch
 
TString fMCTrackBranch
 
int fEventNr
 
int fVerbose
 
double fVertexCut
 
TClonesArray * fHitArray
 
TClonesArray * fHitArray2
 
TClonesArray * fTrackCandArray
 
TClonesArray * fTrackArray
 
TClonesArray * fIdealTrackCandArray
 
TClonesArray * fMCTrackArray
 
TClonesArray * fVertex
 
TClonesArray * fMCVertex
 

Detailed Description

Definition at line 13 of file PndMvdRiemannVertexFinderTask.h.

Constructor & Destructor Documentation

PndMvdRiemannVertexFinderTask::PndMvdRiemannVertexFinderTask ( )

Definition at line 34 of file PndMvdRiemannVertexFinderTask.cxx.

34  : FairTask("MVD Riemann VERTEX Finder"),
35  delta(NULL),
36  wrongV(NULL),
37  eff(),
38  ghosts(),
39  fHitBranch("MVDHitsPixel"),
40  fHitBranch2("MVDHitsStrip"),
41  fTrackBranch("MVDRiemannTrackCand"),
42  fIdealTrackCandBranch("MVDIdealTrackCand"),
43  fMCTrackBranch("MCTrack"),
44  fEventNr(0),
45  fVerbose(0),
46  fVertexCut(0.),
47  fHitArray(NULL),
48  fHitArray2(NULL),
49  fTrackCandArray(NULL),
50  fTrackArray(NULL),
52  fMCTrackArray(NULL),
53  fVertex(NULL),
54  fMCVertex(NULL)
55 {
56 }
PndMvdRiemannVertexFinderTask::~PndMvdRiemannVertexFinderTask ( )
virtual

Definition at line 58 of file PndMvdRiemannVertexFinderTask.cxx.

59 {
60 }
PndMvdRiemannVertexFinderTask::PndMvdRiemannVertexFinderTask ( const PndMvdRiemannVertexFinderTask )
delete

Member Function Documentation

void PndMvdRiemannVertexFinderTask::CalcEfficiency ( std::vector< std::pair< int, int > >  TrueMCCand,
std::vector< std::pair< int, int > >  FalseMCCand,
std::vector< std::pair< int, int > >  MCCand 
)
private

Definition at line 172 of file PndMvdRiemannVertexFinderTask.cxx.

References counter, eff, ghosts, and i.

Referenced by Exec().

174 {
175  int counter=0;
176  for(unsigned int i=0;i<MCCand.size();i++){
177  int MCfirst=MCCand[i].first;
178  int MCsecond=MCCand[i].second;
179  for(unsigned int j=0;j<TrueMCCand.size();j++){
180  int first=TrueMCCand[j].first;
181  int second=TrueMCCand[j].second;
182  if ((MCfirst==first && MCsecond==second) or (MCfirst==second && MCsecond==first)){
183  counter++;
184  break;
185  }
186  }
187  }
188  eff.first+=counter;
189  eff.second+=MCCand.size();
190  ghosts.first+=FalseMCCand.size();
191  ghosts.second+=MCCand.size();
192 // std::cout<<"eff= "<<((double)counter)/((double)MCCand.size())<<std::endl;
193 // std::cout<<"ghosts= "<<((double)FalseMCCand.size())/((double)MCCand.size())<<std::endl;
194 
195 }
Int_t i
Definition: run_full.C:25
int counter
Definition: ZeeAnalysis.C:59
bool PndMvdRiemannVertexFinderTask::CheckRecoTrack ( PndTrackCand cand,
PndMCTrack myTrack 
)
private

Definition at line 377 of file PndMvdRiemannVertexFinderTask.cxx.

References a, count, fHitArray, fHitArray2, fHitBranch, fHitBranch2, PndTrackCandHit::GetDetId(), PndTrackCandHit::GetHitId(), PndMCTrack::GetMotherID(), PndTrackCand::GetNHits(), PndSdsHit::GetPosition(), PndTrackCand::GetSortedHit(), and i.

Referenced by FindVertex().

378 {
379  if (/*myTrack->GetPdgCode()==211 &&*/ (cand->GetNHits()>2) && (myTrack->GetMotherID()==-1)){
380  int count=0;
381  unsigned int detIDi, hitIDi;
382  unsigned int detIDj, hitIDj;
383  for(unsigned int i=0;i<cand->GetNHits();i++){
384  detIDi=cand->GetSortedHit(i).GetDetId();
385  hitIDi=cand->GetSortedHit(i).GetHitId();
386  PndSdsHit *pointI;
387  if ((int)detIDi == FairRootManager::Instance()->GetBranchId(fHitBranch))
388  pointI = (PndSdsHit*)fHitArray->At(hitIDi);
389  else if ((int)detIDi == FairRootManager::Instance()->GetBranchId(fHitBranch2))
390  pointI = (PndSdsHit*)fHitArray2->At(hitIDi);
391  else pointI = 0;
392 
393  for(unsigned int j=0;j<cand->GetNHits();j++){
394  detIDj=cand->GetSortedHit(j).GetDetId();
395  hitIDj=cand->GetSortedHit(j).GetHitId();
396  PndSdsHit *pointJ;
397  if ((int)detIDj == FairRootManager::Instance()->GetBranchId(fHitBranch))
398  pointJ = (PndSdsHit*)fHitArray->At(hitIDj);
399  else if ((int)detIDj == FairRootManager::Instance()->GetBranchId(fHitBranch2))
400  pointJ = (PndSdsHit*)fHitArray2->At(hitIDj);
401  else pointJ = 0;
402 
403  if ((pointI!=0) && (pointJ!=0) && (i!=j)){
404  TVector3 a=pointI->GetPosition() - pointJ->GetPosition();
405  if (a.Mag()<1){
406  count++;
407  }
408  }
409  }
410  }
411  if ((cand->GetNHits()-count/2)<3){
412 // std::cout <<"Less then 3 base points in RecoTrack"<<" "<<cand->getNHits()<<" "<<count <<std::endl;
413  return false;
414  }
415  else return true;
416 
417 
418  }
419  else
420  return false;
421 }
Int_t i
Definition: run_full.C:25
TVector3 GetPosition() const
Definition: PndSdsHit.h:93
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
Int_t a
Definition: anaLmdDigi.C:126
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
int count
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
Int_t GetHitId() const
Int_t GetDetId() const
bool PndMvdRiemannVertexFinderTask::CheckTwoCands ( int  first,
int  second 
)
private

Definition at line 353 of file PndMvdRiemannVertexFinderTask.cxx.

References counter, fTrackCandArray, PndTrackCandHit::GetDetId(), PndTrackCandHit::GetHitId(), PndTrackCand::GetNHits(), PndTrackCand::GetSortedHit(), h1, h2, and i.

Referenced by FindVertex().

354 {
355  PndTrackCand* cand1=(PndTrackCand*)fTrackCandArray->At(first);
356  PndTrackCand* cand2=(PndTrackCand*)fTrackCandArray->At(second);
357  int counter=0;
358  for(unsigned int i=0;i<cand1->GetNHits();i++){
359  unsigned int d1,h1;
360  d1=cand1->GetSortedHit(i).GetDetId();
361  h1=cand1->GetSortedHit(i).GetHitId();
362  for(unsigned int j=0;j<cand2->GetNHits();j++){
363  unsigned int d2,h2;
364  d2=cand2->GetSortedHit(j).GetDetId();
365  h2=cand2->GetSortedHit(j).GetHitId();
366  if(d1==d2 && h1==h2){
367  counter++;
368  }
369  }
370  }
371  if (counter<3)
372  return true;
373  else
374  return false;
375 }
Int_t i
Definition: run_full.C:25
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
int counter
Definition: ZeeAnalysis.C:59
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
Int_t GetHitId() const
Int_t GetDetId() const
bool PndMvdRiemannVertexFinderTask::CheckVertex ( std::vector< int >  Combination,
std::vector< std::pair< int, int > >  PairCand 
)
private

Definition at line 335 of file PndMvdRiemannVertexFinderTask.cxx.

References i.

336 {
337  int count1=0;
338  for(unsigned int i=0;i<PairCand.size();i++){
339  int ind1=PairCand[i].first,ind2=PairCand[i].second;
340  int count2=0;
341  for(unsigned int j=0;j<Combination.size();j++){
342  if ((ind1==Combination[j]) or (ind2==Combination[j]))
343  count2++;
344  }
345  if (count2==2)
346  count1++;
347  }
348 
349  if (count1==3)
350  return true;
351  else return false;
352 }
Int_t i
Definition: run_full.C:25
PndMvdRiemannVertexFinderTask::ClassDef ( PndMvdRiemannVertexFinderTask  ,
 
)
private
void PndMvdRiemannVertexFinderTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 152 of file PndMvdRiemannVertexFinderTask.cxx.

References CalcEfficiency(), fEventNr, FindVertex(), fTrackCandArray, and refit().

153 {
154  // Reset output array
155  if ( ! fTrackCandArray )
156  Fatal("Exec", "No trackCandArray");
157  std::vector<int> CheckedCand;
158  refit(CheckedCand);
159 
160 
161 
162  std::vector< std::pair <int,int> > PairCand;
163  std::vector< std::pair <int,int> > TrueMCCand;
164  std::vector< std::pair <int,int> > FalseMCCand;
165  std::vector< std::pair <int,int> > MCCand;
166  int MaxIndex=0;
167  FindVertex(CheckedCand,PairCand,TrueMCCand,FalseMCCand,MCCand ,MaxIndex);
168  CalcEfficiency(TrueMCCand, FalseMCCand, MCCand);
169 
170  std::cout << "Done event : " << fEventNr++ << std::endl;
171 }
void refit(std::vector< int > &CheckedCand)
void FindVertex(std::vector< int > CheckedCand, std::vector< std::pair< int, int > > &PairCand, std::vector< std::pair< int, int > > &TrueMCCand, std::vector< std::pair< int, int > > &FalseMCCand, std::vector< std::pair< int, int > > &MCCand, int &MaxIndex)
void CalcEfficiency(std::vector< std::pair< int, int > > TrueMCCand, std::vector< std::pair< int, int > > FalseMCCand, std::vector< std::pair< int, int > > MCCand)
void PndMvdRiemannVertexFinderTask::FindVertex ( std::vector< int >  CheckedCand,
std::vector< std::pair< int, int > > &  PairCand,
std::vector< std::pair< int, int > > &  TrueMCCand,
std::vector< std::pair< int, int > > &  FalseMCCand,
std::vector< std::pair< int, int > > &  MCCand,
int &  MaxIndex 
)
private

for D+D-

Definition at line 197 of file PndMvdRiemannVertexFinderTask.cxx.

References c1, c2, PndRiemannTrack::calcIntersection(), CheckRecoTrack(), CheckTwoCands(), delta, fabs(), fIdealTrackCandArray, fMCTrackArray, FoundCandInMCCands(), fTrackArray, fVertex, PndTrackCand::getMcTrackId(), PndMCTrack::GetPdgCode(), PndMCTrack::GetStartVertex(), i, p1, p2, PndRiemannTrack::SetVerbose(), and wrongV.

Referenced by Exec().

199 {
200  for(int i=0;i<fTrackArray->GetEntriesFast();i++){
202  track1.SetVerbose(-1);
203  for(int j=i+1;j<fTrackArray->GetEntriesFast();j++){
205  track2.SetVerbose(-1);
206  TVector3 p1,p2,av;
207  if ((track2.calcIntersection(track1,p1,p2)==1) && CheckTwoCands(CheckedCand[i],CheckedCand[j] ) ){
208  // std::cout<<"Vertex: "<<i<<" : "<<j<<std::endl;
209  pair <int,int> Pair(i,j);
210  PairCand.push_back(Pair);
211  if (j>MaxIndex) MaxIndex=j;
212  av=p1+p2;
213  av*=0.5;
214  PndSdsMCPoint* Vertex1 = new PndSdsMCPoint(0,0,-1,p1,p1,p1,p1,0,0,0);
215  PndSdsMCPoint* Vertex2 = new PndSdsMCPoint(0,0,-1,p2,p2,p2,p2,0,0,0);
216  new((*fVertex)[fVertex->GetEntriesFast()])PndSdsMCPoint(*Vertex1);
217  new((*fVertex)[fVertex->GetEntriesFast()])PndSdsMCPoint(*Vertex2);
218 
219  int c1=FoundCandInMCCands(CheckedCand[i]);
220  int c2=FoundCandInMCCands(CheckedCand[j]);
221  if (c1>-1 && c2>-1){
222  // std::cout<<c1<<" "<<c2<<std::endl;
224  PndMCTrack* myTrack1 = (PndMCTrack*)fMCTrackArray->At(Cand1->getMcTrackId());
226  PndMCTrack* myTrack2 = (PndMCTrack*)fMCTrackArray->At(Cand2->getMcTrackId());
227  if (((myTrack1->GetPdgCode()==myTrack2->GetPdgCode()) or ((myTrack1->GetPdgCode()*myTrack2->GetPdgCode()<0) && (fabs(myTrack1->GetPdgCode())!=fabs(myTrack2->GetPdgCode()))))){
228  delta->Fill((myTrack1->GetStartVertex()-p1).Mag());
229  delta->Fill((myTrack1->GetStartVertex()-p2).Mag());
230  pair <int,int> Pair2(c1,c2);
231  TrueMCCand.push_back(Pair2);
232  }
233  else{
234  wrongV->Fill((myTrack1->GetStartVertex()-p1).Mag());
235  wrongV->Fill((myTrack1->GetStartVertex()-p2).Mag());
236  pair <int,int> Pair3(c1,c2);
237  FalseMCCand.push_back(Pair3);
238  }
239 
240  }
241  }
242  }
243  }
244 
245  for(int i=0;i<fIdealTrackCandArray->GetEntriesFast();i++){
247  PndMCTrack* myTrack1 = (PndMCTrack*)fMCTrackArray->At(Cand1->getMcTrackId());
248  for(int j=i+1;j<fIdealTrackCandArray->GetEntriesFast();j++){
250  PndMCTrack* myTrack2 = (PndMCTrack*)fMCTrackArray->At(Cand2->getMcTrackId());
251  if (((myTrack1->GetPdgCode()==myTrack2->GetPdgCode()) or ((myTrack1->GetPdgCode()*myTrack2->GetPdgCode()<0) && (fabs(myTrack1->GetPdgCode())!=fabs(myTrack2->GetPdgCode()))))){
252  if(CheckRecoTrack(Cand1,myTrack1) && CheckRecoTrack(Cand2,myTrack2) ){
253  pair <int,int> Pair(i,j);
254  MCCand.push_back(Pair);
255  // std::cout<<"IdealVertex: "<<i<<" : "<<j<<std::endl;
256  }
257  }
258 
259  }
260 
261  }
262 
263 
264 }
int getMcTrackId() const
Definition: PndTrackCand.h:60
Int_t i
Definition: run_full.C:25
bool CheckRecoTrack(PndTrackCand *cand, PndMCTrack *myTrack)
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
c2
Definition: plot_dirc.C:39
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
void SetVerbose(int i)
c1
Definition: plot_dirc.C:35
TPad * p2
Definition: hist-t7.C:117
int calcIntersection(PndRiemannTrack &track, TVector3 &p1, TVector3 &p2)
TPad * p1
Definition: hist-t7.C:116
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
void PndMvdRiemannVertexFinderTask::FinishEvent ( )
virtual
int PndMvdRiemannVertexFinderTask::FoundCandInMCCands ( int  candN)
private

Definition at line 302 of file PndMvdRiemannVertexFinderTask.cxx.

References counter, fIdealTrackCandArray, fTrackCandArray, PndTrackCandHit::GetDetId(), PndTrackCandHit::GetHitId(), PndTrackCand::GetNHits(), PndTrackCand::GetSortedHit(), h1, h2, and i.

Referenced by FindVertex().

303 {
304  PndTrackCand* FoundCand=(PndTrackCand*)fTrackCandArray->At(candN);
305  for(int k=0;k<fIdealTrackCandArray->GetEntriesFast();k++){
307  unsigned int counter=0;
308  bool ZeroPos=false;
309  for(unsigned int i=0;i<FoundCand->GetNHits();i++){
310  unsigned int d1,h1;
311  d1=FoundCand->GetSortedHit(i).GetDetId();
312  h1=FoundCand->GetSortedHit(i).GetHitId();
313  if (d1==0){
314  ZeroPos=true;
315  continue;
316  }
317  for(unsigned int j=0;j<TestCand->GetNHits();j++){
318  unsigned int d2,h2;
319  d2=TestCand->GetSortedHit(j).GetDetId();
320  h2=TestCand->GetSortedHit(j).GetHitId();
321  if(d1==d2 && h1==h2){
322  counter++;
323  }
324  }
325  }
326  if ((counter==FoundCand->GetNHits() && !ZeroPos) or (counter==FoundCand->GetNHits()-1 && ZeroPos)){
327  return k;
328  }
329  }
330  return -1;
331 }
Int_t i
Definition: run_full.C:25
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
int counter
Definition: ZeeAnalysis.C:59
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
Int_t GetHitId() const
Int_t GetDetId() const
InitStatus PndMvdRiemannVertexFinderTask::Init ( )
virtual

Definition at line 88 of file PndMvdRiemannVertexFinderTask.cxx.

References delta, fHitArray, fHitArray2, fHitBranch, fHitBranch2, fIdealTrackCandArray, fIdealTrackCandBranch, fMCTrackArray, fMCTrackBranch, fMCVertex, fTrackArray, fTrackBranch, fTrackCandArray, fVertex, and wrongV.

89 {
90 
91  FairRootManager* ioman = FairRootManager::Instance();
92 
93  if ( ! ioman )
94  {
95  std::cout << "-E- PndMvdRiemannVertexFinderTask::Init: "
96  << "RootManager not instantiated!" << std::endl;
97  return kFATAL;
98  }
99 
100  // Get input array
101  fHitArray = (TClonesArray*) ioman->GetObject(fHitBranch);
102  if ( !fHitArray){
103  std::cout << "-W- PndMvdRiemannVertexFinderTask::Init: " << "No hitArray!" << std::endl;
104  return kERROR;
105  }
106 
107  fHitArray2 = (TClonesArray*) ioman->GetObject(fHitBranch2);
108  if ( !fHitArray2){
109  std::cout << "-W- PndMvdRiemannVertexFinderTask::Init: " << "No hitArray2!" << std::endl;
110  return kERROR;
111  }
112 
113 
114  fTrackCandArray = (TClonesArray*) ioman->GetObject(fTrackBranch);
115  if ( !fTrackCandArray){
116  std::cout << "-W- PndMvdRiemannVertexFinderTask::Init: " << "No tracks!" << std::endl;
117  return kERROR;
118  }
119 
120  fIdealTrackCandArray= (TClonesArray*) ioman->GetObject(fIdealTrackCandBranch);
121  if ( !fIdealTrackCandArray){
122  std::cout << "-W- PndMvdRiemannVertexFinderTask::Init: " << "No ideal tracks!" << std::endl;
123  return kERROR;
124  }
125 
126 
127  fMCTrackArray = (TClonesArray*) ioman->GetObject(fMCTrackBranch);
128  if ( !fMCTrackArray){
129  std::cout << "-W- PndMvdRiemannVertexFinderTask::Init: " << "No MC tracks!" << std::endl;
130  return kERROR;
131  }
132 
133  fTrackArray = new TClonesArray("PndRiemannTrack");
134 
135  delta = new TH1F("delta","delta",1000,0,10);
136  wrongV= new TH1F("wrong","wrong",1000,0,10);
137 
138  fVertex = new TClonesArray("PndSdsMCPoint");
139  ioman->Register("Vertex", "MVD", fVertex, kTRUE);
140 
141  fMCVertex = new TClonesArray("PndSdsMCPoint");
142  ioman->Register("MCVertex", "MVD", fMCVertex, kTRUE);
143 
144 // fRiemannTrackArray = new TClonesArray("PndRiemannTrack");
145 // ioman->Register("MVDRiemannTrack","MVD",fRiemannTrackArray, kTRUE);
146 
147  std::cout << "-I- PndMvdRiemannVertexFinderTask: Initialisation successfull" << std::endl;
148  return kSUCCESS;
149 }
PndMvdRiemannVertexFinderTask& PndMvdRiemannVertexFinderTask::operator= ( const PndMvdRiemannVertexFinderTask )
delete
void PndMvdRiemannVertexFinderTask::ProduceHits ( )
private
void PndMvdRiemannVertexFinderTask::refit ( std::vector< int > &  CheckedCand)
private

Definition at line 267 of file PndMvdRiemannVertexFinderTask.cxx.

References PndRiemannTrack::addHit(), fHitArray, fHitArray2, fHitBranch, fHitBranch2, fTrackArray, fTrackCandArray, fVertexCut, PndTrackCandHit::GetDetId(), PndTrackCandHit::GetHitId(), PndTrackCand::GetNHits(), PndTrackCand::GetSortedHit(), hit, i, point, PndRiemannTrack::refit(), PndRiemannHit::setDXYZ(), PndRiemannTrack::SetVertexCut(), PndRiemannHit::setXYZ(), PndRiemannTrack::szFit(), and track.

Referenced by Exec().

268 {
269  for(int i=0;i<fTrackCandArray->GetEntriesFast();i++){
271  unsigned int detId,hitId;
273  for(unsigned int j=0;j<Cand->GetNHits();j++){
274  detId=Cand->GetSortedHit(j).GetDetId();
275  hitId=Cand->GetSortedHit(j).GetHitId();
276  PndSdsHit* point = new PndSdsHit();
277  // std::cout<<"detId= "<<detId<<" "<<"hitId= "<<hitId<<std::endl;
278  if ((int)detId==FairRootManager::Instance()->GetBranchId(fHitBranch))
279  point=(PndSdsHit*)fHitArray->At(hitId);
280  else if ((int)detId==FairRootManager::Instance()->GetBranchId(fHitBranch2))
281  point=(PndSdsHit*)fHitArray2->At(hitId);
282  else point=0;
283  if (point!=0){
285  hit.setXYZ(point->GetX(),point->GetY(),point->GetZ());
286  hit.setDXYZ(point->GetDx(),point->GetDy(),point->GetDz());
287  track.addHit(hit);
288  }
289  }
290  track.refit(true);
291  track.szFit(true);
292  track.SetVertexCut(fVertexCut);
293  int size=fTrackArray->GetEntriesFast();
294  new ((*fTrackArray)[size])PndRiemannTrack(track);
295  CheckedCand.push_back(i);
296  }
297  std::cout<<"--------------------------------------------------"<<fTrackCandArray->GetEntriesFast()<<std::endl;
298 }
Int_t i
Definition: run_full.C:25
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
void setDXYZ(double dx, double dy, double dz)
void SetVertexCut(double cut)
void setXYZ(double x, double y, double z)
PndMCTrack * track
Definition: anaLmdCluster.C:89
void refit(bool withErrorCalc=true)
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
void szFit(bool withErrorCalc=true)
void addHit(PndRiemannHit &hit)
PndSdsMCPoint * hit
Definition: anasim.C:70
Int_t GetHitId() const
Int_t GetDetId() const
PndSdsMCPoint * point
Definition: anaLmdCluster.C:72
void PndMvdRiemannVertexFinderTask::Register ( )
private
InitStatus PndMvdRiemannVertexFinderTask::ReInit ( )
virtual

Definition at line 72 of file PndMvdRiemannVertexFinderTask.cxx.

73 {
74 
75  InitStatus stat=kSUCCESS;
76  return stat;
77 
78  /*
79  FairRun* ana = FairRun::Instance();
80  FairRuntimeDb* rtdb=ana->GetRuntimeDb();
81  fGeoPar=(PndMvdGeoPar*)(rtdb->getContainer("PndMvdGeoPar"));
82 
83  return kSUCCESS;
84  */
85 }
void PndMvdRiemannVertexFinderTask::Reset ( )
private
void PndMvdRiemannVertexFinderTask::SetParContainers ( )
virtual

Virtual method Init

Definition at line 62 of file PndMvdRiemannVertexFinderTask.cxx.

63 {
64  // Get Base Container
65 /*
66  FairRun* ana = FairRun::Instance();
67  FairRuntimeDb* rtdb=ana->GetRuntimeDb();
68  fGeoPar = (PndMvdGeoPar*)(rtdb->getContainer("PndMvdGeoPar"));
69 */
70 }
void PndMvdRiemannVertexFinderTask::SetVerbose ( Int_t  verbose)
inline

Definition at line 31 of file PndMvdRiemannVertexFinderTask.h.

References fVerbose, and verbose.

void PndMvdRiemannVertexFinderTask::SetVertexCut ( double  cut)
inline

Definition at line 32 of file PndMvdRiemannVertexFinderTask.h.

References cut, and fVertexCut.

32 { fVertexCut =cut;};
double cut[MAX]
Definition: autocutx.C:36

Member Data Documentation

TH1F* PndMvdRiemannVertexFinderTask::delta

Definition at line 32 of file PndMvdRiemannVertexFinderTask.h.

Referenced by FindVertex(), and Init().

std::pair<double,double> PndMvdRiemannVertexFinderTask::eff

Definition at line 36 of file PndMvdRiemannVertexFinderTask.h.

Referenced by CalcEfficiency().

int PndMvdRiemannVertexFinderTask::fEventNr
private

Definition at line 46 of file PndMvdRiemannVertexFinderTask.h.

Referenced by Exec().

TClonesArray* PndMvdRiemannVertexFinderTask::fHitArray
private

Definition at line 50 of file PndMvdRiemannVertexFinderTask.h.

Referenced by CheckRecoTrack(), Init(), and refit().

TClonesArray* PndMvdRiemannVertexFinderTask::fHitArray2
private

Definition at line 51 of file PndMvdRiemannVertexFinderTask.h.

Referenced by CheckRecoTrack(), Init(), and refit().

TString PndMvdRiemannVertexFinderTask::fHitBranch
private

Definition at line 40 of file PndMvdRiemannVertexFinderTask.h.

Referenced by CheckRecoTrack(), Init(), and refit().

TString PndMvdRiemannVertexFinderTask::fHitBranch2
private

Definition at line 41 of file PndMvdRiemannVertexFinderTask.h.

Referenced by CheckRecoTrack(), Init(), and refit().

TClonesArray* PndMvdRiemannVertexFinderTask::fIdealTrackCandArray
private

Definition at line 55 of file PndMvdRiemannVertexFinderTask.h.

Referenced by FindVertex(), FinishEvent(), FoundCandInMCCands(), and Init().

TString PndMvdRiemannVertexFinderTask::fIdealTrackCandBranch
private

Definition at line 43 of file PndMvdRiemannVertexFinderTask.h.

Referenced by Init().

TClonesArray* PndMvdRiemannVertexFinderTask::fMCTrackArray
private

Definition at line 56 of file PndMvdRiemannVertexFinderTask.h.

Referenced by FindVertex(), and Init().

TString PndMvdRiemannVertexFinderTask::fMCTrackBranch
private

Definition at line 44 of file PndMvdRiemannVertexFinderTask.h.

Referenced by Init().

TClonesArray* PndMvdRiemannVertexFinderTask::fMCVertex
private

Definition at line 59 of file PndMvdRiemannVertexFinderTask.h.

Referenced by FinishEvent(), and Init().

TClonesArray* PndMvdRiemannVertexFinderTask::fTrackArray
private

Definition at line 54 of file PndMvdRiemannVertexFinderTask.h.

Referenced by FindVertex(), FinishEvent(), Init(), and refit().

TString PndMvdRiemannVertexFinderTask::fTrackBranch
private

Definition at line 42 of file PndMvdRiemannVertexFinderTask.h.

Referenced by Init().

TClonesArray* PndMvdRiemannVertexFinderTask::fTrackCandArray
private
int PndMvdRiemannVertexFinderTask::fVerbose
private

Definition at line 47 of file PndMvdRiemannVertexFinderTask.h.

Referenced by SetVerbose().

TClonesArray* PndMvdRiemannVertexFinderTask::fVertex
private

Definition at line 58 of file PndMvdRiemannVertexFinderTask.h.

Referenced by FindVertex(), FinishEvent(), and Init().

double PndMvdRiemannVertexFinderTask::fVertexCut
private

Definition at line 48 of file PndMvdRiemannVertexFinderTask.h.

Referenced by refit(), and SetVertexCut().

std::pair<double,double> PndMvdRiemannVertexFinderTask::ghosts

Definition at line 37 of file PndMvdRiemannVertexFinderTask.h.

Referenced by CalcEfficiency().

TH1F* PndMvdRiemannVertexFinderTask::wrongV

Definition at line 34 of file PndMvdRiemannVertexFinderTask.h.

Referenced by FindVertex(), and Init().


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