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

#include <PndHypDPatternRecoTask.h>

Inheritance diagram for PndHypDPatternRecoTask:

Public Member Functions

 PndHypDPatternRecoTask ()
 
 ~PndHypDPatternRecoTask ()
 
void SetVtxAbsName (TString name)
 
void AddHitBranch (unsigned int detId, const TString &m)
 
void SetPersistence (Bool_t opt=kTRUE)
 
void SetField (FairField *f)
 
void SetHitFL (Bool_t opt)
 
void UseGeane (bool f=true)
 
void UseMVD (bool mvd)
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 

Private Member Functions

Int_t GetChargeIon (Int_t ion)
 

Private Attributes

TClonesArray * fPointArray
 
TClonesArray * fSdsArray
 
TClonesArray * fTrackArray
 
TClonesArray * fHitArray
 
TClonesArray * fMcArray
 
std::map< unsigned int, TStringfHitBranchNameMap
 
std::map< unsigned int,
TClonesArray * > 
fHitBranchMap
 
int fEventNr
 
TString fVtxName
 
Bool_t fPersistence
 
Bool_t fUseGeane
 
bool fUseMVD
 
Bool_t fMCvalue
 
Bool_t fVtx
 
GFRecoHitFactoryfTheRecoHitFactory
 
FairField * fField
 
FairGeanePro * fGeanePro
 

Detailed Description

Definition at line 35 of file PndHypDPatternRecoTask.h.

Constructor & Destructor Documentation

PndHypDPatternRecoTask::PndHypDPatternRecoTask ( )

Definition at line 54 of file PndHypDPatternRecoTask.cxx.

PndHypDPatternRecoTask::~PndHypDPatternRecoTask ( )

Definition at line 61 of file PndHypDPatternRecoTask.cxx.

62 {
63 }

Member Function Documentation

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

Definition at line 53 of file PndHypDPatternRecoTask.h.

References fHitBranchNameMap, and m.

53 {fHitBranchNameMap[detId]=m;};
__m128 m
Definition: P4_F32vec4.h:28
std::map< unsigned int, TString > fHitBranchNameMap
void PndHypDPatternRecoTask::Exec ( Option_t *  opt)
virtual

Definition at line 131 of file PndHypDPatternRecoTask.cxx.

References GFTrackCand::addHit(), b, fb, fEventNr, fField, fGeanePro, fHitArray, fHitBranchMap, fMcArray, fMCvalue, fPointArray, fSdsArray, fTrackArray, fUseGeane, fUseMVD, fVtx, fVtxName, GetChargeIon(), GFTrackCand::getHit(), PndMCTrack::GetMomentum(), PndMCTrack::GetMotherID(), GFTrackCand::getNHits(), PndMCTrack::GetPdgCode(), PndMCTrack::GetStartVertex(), gGeoManager, hit, mom, point, pos, GFTrack::setCandidate(), trk, TString, and v.

132 {
133  std::cout<<"PndHypDPatternRecoTask::Exec"<<std::endl;
134  // Reset output Array
135  if(fTrackArray==0) Fatal("PndHypDPatternReco::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  FairMCPoint* point;
142  FairHit* hit;
143 
144  std::map<unsigned int,TClonesArray*>::iterator iter=fHitBranchMap.begin();
145  while(iter!=fHitBranchMap.end()){
146  //TClonesArray* pointArray = iter->second;//_hitBranchMap[iter->first];
147  fHitArray = iter->second;
148  //loop over points
149  Int_t np=fHitArray->GetEntriesFast();
150 
151  for(Int_t ip=0; ip<np; ++ip){
152  //PndHypPoint* point=(PndHypPoint*)hitArray->At(ip);
153 
154  if(iter->first==2){
155  hit =(PndHypHit*)fHitArray->At(ip);
156  point=(PndHypPoint*)fPointArray->At(hit->GetRefIndex());
157  }
158 
159  if(fUseMVD==true){
160  if(iter->first==3){
161  hit =(PndSdsHit*)fHitArray->At(ip);
162  point=(PndSdsMCPoint*)fSdsArray->At(hit->GetRefIndex());
163  }
164  }
165 
166  if(point==0)continue;
167 
168  unsigned int id=point->GetTrackID();
169  // cut on insane ids
170  if(id>10000)continue;
171  // look for track in mc map
172  if(candmap[id]==NULL){
173  // create new track
174  //std::cout<<"Creating new track candidate with id="<<id<<std::endl;
175  candmap[id]=new GFTrackCand;
176  }
177  // add hit to track
178  // set detectorid=2;
179  candmap[id]->addHit(iter->first,ip);
180  }// end loop over cluster
181  ++iter;
182  }
183 
184  // ------------------------------------------------------------------------
185  // try to find some starting values
186  // loop over tracks
187 
188  std::map<unsigned int,GFTrackCand*>::iterator candIter=candmap.begin();
189  while(candIter!=candmap.end()){
190  GFTrackCand* cand=candIter->second;
191  if(cand->getNHits()<3){
192  ++candIter;
193  continue;
194  }
195 
196 
197 
198  // Get MCTrack
199  PndMCTrack* mc=(PndMCTrack*)fMcArray->At(candIter->first);
200  if(mc==0){
201  Error("PndHypDPRTask::Exec","MCTrack Id=&i not found!",candIter->first);
202  ++candIter;
203  continue;
204  }
205 
206  int pdg=mc->GetPdgCode();
207  double q;int PDG;
208 
209  std::cout<<" particle "<<pdg<<std::endl;
210 
211 
212  if(pdg>5000){
213 
214  q = GetChargeIon(pdg);
215  //std::cout<<" particle "<<q<<std::endl;
216  }
217  else{ q=TDatabasePDG::Instance()->GetParticle(pdg)->Charge()/3.;
218  //std::cout<<" particle "<<q<<std::endl;
219  }
220 
221  TVector3 vtx;
222  TString VolName;
223 
224  if(fVtx){
225  if(mc->GetMotherID()!=-1){
226  ++candIter;
227  continue;
228  }
229 
230  vtx = mc->GetStartVertex();
231  std::cout<<" pion prim "<<mc->GetMotherID()<<std::endl;
232  std::cout<<" pion vertex "<<vtx.X()<<" "<<vtx.Y()<<" "<<vtx.Z()<<std::endl;
233 
234  VolName = gGeoManager->FindNode(vtx.X(),vtx.Y(),vtx.Z())->GetName();
235 
236  std::cout<<" vol name "<<VolName.Data()<<std::endl;
237 
238  if(!(VolName.Contains(fVtxName.Data())))
239  {
240  ++candIter;
241  continue;
242  }
243 
244  std::cout<<" vol name "<<VolName.Data()<<" "<<fVtxName.Data()<<std::endl;
245  }
246 
247  TVector3 pos;
248  double pre;
249  TVector3 mom;
250 
251  if(fMCvalue){
252  unsigned int detId, hitId;
253  //getting the first or last hit in the track candidate
254  cand->getHit(0, detId, hitId);
255  PndHypHit* myHit;
256  myHit = (PndHypHit*)fHitArray->At(hitId);
257  FairMCPoint* pointF=(FairMCPoint*)fPointArray->At(myHit->GetRefIndex());
258 
259  if(pointF==0){
260  Error("PndHypDPRTask::Exec","MCTrack Id=&i not found!",candIter->first);
261  ++candIter;
262  continue;
263  }
264  std::cout << "Detector no. " << detId <<std::endl;
265 
266 
267  pointF->Position(pos);
268  pointF->Momentum(mom);
269 
270  }else {
271 
272  pos=mc->GetStartVertex()+TVector3(0.001,0.001,0.001);
273  pre=mc->GetMomentum().Mag()*0.02;
274  mom=mc->GetMomentum()+TVector3(gRandom->Gaus(0,pre),
275  gRandom->Gaus(0,pre),
276  gRandom->Gaus(0,pre));
277 
278  std::cout << "MC starting values " << mom.Mag()<<std::endl;
279  }
280 
281  if(mom.Mag()==0) {
282  Error("PndHypDPRTask::Exec","Momentum is zero!",candIter->first);
283  ++candIter;
284  continue;
285  }
286 
287  TVector3 poserr(0.01,0.01,0.005);
288  TVector3 momerr=0.3*mom; // 10% error
289  momerr+=TVector3(0.1,0.1,0.1);
290 
291  TVector3 u(1.,0.,0.);//=mom.Orthogonal();
292  //u.SetMag(1.);
293  TVector3 v(0.,1.,0.);//=mom.Cross(u);
294  //v.SetMag(1.);
295 
296  // create track-representation object and initialize with start values
297  GFAbsTrackRep* rep=0;
298  // setting the magnetic field properly
299  fField= FairRunAna::Instance()->GetField();
301  GFFieldManager * fb;
302  //fb->getInstance();
303  //fb->init(b);
304 
305  if(fUseGeane){
306  GFDetPlane pl(pos,u,v);
307  if(q>0)PDG=211;
308  if(q<0)PDG=-211;
309 
310  GeaneTrackRep* grep=new GeaneTrackRep(fGeanePro,pl,mom,poserr,momerr,q,PDG);
311 
312 
313  //grep->setPropDir(1); // propagate in flight direction!
314  rep=grep;
315  }
316  else { // use LSLTrackRep
317  // calc momentum projections
318  TVector3 dir=mom.Unit();
319  double dxdz=dir.X()/dir.Z();
320  double dydz=dir.Y()/dir.Z();
321  double qp=q/mom.Mag();
322  rep=new LSLTrackRep(pos.Z(),pos.X(),pos.Y(),dxdz,dydz,qp,
323  poserr.X(),poserr.Y(),0.1,0.1,0.1,b);//NULL instead of b field
324 
325  // use RKTrackRep
326  //rep=new RKtrackRep(pos,mom,poserr,momerr,pdg);
327 
328  }
329 
330 
331 
332 
333  // create track object
334  GFTrack* trk=new((*fTrackArray)[fTrackArray->GetEntriesFast()]) GFTrack(rep);
335  //std::cout<<" particle "<<cand->getNHits()<<std::endl;
336  if(cand==NULL) {
337  Error("PndHypDPRTask::Exec","Momentum is zero!",candIter->first);
338  //std::cout<<_trackArray->GetEntriesFast()<<" tracks created"<<std::endl;
339  ++candIter;
340  continue;
341  }
342  if(trk!=NULL)trk->setCandidate(*cand); // here the candidate is copied!
343  //std::cout<<" particle "<<(trk->getNumHits())<<std::endl;
344 
345  ++candIter;
346  }// end loop over tracks
347 
348  std::cout<<fTrackArray->GetEntriesFast()<<" tracks created"<<std::endl;
349 
350 
351  //return;
352 }
TVector3 pos
Base Class for genfit track representations. Defines interface for track parameterizations.
Definition: GFAbsTrackRep.h:80
TTree * b
Detector plane genfit geometry class.
Definition: GFDetPlane.h:59
Track object for genfit. genfit algorithms work on these objects.
Definition: GFTrack.h:60
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
std::map< unsigned int, TClonesArray * > fHitBranchMap
unsigned int getNHits() const
Definition: GFTrackCand.h:113
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
void setCandidate(const GFTrackCand &cand, bool doreset=false)
set track candidate
Definition: GFTrack.cxx:148
Double_t mom
Definition: plot_dirc.C:14
TGeoManager * gGeoManager
__m128 v
Definition: P4_F32vec4.h:4
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
Track candidate – a list of cluster indices.
Definition: GFTrackCand.h:55
GFTrack * trk
Definition: checkgenfit.C:13
void addHit(unsigned int detId, unsigned int hitId, double rho=0., unsigned int planeId=0)
Definition: GFTrackCand.cxx:44
TFile * fb
PndSdsMCPoint * hit
Definition: anasim.C:70
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
Singleton which provides access to magnetic field for track representations.
PndSdsMCPoint * point
Definition: anaLmdCluster.C:72
Int_t PndHypDPatternRecoTask::GetChargeIon ( Int_t  ion)
private

Definition at line 354 of file PndHypDPatternRecoTask.cxx.

References L, and Z.

Referenced by Exec().

355 {
356  Int_t A,Z,L;
357 
358  if(ion>1000000000&&(ion<1010000000))
359  { ion -= 1000000000;
360  Z = ion/10000;
361  ion -= 10000*Z;
362  A = ion/10;
363  cout<<" ion charge "<<Z<<endl;
364 
365  return Z;
366 
367  }
368 if((ion>1010000000||ion>1020000000))
369  {
370  ion -= 1000000000;
371  L = ion/10000000;
372  ion -= 10000000*L;
373  Z = ion/10000;
374  ion -= 10000*Z;
375  A = ion/10;
376  cout<<L<<" hypernuclei charge "<<Z<<endl;
377 
378  return Z;
379 
380  }
381 }
double Z
Definition: anaLmdDigi.C:68
InitStatus PndHypDPatternRecoTask::Init ( )
virtual

Definition at line 66 of file PndHypDPatternRecoTask.cxx.

References fGeanePro, fHitBranchMap, fHitBranchNameMap, fMcArray, fPersistence, fPointArray, fSdsArray, fTrackArray, fUseMVD, and gGeoManager.

67 {
68  //Get ROOT Manager
69  FairRootManager* ioman= FairRootManager::Instance();
70 
71  if(ioman==0)
72  {
73  Error("PndHypDPatternRecoTask::Init","RootManager not instantiated!");
74  return kERROR;
75  }
76 
77  // open hit arrays
78  std::map<unsigned int,TString>::iterator iter=fHitBranchNameMap.begin();
79  while(iter!=fHitBranchNameMap.end()){
80 
81  TClonesArray* ar=(TClonesArray*) ioman->GetObject(iter->second);
82 
83  if(ar==0){
84  Error("PndHypDPRTask::Init","point-array %s not found!",iter->second.Data());
85  }
86  else{
87  std::cout<<" hit array "<<iter->second<<std::endl;
88  fHitBranchMap[iter->first] = ar;
89  }
90 
91  ++iter;
92  }//end loops over hit types
93 
94 
95  // open MCTruth array
96  fMcArray=(TClonesArray*) ioman->GetObject("MCTrack");
97  if(fMcArray==0){
98  Error("PndHypDPRTask::Init","mctrack-array not found!");
99  return kERROR;
100  }
101 
102  fPointArray=(TClonesArray*) ioman->GetObject("HypPoint");
103  if(fPointArray==0){
104  Error("PndHypDPRTask::Init","hit-array not found!");
105  return kERROR;
106  }
107 
108  if(fUseMVD==true){
109  fSdsArray=(TClonesArray*) ioman->GetObject("MVDPoint");
110  if(fSdsArray==0){
111  Error("PndHypDPRTask::Init","hit-array not found!");
112  return kERROR;
113  }
114  }
115 
116  // create and register output array
117  fTrackArray = new TClonesArray("GFTrack");
118  ioman->Register("Track","GenFit",fTrackArray,fPersistence);
119 
120 
121  // GeanePro will get Geometry and BField from the Run
122  fGeanePro=new FairGeanePro();
123 
124  std::cout << "-I- gGeoManager = "<<gGeoManager << std::endl;
125 
126  return kSUCCESS;
127 }
std::map< unsigned int, TClonesArray * > fHitBranchMap
TGeoManager * gGeoManager
std::map< unsigned int, TString > fHitBranchNameMap
void PndHypDPatternRecoTask::SetField ( FairField *  f)
inline

Definition at line 55 of file PndHypDPatternRecoTask.h.

References f, and fField.

55 {fField=f;}
TFile * f
Definition: bump_analys.C:12
void PndHypDPatternRecoTask::SetHitFL ( Bool_t  opt)
inline

Definition at line 56 of file PndHypDPatternRecoTask.h.

References fMCvalue.

void PndHypDPatternRecoTask::SetPersistence ( Bool_t  opt = kTRUE)
inline

Definition at line 54 of file PndHypDPatternRecoTask.h.

References fPersistence.

void PndHypDPatternRecoTask::SetVtxAbsName ( TString  name)
inline

Definition at line 49 of file PndHypDPatternRecoTask.h.

References fVtx, fVtxName, and name.

49  {
50  fVtx = true;
51  fVtxName=name;
52  }
TString name
void PndHypDPatternRecoTask::UseGeane ( bool  f = true)
inline

Definition at line 57 of file PndHypDPatternRecoTask.h.

References f, and fUseGeane.

57 {fUseGeane=f;}
TFile * f
Definition: bump_analys.C:12
void PndHypDPatternRecoTask::UseMVD ( bool  mvd)
inline

Definition at line 58 of file PndHypDPatternRecoTask.h.

References fUseMVD.

Member Data Documentation

int PndHypDPatternRecoTask::fEventNr
private

Definition at line 81 of file PndHypDPatternRecoTask.h.

Referenced by Exec().

FairField* PndHypDPatternRecoTask::fField
private

Definition at line 92 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and SetField().

FairGeanePro* PndHypDPatternRecoTask::fGeanePro
private

Definition at line 93 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndHypDPatternRecoTask::fHitArray
private

Definition at line 74 of file PndHypDPatternRecoTask.h.

Referenced by Exec().

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

Definition at line 80 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and Init().

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

Definition at line 79 of file PndHypDPatternRecoTask.h.

Referenced by AddHitBranch(), and Init().

TClonesArray* PndHypDPatternRecoTask::fMcArray
private

Definition at line 75 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and Init().

Bool_t PndHypDPatternRecoTask::fMCvalue
private

Definition at line 87 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and SetHitFL().

Bool_t PndHypDPatternRecoTask::fPersistence
private

Definition at line 84 of file PndHypDPatternRecoTask.h.

Referenced by Init(), and SetPersistence().

TClonesArray* PndHypDPatternRecoTask::fPointArray
private

Definition at line 71 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndHypDPatternRecoTask::fSdsArray
private

Definition at line 72 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and Init().

GFRecoHitFactory* PndHypDPatternRecoTask::fTheRecoHitFactory
private

Definition at line 90 of file PndHypDPatternRecoTask.h.

TClonesArray* PndHypDPatternRecoTask::fTrackArray
private

Definition at line 73 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and Init().

Bool_t PndHypDPatternRecoTask::fUseGeane
private

Definition at line 85 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and UseGeane().

bool PndHypDPatternRecoTask::fUseMVD
private

Definition at line 86 of file PndHypDPatternRecoTask.h.

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

Bool_t PndHypDPatternRecoTask::fVtx
private

Definition at line 88 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and SetVtxAbsName().

TString PndHypDPatternRecoTask::fVtxName
private

Definition at line 83 of file PndHypDPatternRecoTask.h.

Referenced by Exec(), and SetVtxAbsName().


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