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

#include <PndMvdIdealTrackFinderTask.h>

Inheritance diagram for PndMvdIdealTrackFinderTask:

Public Member Functions

 PndMvdIdealTrackFinderTask ()
 
virtual ~PndMvdIdealTrackFinderTask ()
 
 PndMvdIdealTrackFinderTask (const PndMvdIdealTrackFinderTask &)=delete
 
PndMvdIdealTrackFinderTaskoperator= (const PndMvdIdealTrackFinderTask &)=delete
 
virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
virtual void Exec (Option_t *opt)
 
void PrintResult ()
 
void SetVerbose (Int_t verbose)
 

Private Member Functions

void ClearTrackCandMap ()
 
void AddAndExpand (Int_t trackID, Int_t detnum, Int_t iHit, PndSdsHit *theHit)
 
Double_t GetTrackDip (PndMCTrack *myTrack)
 
Double_t GetTrackCurvature (PndMCTrack *myTrack)
 
void Register ()
 
void Reset ()
 
void ProduceHits ()
 
 ClassDef (PndMvdIdealTrackFinderTask, 2)
 

Private Attributes

TString fHitBranchStrip
 
TString fHitBranchPixel
 
TString fClusterBranchStrip
 
TString fClusterBranchPixel
 
TString fDigiBranchStrip
 
TString fDigiBranchPixel
 
TString fMcBranch
 
TString fTrackBranch
 
TClonesArray * fStripHitArray
 
TClonesArray * fPixelHitArray
 
TClonesArray * fStripClusterArray
 
TClonesArray * fPixelClusterArray
 
TClonesArray * fStripDigiArray
 
TClonesArray * fPixelDigiArray
 
TClonesArray * fMcArray
 
TClonesArray * fTrackArray
 
TClonesArray * fTrackCandArray
 
std::map< Int_t, PndTrackCand * > fTrackCandMap
 

Detailed Description

Definition at line 27 of file PndMvdIdealTrackFinderTask.h.

Constructor & Destructor Documentation

PndMvdIdealTrackFinderTask::PndMvdIdealTrackFinderTask ( )

Default constructor

Definition at line 22 of file PndMvdIdealTrackFinderTask.cxx.

22  :
23  FairTask("MVD Ideal Track Finding Task"),
24  fHitBranchStrip("MVDHitsStrip"),
25  fHitBranchPixel("MVDHitsPixel"),
26  fClusterBranchStrip("MVDStripClusterCand"),
27  fClusterBranchPixel("MVDPixelClusterCand"),
28  fDigiBranchStrip("MVDStripDigis"),
29  fDigiBranchPixel("MVDPixelDigis"),
30  fMcBranch("MVDPoint"),
31  fTrackBranch("MCTrack"),
32  fStripHitArray(NULL),
33  fPixelHitArray(NULL),
34  fStripClusterArray(NULL),
35  fPixelClusterArray(NULL),
36  fStripDigiArray(NULL),
37  fPixelDigiArray(NULL),
38  fMcArray(NULL),
39  fTrackArray(NULL),
40  fTrackCandArray(NULL),
42 {
43 }
std::map< Int_t, PndTrackCand * > fTrackCandMap
PndMvdIdealTrackFinderTask::~PndMvdIdealTrackFinderTask ( )
virtual

Destructor

Definition at line 50 of file PndMvdIdealTrackFinderTask.cxx.

51 {
52 }
PndMvdIdealTrackFinderTask::PndMvdIdealTrackFinderTask ( const PndMvdIdealTrackFinderTask )
delete

Member Function Documentation

void PndMvdIdealTrackFinderTask::AddAndExpand ( Int_t  trackID,
Int_t  detnum,
Int_t  iHit,
PndSdsHit theHit 
)
private

Definition at line 198 of file PndMvdIdealTrackFinderTask.cxx.

References fTrackCandMap, PndSdsHit::GetPosition(), and PndTrackCand::setMcTrackId().

Referenced by Exec().

198  {
199  if (fTrackCandMap[trackID] == 0){
200  PndTrackCand *myTCand = new PndTrackCand();
201  //PndMCTrack* myMCTrack = (PndMCTrack*)fTrackArray->At(trackID); //[R.K.03/2017] unused variable
202 // myTCand->setCurv(GetTrackCurvature(myMCTrack));
203 // myTCand->setDip(GetTrackDip(myMCTrack));
204 // myTCand->setInverted(false);
205  //int pdg = myMCTrack->GetPdgCode(); //[R.K.02/2017] Unused variable?
206  //double charge; //[R.K.02/2017] Unused variable?
207  //if(pdg<100000000){ //[R.K.02/2017] Unused variable?
208  //charge = TDatabasePDG::Instance()->GetParticle(pdg)->Charge()/3.; //[R.K.02/2017] Unused variable?
209  //} //[R.K.02/2017] Unused variable?
210  //else{ //[R.K.02/2017] Unused variable?
211  //charge = 0.; //[R.K.02/2017] Unused variable?
212  //} //[R.K.02/2017] Unused variable?
213 // myTCand->setTrackSeed(
214 // myMCTrack->GetStartVertex(),
215 // myMCTrack->GetMomentum(),
216 // charge/myMCTrack->GetMomentum().Mag()
217 // );
218  myTCand->setMcTrackId(trackID);
219  fTrackCandMap[trackID] = myTCand;
220  }
221 
222  fTrackCandMap[trackID]->AddHit(detnum,iHit,theHit->GetPosition().Mag());
223 }
TVector3 GetPosition() const
Definition: PndSdsHit.h:93
void setMcTrackId(int i)
Definition: PndTrackCand.h:72
std::map< Int_t, PndTrackCand * > fTrackCandMap
PndMvdIdealTrackFinderTask::ClassDef ( PndMvdIdealTrackFinderTask  ,
 
)
private
void PndMvdIdealTrackFinderTask::ClearTrackCandMap ( )
private

Definition at line 258 of file PndMvdIdealTrackFinderTask.cxx.

References fTrackCandMap.

Referenced by Exec().

259 {
260  for (std::map<Int_t,PndTrackCand*>::const_iterator kIt=fTrackCandMap.begin();
261  kIt != fTrackCandMap.end(); kIt++){
262  delete(kIt->second);
263  }
264  fTrackCandMap.clear();
265 }
std::map< Int_t, PndTrackCand * > fTrackCandMap
void PndMvdIdealTrackFinderTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 154 of file PndMvdIdealTrackFinderTask.cxx.

References AddAndExpand(), ClearTrackCandMap(), Double_t, fMcArray, fPixelClusterArray, fPixelDigiArray, fPixelHitArray, fStripClusterArray, fStripDigiArray, fStripHitArray, fTrackCandArray, fTrackCandMap, fVerbose, PndSdsHit::GetClusterIndex(), PndSdsDigi::GetDetID(), PndSdsCluster::GetDigiIndex(), PndSdsDigi::GetIndex(), i, and PrintResult().

155 {
156 
157  // Reset output array
158  if ( ! fTrackCandArray )
159  Fatal("Exec", "No trackCandArray");
160  fTrackCandArray->Delete();
161 
162  Int_t nStripHits = fStripHitArray->GetEntriesFast();
163  Int_t nPixelHits = fPixelHitArray->GetEntriesFast();
164 //TODO: Thi becomes easier now, with the FairHit::GetRefIndex()
165  //pixel part
166  for (Int_t iHit = 0; iHit < nPixelHits; iHit++){
167  PndSdsHit* myHit = (PndSdsHit*)(fPixelHitArray->At(iHit));
168  PndSdsCluster* myCluster = (PndSdsCluster*)(fPixelClusterArray->At(myHit->GetClusterIndex()));
169  PndSdsDigiPixel* apixeldigi = (PndSdsDigiPixel*)fPixelDigiArray->At(myCluster->GetDigiIndex(0));
170  if (apixeldigi->GetIndex(0) == -1) continue; // sort out noise
171  PndSdsMCPoint* myPoint = (PndSdsMCPoint*)(fMcArray->At(apixeldigi->GetIndex(0)));
172 
173  AddAndExpand(myPoint->GetTrackID(),apixeldigi->GetDetID(),iHit, myHit);
174 
175  }
176  //strip part
177  for (Int_t iHit = 0; iHit < nStripHits; iHit++){
178  PndSdsHit* myHit = (PndSdsHit*)(fStripHitArray->At(iHit));
179  PndSdsCluster* myCluster = (PndSdsCluster*)(fStripClusterArray->At(myHit->GetClusterIndex()));
180  PndSdsDigiStrip* astripdigi = (PndSdsDigiStrip*)fStripDigiArray->At(myCluster->GetDigiIndex(0));
181  if (astripdigi->GetIndex(0) == -1) continue; // sort out noise
182  PndSdsMCPoint* myPoint = (PndSdsMCPoint*)(fMcArray->At(astripdigi->GetIndex(0)));
183 
184  AddAndExpand(myPoint->GetTrackID(),astripdigi->GetDetID(),iHit, myHit);
185  }
186 
187  if(fVerbose>0) PrintResult();
188 
189  Int_t i = 0;
190  for (std::map<Int_t,PndTrackCand*>::const_iterator kIt=fTrackCandMap.begin(); kIt != fTrackCandMap.end(); kIt++){
191  Double_t eventTime = FairRootManager::Instance()->GetEventTime();
192  kIt->second->SetTimeStamp(eventTime);
193  new((*fTrackCandArray)[i]) PndTrackCand(*(kIt->second));
194  }
196 }
void AddAndExpand(Int_t trackID, Int_t detnum, Int_t iHit, PndSdsHit *theHit)
int fVerbose
Definition: poormantracks.C:24
Int_t i
Definition: run_full.C:25
Class to store the Digis which belong to one cluster This class holds the information which Digi belo...
Definition: PndSdsCluster.h:19
Int_t GetIndex(int i=0) const
Definition: PndSdsDigi.h:63
Class for digitised strip hits.
std::map< Int_t, PndTrackCand * > fTrackCandMap
Int_t GetDigiIndex(Int_t i) const
Definition: PndSdsCluster.h:40
Double_t
Int_t GetDetID() const
Definition: PndSdsDigi.h:61
Data class to store the digi output of a pixel module.
Int_t GetClusterIndex() const
Definition: PndSdsHit.h:94
Double_t PndMvdIdealTrackFinderTask::GetTrackCurvature ( PndMCTrack myTrack)
private

Definition at line 225 of file PndMvdIdealTrackFinderTask.cxx.

References PndMCTrack::GetMomentum(), p, and CAMath::Sqrt().

226 {
227  TVector3 p = myTrack->GetMomentum();
228  return (2/TMath::Sqrt(p.Px()*p.Px() + p.Py()*p.Py()));
229 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
Double_t p
Definition: anasim.C:58
Double_t PndMvdIdealTrackFinderTask::GetTrackDip ( PndMCTrack myTrack)
private

Definition at line 231 of file PndMvdIdealTrackFinderTask.cxx.

References PndMCTrack::GetMomentum(), p, and CAMath::Sqrt().

232 {
233  TVector3 p= myTrack->GetMomentum();
234  return (p.Mag()/TMath::Sqrt(p.Px()*p.Px() + p.Py()*p.Py()));
235 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
Double_t p
Definition: anasim.C:58
InitStatus PndMvdIdealTrackFinderTask::Init ( )
virtual

Definition at line 82 of file PndMvdIdealTrackFinderTask.cxx.

References fClusterBranchPixel, fClusterBranchStrip, fDigiBranchPixel, fDigiBranchStrip, fHitBranchPixel, fHitBranchStrip, fMcArray, fMcBranch, fPixelClusterArray, fPixelDigiArray, fPixelHitArray, fStripClusterArray, fStripDigiArray, fStripHitArray, fTrackArray, fTrackBranch, and fTrackCandArray.

83 {
84 
85  FairRootManager* ioman = FairRootManager::Instance();
86 
87  if ( ! ioman )
88  {
89  std::cout << "-E- PndMvdIdealTrackFinderTask::Init: "
90  << "RootManager not instantiated!" << std::endl;
91  return kFATAL;
92  }
93 
94  // Get input array
95  fStripHitArray = (TClonesArray*) ioman->GetObject(fHitBranchStrip);
96  if ( !fStripHitArray){
97  std::cout << "-W- PndMvdIdealTrackFinderTask::Init: " << "No fStripHitArray!" << std::endl;
98  return kERROR;
99  }
100 
101  fPixelHitArray = (TClonesArray*) ioman->GetObject(fHitBranchPixel);
102  if ( !fPixelHitArray){
103  std::cout << "-W- PndMvdIdealTrackFinderTask::Init: " << "No fPixelHitArray!" << std::endl;
104  return kERROR;
105  }
106 
107  fStripClusterArray = (TClonesArray*) ioman->GetObject(fClusterBranchStrip);
108  if ( !fStripClusterArray){
109  std::cout << "-W- PndMvdIdealTrackFinderTask::Init: " << "No StripclusterArray!" << std::endl;
110  return kERROR;
111  }
112  fPixelClusterArray = (TClonesArray*) ioman->GetObject(fClusterBranchPixel);
113  if ( !fPixelClusterArray){
114  std::cout << "-W- PndMvdIdealTrackFinderTask::Init: " << "No PixelclusterArray!" << std::endl;
115  return kERROR;
116  }
117 
118 
119  fStripDigiArray = (TClonesArray*) ioman->GetObject(fDigiBranchStrip);
120  if ( !fStripDigiArray){
121  std::cout << "-W- PndMvdIdealTrackFinderTask::Init: " << "No StripdigiArray!" << std::endl;
122  return kERROR;
123  }
124 
125  fPixelDigiArray = (TClonesArray*) ioman->GetObject(fDigiBranchPixel);
126  if ( !fPixelDigiArray){
127  std::cout << "-W- PndMvdIdealTrackFinderTask::Init: " << "No PixeldigiArray!" << std::endl;
128  return kERROR;
129  }
130 
131  fMcArray = (TClonesArray*) ioman->GetObject(fMcBranch);
132  if ( !fMcArray){
133  std::cout << "-W- PndMvdIdealTrackFinderTask::Init: " << "No mcArray!" << std::endl;
134  return kERROR;
135  }
136 
137  fTrackArray = (TClonesArray*) ioman->GetObject(fTrackBranch);
138  if ( !fTrackArray){
139  std::cout << "-W- PndMvdIdealTrackFinderTask::Init: " << "No trackArray!" << std::endl;
140  return kERROR;
141  }
142 
143  fTrackCandArray = new TClonesArray("PndTrackCand");
144  ioman->Register("MVDIdealTrackCand", "MVD", fTrackCandArray, kTRUE);
145 
146  std::cout << "-I- PndMvdIdealTrackFinderTask: Initialisation successfull" << std::endl;
147  return kSUCCESS;
148 }
PndMvdIdealTrackFinderTask& PndMvdIdealTrackFinderTask::operator= ( const PndMvdIdealTrackFinderTask )
delete
void PndMvdIdealTrackFinderTask::PrintResult ( )

Definition at line 237 of file PndMvdIdealTrackFinderTask.cxx.

References fPixelHitArray, fStripHitArray, fTrackCandMap, PndTrackCandHit::GetDetId(), PndTrackCandHit::GetHitId(), PndTrackCand::GetNHits(), PndTrackCand::GetSortedHit(), and i.

Referenced by Exec().

238 {
239  std::cout << "**** TrackFinding *****" << std::endl;
240  Int_t nStripHits = fStripHitArray->GetEntriesFast();
241  for (std::map<Int_t, PndTrackCand*>::const_iterator kIt=fTrackCandMap.begin();
242  kIt != fTrackCandMap.end(); kIt++){
243  std::cout << "TrackID: " << kIt->first << std::endl;
244  PndTrackCand* trackCand = kIt->second;
245  unsigned int detId, hitId;
246  for (unsigned int i = 0; i < trackCand->GetNHits(); i++){
247  detId=trackCand->GetSortedHit(i).GetDetId();
248  hitId=trackCand->GetSortedHit(i).GetHitId();
249  PndSdsHit* myHit;
250  if(hitId<(UInt_t)nStripHits) myHit = (PndSdsHit*)(fStripHitArray->At(hitId));
251  else myHit = (PndSdsHit*)(fPixelHitArray->At(hitId - nStripHits));
252  std::cout << "Detector no. " << detId <<": "<< *myHit;
253  }
254  }
255  std::cout << std::endl;
256 }
Int_t i
Definition: run_full.C:25
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
std::map< Int_t, PndTrackCand * > fTrackCandMap
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
Int_t GetHitId() const
Int_t GetDetId() const
void PndMvdIdealTrackFinderTask::ProduceHits ( )
private
void PndMvdIdealTrackFinderTask::Register ( )
private
InitStatus PndMvdIdealTrackFinderTask::ReInit ( )
virtual

Definition at line 66 of file PndMvdIdealTrackFinderTask.cxx.

67 {
68 
69  InitStatus stat=kSUCCESS;
70  return stat;
71 
72  /*
73  FairRun* ana = FairRun::Instance();
74  FairRuntimeDb* rtdb=ana->GetRuntimeDb();
75  fGeoPar=(PndMvdGeoPar*)(rtdb->getContainer("PndMvdGeoPar"));
76 
77  return kSUCCESS;
78  */
79 }
void PndMvdIdealTrackFinderTask::Reset ( )
private
void PndMvdIdealTrackFinderTask::SetParContainers ( )
virtual

Virtual method Init

Definition at line 56 of file PndMvdIdealTrackFinderTask.cxx.

57 {
58  // Get Base Container
59 /*
60  FairRun* ana = FairRun::Instance();
61  FairRuntimeDb* rtdb=ana->GetRuntimeDb();
62  fGeoPar = (PndMvdGeoPar*)(rtdb->getContainer("PndMvdGeoPar"));
63 */
64 }
void PndMvdIdealTrackFinderTask::SetVerbose ( Int_t  verbose)
inline

Definition at line 47 of file PndMvdIdealTrackFinderTask.h.

References fVerbose, and verbose.

Referenced by locT_all(), locT_theta(), and reco_mvd().

47 { fVerbose = verbose;};
int fVerbose
Definition: poormantracks.C:24
#define verbose

Member Data Documentation

TString PndMvdIdealTrackFinderTask::fClusterBranchPixel
private

Definition at line 61 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TString PndMvdIdealTrackFinderTask::fClusterBranchStrip
private

Definition at line 60 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TString PndMvdIdealTrackFinderTask::fDigiBranchPixel
private

Definition at line 63 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TString PndMvdIdealTrackFinderTask::fDigiBranchStrip
private

Definition at line 62 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TString PndMvdIdealTrackFinderTask::fHitBranchPixel
private

Definition at line 59 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TString PndMvdIdealTrackFinderTask::fHitBranchStrip
private

Definition at line 58 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TClonesArray* PndMvdIdealTrackFinderTask::fMcArray
private

Definition at line 74 of file PndMvdIdealTrackFinderTask.h.

Referenced by Exec(), and Init().

TString PndMvdIdealTrackFinderTask::fMcBranch
private

Definition at line 64 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TClonesArray* PndMvdIdealTrackFinderTask::fPixelClusterArray
private

Definition at line 71 of file PndMvdIdealTrackFinderTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndMvdIdealTrackFinderTask::fPixelDigiArray
private

Definition at line 73 of file PndMvdIdealTrackFinderTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndMvdIdealTrackFinderTask::fPixelHitArray
private

Definition at line 69 of file PndMvdIdealTrackFinderTask.h.

Referenced by Exec(), Init(), and PrintResult().

TClonesArray* PndMvdIdealTrackFinderTask::fStripClusterArray
private

Definition at line 70 of file PndMvdIdealTrackFinderTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndMvdIdealTrackFinderTask::fStripDigiArray
private

Definition at line 72 of file PndMvdIdealTrackFinderTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndMvdIdealTrackFinderTask::fStripHitArray
private

Input array of PndSdsDigis

Definition at line 68 of file PndMvdIdealTrackFinderTask.h.

Referenced by Exec(), Init(), and PrintResult().

TClonesArray* PndMvdIdealTrackFinderTask::fTrackArray
private

Definition at line 75 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TString PndMvdIdealTrackFinderTask::fTrackBranch
private

Definition at line 65 of file PndMvdIdealTrackFinderTask.h.

Referenced by Init().

TClonesArray* PndMvdIdealTrackFinderTask::fTrackCandArray
private

Output array of PndMvdHits

Definition at line 78 of file PndMvdIdealTrackFinderTask.h.

Referenced by Exec(), and Init().

std::map<Int_t, PndTrackCand*> PndMvdIdealTrackFinderTask::fTrackCandMap
private

Definition at line 79 of file PndMvdIdealTrackFinderTask.h.

Referenced by AddAndExpand(), ClearTrackCandMap(), Exec(), and PrintResult().


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