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

#include <PndHypIdealPRTask.h>

Inheritance diagram for PndHypIdealPRTask:

Public Member Functions

 PndHypIdealPRTask ()
 
 ~PndHypIdealPRTask ()
 
void AddHitBranch (unsigned int detId, const TString &m)
 
void SetPersistence (Bool_t opt=kTRUE)
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
virtual void Exec (Option_t *opt)
 
void WriteHistograms ()
 

Private Member Functions

Int_t GetChargeIon (Int_t ion)
 
void AddHitToTrack (Int_t trackID, Int_t detnum, Int_t iHit)
 
void ClearTrackCandMap ()
 

Private Attributes

TClonesArray * fHitArray
 
TClonesArray * fPointArray
 
TClonesArray * fSdsArray
 
TClonesArray * fTrackArray
 
TClonesArray * fMcArray
 
Bool_t fPersistence
 
std::map< unsigned int, TStringfHitBranchNameMap
 
std::map< unsigned int,
TClonesArray * > 
fHitBranchMap
 
std::map< Int_t, GFTrackCand * > fTrackCandMap
 
int fEventNr
 
TH1D * fPH
 
GFRecoHitFactoryfTheRecoHitFactory
 

Detailed Description

Definition at line 47 of file PndHypIdealPRTask.h.

Constructor & Destructor Documentation

PndHypIdealPRTask::PndHypIdealPRTask ( )

Definition at line 49 of file PndHypIdealPRTask.cxx.

50  : FairTask("Ideal Pattern Reco"), fPersistence(kFALSE),fEventNr(0)
51 {
52 }
PndHypIdealPRTask::~PndHypIdealPRTask ( )

Definition at line 55 of file PndHypIdealPRTask.cxx.

References fPH.

56 {
57  if(fPH!=NULL)delete fPH;
58 }

Member Function Documentation

void PndHypIdealPRTask::AddHitBranch ( unsigned int  detId,
const TString m 
)
inline

Definition at line 59 of file PndHypIdealPRTask.h.

References fHitBranchNameMap, and m.

59 {fHitBranchNameMap[detId]=m;};
std::map< unsigned int, TString > fHitBranchNameMap
__m128 m
Definition: P4_F32vec4.h:28
void PndHypIdealPRTask::AddHitToTrack ( Int_t  trackID,
Int_t  detnum,
Int_t  iHit 
)
private

Definition at line 228 of file PndHypIdealPRTask.cxx.

References fMcArray, fTrackCandMap, GetChargeIon(), PndMCTrack::GetMomentum(), PndMCTrack::GetPdgCode(), PndMCTrack::GetStartVertex(), GFTrackCand::setMcTrackId(), GFTrackCand::setPdgCode(), and GFTrackCand::setTrackSeed().

Referenced by Exec().

229  {
230 
231  if(fTrackCandMap[trackID]==NULL){
232  // create new track
233  //std::cout<<"Creating new track candidate with id="<<id<<std::endl;
234 
235  GFTrackCand *myTCand=new GFTrackCand();
236 
237  // Get MCTrack
238  PndMCTrack* mc=(PndMCTrack*)fMcArray->At(trackID);
239 
240 
241  int pdg=mc->GetPdgCode();
242  double q;
243  int PDG;
244 
245  std::cout<<" particle "<<pdg<<std::endl;
246 
247 
248  if(pdg>5000){
249  q = GetChargeIon(pdg);
250  //std::cout<<" particle "<<q<<std::endl;
251  }
252  else { q=TDatabasePDG::Instance()->GetParticle(pdg)->Charge()/3.;
253  //std::cout<<" particle "<<q<<std::endl;
254  }
255 
256  myTCand->setTrackSeed(
257  mc->GetStartVertex(),
258  mc->GetMomentum(),
259  q/mc->GetMomentum().Mag()
260  );
261 
262 
263  myTCand->setMcTrackId(trackID);
264  myTCand->setPdgCode(pdg);
265  fTrackCandMap[trackID]= myTCand;
266 
267  }
268  // add hit to track
269  // set detectorid=2;
270  fTrackCandMap[trackID]->addHit(detnum,iHit);
271 
272 
273 }
Int_t GetChargeIon(Int_t ion)
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
void setMcTrackId(int i)
set the MCT track id, for MC simulations
Definition: GFTrackCand.h:149
void setTrackSeed(const TVector3 &pos, const TVector3 &direction, const double qop)
set the seed values for track: pos, direction, q/p
Definition: GFTrackCand.h:155
Track candidate – a list of cluster indices.
Definition: GFTrackCand.h:55
TClonesArray * fMcArray
void setPdgCode(int pdgCode)
set a particle hypothesis in form of a PDG code
Definition: GFTrackCand.h:163
std::map< Int_t, GFTrackCand * > fTrackCandMap
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
void PndHypIdealPRTask::ClearTrackCandMap ( )
private

Definition at line 304 of file PndHypIdealPRTask.cxx.

References fTrackCandMap.

Referenced by Exec().

305 {
306  for (std::map<Int_t,GFTrackCand*>::const_iterator ci=fTrackCandMap.begin();
307  ci != fTrackCandMap.end(); ci++){
308  delete(ci->second);
309  }
310  fTrackCandMap.clear();
311 }
std::map< Int_t, GFTrackCand * > fTrackCandMap
void PndHypIdealPRTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 131 of file PndHypIdealPRTask.cxx.

References AddHitToTrack(), ClearTrackCandMap(), fEventNr, fHitArray, fHitBranchMap, fMcArray, fPH, fPointArray, fSdsArray, fTrackArray, fTrackCandMap, GFTrackCand::getHit(), GFTrackCand::getNHits(), hit, and point.

132 {
133  std::cout<<"PndHypIdealPRTask::Exec"<<std::endl;
134  // Reset output Array
135  if(fTrackArray==0) Fatal("PndHypIdealPR::Exec)","No TrackArray");
136  fTrackArray->Delete();
137 
138  // use McId to distinguish data from different tracks
139  //std::map<unsigned int,GFTrackCand*> candmap;
140  std::cout<<"<<<<< Event "<<fEventNr++<<" <<<"<<std::endl;
141 
142 
143  std::map<unsigned int,TClonesArray*>::iterator iter=fHitBranchMap.begin();
144  while(iter!=fHitBranchMap.end()){
145 
146  fHitArray = iter->second;
147  //loop over points
148  Int_t np=fHitArray->GetEntriesFast();
149  FairMCPoint* point;
150  FairHit* hit;
151 
152  for(Int_t ip=0; ip<np; ++ip){
153 
154  if(iter->first==2){
155  hit = (PndHypHit*)fHitArray->At(ip);
156 
157  point = (PndHypPoint*)fPointArray->At(hit->GetRefIndex());
158  }
159  else if(iter->first==3){
160  hit = (PndSdsHit*)fHitArray->At(ip);
161  point = (PndSdsMCPoint*)fSdsArray->At(hit->GetRefIndex());
162  }
163  if(point==0)continue;
164 
165  unsigned int id=point->GetTrackID();
166  // cut on insane ids
167  if(id>10000)continue;
168  // look for track in mc map
169  AddHitToTrack(id,iter->first,ip);
170 
171  /* if(candmap[id]==NULL){
172  // create new track
173  //std::cout<<"Creating new track candidate with id="<<id<<std::endl;
174  candmap[id]=new GFTrackCand;
175  }
176  // add hit to track
177  // set detectorid=2;
178  candmap[id]->addHit(iter->first,ip);
179  */
180  }// end loop over cluster
181 
182  ++iter;
183  }
184 
185  // ------------------------------------------------------------------------
186  // try to find some starting values
187  // loop over tracks
188 
189  std::map<Int_t ,GFTrackCand*>::iterator candIter=fTrackCandMap.begin();
190  while(candIter!=fTrackCandMap.end()){
191  GFTrackCand* cand=candIter->second;
192  if(cand->getNHits()<3){
193  ++candIter;
194  continue;
195  }
196  TVector3 frmom;
197  unsigned int detId, hitId;
198  cand->getHit(0, detId, hitId);
199  PndHypHit* myHit;
200  myHit = (PndHypHit*)fHitArray->At(hitId);
201  FairMCPoint* pointF;
202 
203  if(detId==2) pointF =(FairMCPoint*)fPointArray->At(myHit->GetRefIndex());
204  std::cout << "Detector no. " << detId <<": "<< *pointF;
205  pointF->Momentum(frmom);
206  fPH->Fill(frmom.Mag());
207 
208 
209  // create track object
210 
211  new((*fTrackArray)[fTrackArray->GetEntriesFast()]) GFTrackCand(*(candIter->second));
212 
213 
214  ++candIter;
215  }// end loop over tracks
216 
217  std::cout<<fTrackArray->GetEntriesFast()<<" tracks created "<<std::endl;
218 
220  fHitArray->Delete();
221  fMcArray->Delete();
222 
223 
224 
225  //return;
226 }
TClonesArray * fPointArray
unsigned int getNHits() const
Definition: GFTrackCand.h:113
TClonesArray * fSdsArray
std::map< unsigned int, TClonesArray * > fHitBranchMap
void getHit(unsigned int i, unsigned int &detId, unsigned int &hitId) const
Get detector ID and cluster index (hitId) for hit number i.
Definition: GFTrackCand.h:84
void AddHitToTrack(Int_t trackID, Int_t detnum, Int_t iHit)
Track candidate – a list of cluster indices.
Definition: GFTrackCand.h:55
TClonesArray * fHitArray
TClonesArray * fMcArray
TClonesArray * fTrackArray
std::map< Int_t, GFTrackCand * > fTrackCandMap
PndSdsMCPoint * hit
Definition: anasim.C:70
PndSdsMCPoint * point
Definition: anaLmdCluster.C:72
Int_t PndHypIdealPRTask::GetChargeIon ( Int_t  ion)
private

Definition at line 275 of file PndHypIdealPRTask.cxx.

References L, and Z.

Referenced by AddHitToTrack().

276 {
277  Int_t A,Z,L;
278 
279  if(ion>1000000000&&(ion<1010000000))
280  { ion -= 1000000000;
281  Z = ion/10000;
282  ion -= 10000*Z;
283  A = ion/10;
284  cout<<" ion charge "<<Z<<endl;
285 
286  return Z;
287 
288  }
289 if((ion>1010000000||ion>1020000000))
290  {
291  ion -= 1000000000;
292  L = ion/10000000;
293  ion -= 10000000*L;
294  Z = ion/10000;
295  ion -= 10000*Z;
296  A = ion/10;
297  cout<<L<<" hypernuclei charge "<<Z<<endl;
298 
299  return Z;
300 
301  }
302 }
double Z
Definition: anaLmdDigi.C:68
InitStatus PndHypIdealPRTask::Init ( )
virtual

Virtual method Init

Definition at line 69 of file PndHypIdealPRTask.cxx.

References fHitBranchMap, fHitBranchNameMap, fMcArray, fPersistence, fPH, fPointArray, fSdsArray, and fTrackArray.

70 {
71  //Get ROOT Manager
72  FairRootManager* ioman= FairRootManager::Instance();
73 
74  if(ioman==0)
75  {
76  Error("PndHypIdealPRTask::Init","RootManager not instantiated!");
77  return kERROR;
78  }
79 
80  // open hit arrays
81  std::map<unsigned int,TString>::iterator iter=fHitBranchNameMap.begin();
82  while(iter!=fHitBranchNameMap.end()){
83  TClonesArray* ar=(TClonesArray*) ioman->GetObject(iter->second);
84  if(ar==0){
85  Error("PndHypDPRTask::Init","point-array %s not found!",iter->second.Data());
86  }
87  else{
88  fHitBranchMap[iter->first] = ar;
89  }
90  ++iter;
91  }//end loops over hit types
92 
93 
94  // open MCTruth array
95  fMcArray=(TClonesArray*) ioman->GetObject("MCTrack");
96  if(fMcArray==0){
97  Error("PndHypDPRTask::Init","mctrack-array not found!");
98  return kERROR;
99  }
100 
101  fPointArray=(TClonesArray*) ioman->GetObject("HypPoint");
102  if(fPointArray==0){
103  Error("PndHypDPRTask::Init","hyp hit-array not found!");
104  return kERROR;
105  }
106 
107  fSdsArray=(TClonesArray*) ioman->GetObject("MVDPoint");
108  if(fSdsArray==0){
109  Error("PndHypDPRTask::Init","mvd hit-array not found!");
110  return kERROR;
111  }
112 
113 
114  // create and register output array
115  fTrackArray = new TClonesArray("GFTrackCand");
116  ioman->Register("HypTrackCand","HYP",fTrackArray,fPersistence);
117 
118 
119  // setup histograms
120  fPH=new TH1D("pH","p",1000,0.02,0.3);
121 
122 
123  std::cout << "-I- PndHypIdealPRTask: Initialisation successfull" << std::endl;
124 
125 
126  return kSUCCESS;
127 }
TClonesArray * fPointArray
std::map< unsigned int, TString > fHitBranchNameMap
TClonesArray * fSdsArray
std::map< unsigned int, TClonesArray * > fHitBranchMap
TClonesArray * fMcArray
TClonesArray * fTrackArray
InitStatus PndHypIdealPRTask::ReInit ( )
virtual

Definition at line 60 of file PndHypIdealPRTask.cxx.

61 {
62 
63  InitStatus stat=kERROR;
64  return stat;
65 
66 }
void PndHypIdealPRTask::SetPersistence ( Bool_t  opt = kTRUE)
inline

Definition at line 60 of file PndHypIdealPRTask.h.

References fPersistence.

60 {fPersistence=opt;}
void PndHypIdealPRTask::WriteHistograms ( )

Definition at line 314 of file PndHypIdealPRTask.cxx.

References file, and fPH.

314  {
315  TFile* file = FairRootManager::Instance()->GetOutFile();
316  file->cd();
317  file->mkdir("HypPrefit");
318  file->cd("HypPrefit");
319 
320 
321  fPH->Write();
322  delete fPH;
323  fPH=NULL;
324 
325  file->Close();
326  delete file;
327 
328 }
TFile * file

Member Data Documentation

int PndHypIdealPRTask::fEventNr
private

Definition at line 91 of file PndHypIdealPRTask.h.

Referenced by Exec().

TClonesArray* PndHypIdealPRTask::fHitArray
private

Definition at line 79 of file PndHypIdealPRTask.h.

Referenced by Exec().

std::map<unsigned int,TClonesArray*> PndHypIdealPRTask::fHitBranchMap
private

Definition at line 89 of file PndHypIdealPRTask.h.

Referenced by Exec(), and Init().

std::map<unsigned int,TString> PndHypIdealPRTask::fHitBranchNameMap
private

Definition at line 88 of file PndHypIdealPRTask.h.

Referenced by AddHitBranch(), and Init().

TClonesArray* PndHypIdealPRTask::fMcArray
private

Definition at line 83 of file PndHypIdealPRTask.h.

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

Bool_t PndHypIdealPRTask::fPersistence
private

Definition at line 85 of file PndHypIdealPRTask.h.

Referenced by Init(), and SetPersistence().

TH1D* PndHypIdealPRTask::fPH
private

Definition at line 92 of file PndHypIdealPRTask.h.

Referenced by Exec(), Init(), WriteHistograms(), and ~PndHypIdealPRTask().

TClonesArray* PndHypIdealPRTask::fPointArray
private

Definition at line 80 of file PndHypIdealPRTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndHypIdealPRTask::fSdsArray
private

Definition at line 81 of file PndHypIdealPRTask.h.

Referenced by Exec(), and Init().

GFRecoHitFactory* PndHypIdealPRTask::fTheRecoHitFactory
private

Definition at line 97 of file PndHypIdealPRTask.h.

TClonesArray* PndHypIdealPRTask::fTrackArray
private

Definition at line 82 of file PndHypIdealPRTask.h.

Referenced by Exec(), and Init().

std::map<Int_t, GFTrackCand*> PndHypIdealPRTask::fTrackCandMap
private

Definition at line 90 of file PndHypIdealPRTask.h.

Referenced by AddHitToTrack(), ClearTrackCandMap(), and Exec().


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