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

#include <CbmGeaneTrT.h>

Inheritance diagram for FairGeaneTrT:

Public Member Functions

 FairGeaneTrT ()
 
 ~FairGeaneTrT ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
Bool_t CoordSDToMARS (TVector3 o, TVector3 y, TVector3 z, TMatrixT< double > coor, TVector3 &coordinate)
 

Private Member Functions

 ClassDef (FairGeaneTrT, 1)
 

Private Attributes

TClonesArray * fHitArray
 
TClonesArray * fPointArray
 
TClonesArray * fTrackArray
 
TTree * t
 
TFile * f
 
TClonesArray * fTrackParIni
 
TClonesArray * fTrackParGeane
 
TClonesArray * fTrackParFinal
 
TClonesArray * fTrackParMC
 
TGeant3 * gMC3
 
Int_t fEvent
 
FairGeanePro * fPro
 
FairGeaneUtil * fUtil
 

Detailed Description

Definition at line 21 of file CbmGeaneTrT.h.

Constructor & Destructor Documentation

FairGeaneTrT::FairGeaneTrT ( )

Default constructor

Definition at line 31 of file CbmGeaneTrT.cxx.

31  :
32  FairTask("Test") { }
FairGeaneTrT::~FairGeaneTrT ( )

Destructor

Definition at line 38 of file CbmGeaneTrT.cxx.

38 { }

Member Function Documentation

FairGeaneTrT::ClassDef ( FairGeaneTrT  ,
 
)
private
Bool_t FairGeaneTrT::CoordSDToMARS ( TVector3  o,
TVector3  y,
TVector3  z,
TMatrixT< double >  coor,
TVector3 &  coordinate 
)
void FairGeaneTrT::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 93 of file CbmGeaneTrT.cxx.

References CAMath::Abs(), Bool_t, Double_t, fHitArray, fParticle, fPointArray, fPro, fTrackArray, fTrackParFinal, fTrackParGeane, fTrackParIni, fUtil, PndSttTrack::GetFlag(), PndSttPoint::GetXInLocal(), PndSttPoint::GetXtot(), PndSttPoint::GetYInLocal(), PndSttPoint::GetYtot(), PndSttPoint::GetZtot(), CAMath::Sqrt(), and v.

93  {
94  // cout << "FairGeaneTrT::Exec" << endl;
95  fTrackParGeane->Delete();
96  fTrackParIni->Delete();
97  fTrackParFinal->Delete();
98 
99  if (!fTrackArray ) return;
100 
101 
102  Int_t nTracks = fTrackArray->GetEntriesFast();
103 
104  for (Int_t iTrack = 0; iTrack < nTracks; iTrack++)
105  {
106 
107  PndSttTrack* pTrack = (PndSttTrack*) fTrackArray->At(iTrack);
108  if(!pTrack) continue;
109  if(pTrack->GetFlag() < 3) continue; // per ora voglio solo le tracce perfette
110 
111 
112  Int_t iHit1 = pTrack->GetHitIndex(0);
113  PndSttHit *currenthit1 = (PndSttHit*) fHitArray->At(iHit1);
114  if(!currenthit1) continue;
115  if(currenthit1->GetXint() == -999 || currenthit1->GetYint() == -999) continue;
116  Int_t refindex1 = currenthit1->GetRefIndex();
117  // get point
118  PndSttPoint *fPoint1 = (PndSttPoint*) fPointArray->At(refindex1);
119 
120  TVector3 StartPos = TVector3 (fPoint1->GetXtot(),fPoint1->GetYtot(),fPoint1->GetZtot());
121  // TVector3 StartPos = TVector3(0,0,0);
122  TVector3 StartPosErr = TVector3(0,0,0);
123  TVector3 StartMom = TVector3 (fPoint1->GetPXtot(),fPoint1->GetPYtot(),fPoint1->GetPZtot());
124  // TVector3 StartMom = TVector3(1.,0.1,0.1);
125  TVector3 StartMomErr = TVector3(0,0,0);
126 
127 
128  // tracking from vertex ..........................
129  // using prefit result from vertex
130  // PndSttHelixTrackFitter *fFitter;
131  // TVector3 *vertex = new TVector3(0.,0.,0.); // CHECK to set the vertex in (0.,0.,0.)
132  // TVector3 *recovtx = fFitter->PCAToPoint(pTrack, vertex);
133  // TVector3 *momeatvtx;
134  // momeatvtx = fFitter->MomentumAtPoint(pTrack, recovtx);
135  // // ERRORS ON VERTEX FROM checkMomentumPCA.C
136  // // MOMENTUM 0.018, 0.018, 0.04
137  // // POSITION 0.007, 0.007, 0.46
138  // TVector3 StartPos = *recovtx;
139  // TVector3 StartPosErr = TVector3(0.018,0.018,0.04);
140  // TVector3 StartMom = *momeatvtx;
141  // TVector3 StartMomErr = TVector3(0.007,0.007,0.46);
142 // dstvtx->Fill(StartPos.Mag());
143 // momvtx->Fill(StartMom.Mag() - 1.5);
144 // if(StartPos.Mag() > 1.) continue; // cut on distance from vertex
145 // if((StartMom.Mag() - 1.5) > 0.05) continue; // cut on momentum at vertex
146  // ..................................................
147 
148  //tracking up to the last hit on track on the PCA to the wire
149  Int_t hitcounter = pTrack->GetNofHits();
150  Int_t iHit = pTrack->GetHitIndex(hitcounter - 1);
151  PndSttHit *currenthit = (PndSttHit*) fHitArray->At(iHit);
152  if(!currenthit) continue;
153  if(currenthit->GetXint() == -999 || currenthit->GetYint() == -999) continue;
154  Int_t refindex = currenthit->GetRefIndex();
155  // get point
156  PndSttPoint *iPoint = (PndSttPoint*) fPointArray->At(refindex);
157  if(TMath::Sqrt(iPoint->GetXInLocal()*iPoint->GetXInLocal() + iPoint->GetYInLocal()*iPoint->GetYInLocal()) < 0.4999) continue;
158 
159  // out
160  // let' s take (approx!) the medium point
161  // position --> (posout - posin) /2
162  // momentum --> (momout - momin) / 2
163  TVector3 EndPos = TVector3(iPoint->GetXtot(), iPoint->GetYtot(), iPoint->GetZtot());
164  TVector3 EndPosErr = TVector3(0,0,0);
165  TVector3 EndMom = TVector3(iPoint->GetPXtot(), iPoint->GetPYtot(), iPoint->GetPZtot());
166  TVector3 EndMomErr = TVector3(0,0,0);
167 
168  Int_t PDGCode= 13;
169  TDatabasePDG *fdbPDG= TDatabasePDG::Instance();
170  TParticlePDG *fParticle= fdbPDG->GetParticle(PDGCode);
171  Double_t fCharge= fParticle->Charge();
172 
173  TClonesArray& clref1 = *fTrackParIni;
174  Int_t size1 = clref1.GetEntriesFast();
175  // FairTrackParH *fStart = new (clref1[size1]) FairTrackParH(StartPos, StartMom, StartPosErr, StartMomErr, fCharge);
176  FairTrackParP *fStart = new (clref1[size1]) FairTrackParP(StartPos, StartMom, StartPosErr, StartMomErr, fCharge, TVector3(0.,0.,0.), TVector3(0.,1.,0.), TVector3(0.,0.,1.));
177 
178  TClonesArray& clref = *fTrackParGeane;
179  Int_t size = clref.GetEntriesFast();
180  FairTrackParP *fRes = new(clref[size]) FairTrackParP();
181 
182  TClonesArray& clref2 = *fTrackParFinal;
183  Int_t size2 = clref2.GetEntriesFast();
184  FairTrackParP *fFinal = new(clref2[size2]) FairTrackParP();
185  FairTrackParH *fFinalH = new FairTrackParH(EndPos, EndMom, EndPosErr, EndMomErr, fCharge);
186 
187  // ----- propagation: using method propagate to closest ------
188  TVector3 v0 = TVector3(0, 0, 0);
189  fPro->SetPoint(v0);
190 
191  TVector3 center = TVector3(currenthit->GetX(), currenthit->GetY(), currenthit->GetZ());
192  TVector3 wiredir = currenthit->GetWireDirection();
193  wiredir *= 75.;
194  TVector3 second = center + wiredir;
195  TVector3 first = center - wiredir;
196 
197  fPro->SetWire(first, second);
198  fPro->PropagateToVirtualPlaneAtPCA(2); // 1 if point; 2 if wire
199 
200  Bool_t rc = fPro->Propagate(fStart, fRes, PDGCode);
201 
202  if(rc == kTRUE) {
203  Double_t PD[3], RD[15], H[3], CH, SPU, DJ[3], DK[3], PC[3], RC[15];
204  Int_t IERR;
205 
206  PC[0] = fFinalH->GetQp()/fFinalH->GetQ();
207  PC[1] = fFinalH->GetLambda();
208  PC[2] = fFinalH->GetPhi();
209 
210  fFinalH->GetCov(RC);
211 
212  H[0] = 0; H[1] = 0; H[2] = 20;
213  CH = fRes->GetQ()/TMath::Abs(fRes->GetQ());
214 
215  DJ[0] = fRes->GetJVer().X();
216  DJ[1] = fRes->GetJVer().Y();
217  DJ[2] = fRes->GetJVer().Z();
218  DK[0] = fRes->GetKVer().X();
219  DK[1] = fRes->GetKVer().Y();
220  DK[2] = fRes->GetKVer().Z();
221 
222  fUtil->FromSCToSD(PC, RC, H, CH, DJ, DK, IERR, SPU, PD, RD);
223 
224  TVector3 positionsd = fUtil->FromMARSToSDCoord(TVector3(fFinalH->GetX(), fFinalH->GetY(), fFinalH->GetZ()), fRes->GetOrigin(), fRes->GetIVer(), fRes->GetJVer(), fRes->GetKVer());
225  Double_t u = positionsd.X();
226  Double_t v = positionsd.Y();
227  Double_t w = positionsd.Z();
228 
229  fFinal->SetTrackPar(v, w, PD[1], PD[2], CH*PD[0], RD,
230  fRes->GetOrigin(), fRes->GetIVer(), fRes->GetJVer(),
231  fRes->GetKVer(),SPU);
232  }
233 
234 
235 // -----------------------------------------------------
236 
237 
238 
239  }
240  fTrackArray->Delete();
241  fHitArray->Delete();
242  fPointArray->Delete();
243 
244  // can2->cd(1);
245 // dstvtx->Draw();
246 // can2->cd(2);
247 // momvtx->Draw();
248 
249 }
FairGeaneUtil * fUtil
Definition: CbmGeaneTrT.h:57
TClonesArray * fTrackParFinal
Definition: CbmGeaneTrT.h:50
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t GetXInLocal() const
Definition: PndSttPoint.h:49
__m128 v
Definition: P4_F32vec4.h:4
Double_t GetYInLocal() const
Definition: PndSttPoint.h:50
static T Abs(const T &x)
Definition: PndCAMath.h:39
TClonesArray * fPointArray
Definition: CbmGeaneTrT.h:42
Int_t * fParticle
Definition: run_full.C:24
Double_t GetYtot() const
Definition: PndSttPoint.h:55
Double_t GetZtot() const
Definition: PndSttPoint.h:56
Double_t
TClonesArray * fTrackParIni
Definition: CbmGeaneTrT.h:48
FairGeanePro * fPro
Definition: CbmGeaneTrT.h:56
TClonesArray * fHitArray
Definition: CbmGeaneTrT.h:41
Double_t GetXtot() const
Definition: PndSttPoint.h:54
TClonesArray * fTrackParGeane
Definition: CbmGeaneTrT.h:49
Int_t GetFlag() const
Definition: PndSttTrack.h:51
TClonesArray * fTrackArray
Definition: CbmGeaneTrT.h:43
InitStatus FairGeaneTrT::Init ( )
virtual

Virtual method Init

Definition at line 43 of file CbmGeaneTrT.cxx.

References fHitArray, fPointArray, fPro, fTrackArray, fTrackParFinal, fTrackParGeane, fTrackParIni, and fUtil.

43  {
44 
45  // Get RootManager
46  FairRootManager* ioman = FairRootManager::Instance();
47  if ( ! ioman ) {
48  cout << "-E- FairGeaneTrT::Init: "
49  << "RootManager not instantised!" << endl;
50  return kFATAL;
51  }
52 
53  // can2 = new TCanvas();
54 // can2->Divide(1,2);
55 // dstvtx = new TH1F("dstvtx","dstvtx",100,0,10);
56 // momvtx = new TH1F("momvtx","momvtx",100,0,0.5);
57 
58 
59  fTrackParGeane = new TClonesArray("FairTrackParP");
60  ioman->Register("GeaneTrackPar","Geane", fTrackParGeane, kTRUE);
61 
62  fTrackParIni = new TClonesArray("FairTrackParP");
63  ioman->Register("GeaneTrackIni","Geane", fTrackParIni, kTRUE);
64 
65  fTrackParFinal = new TClonesArray("FairTrackParP");
66  ioman->Register("GeaneTrackFinal","Geane", fTrackParFinal, kTRUE);
67 
68  // Get hit Array
69  fHitArray = (TClonesArray*) ioman->GetObject("STTHit");
70  if (!fHitArray) cout << "-W- FairGeaneTrT::Init: No Hit array!" << endl;
71  // Get point Array
72  fPointArray = (TClonesArray*) ioman->GetObject("STTPoint");
73  if (!fPointArray) cout << "-W- FairGeaneTrT::Init: No Point array!" << endl;
74 
75  // Get SttTrack array
76  fTrackArray = (TClonesArray*) ioman->GetObject("STTTrack");
77  if (!fTrackArray) {
78  cout << "-E- FairGeaneTrT::Init: No SttTrack array!" << endl;
79  return kERROR;
80  }
81 
82  fPro = new FairGeanePro();
83  fUtil = new FairGeaneUtil();
84  // fPro->PropagateToVolume("mL3",0,1);
85 
86  return kSUCCESS;
87 
88 }
FairGeaneUtil * fUtil
Definition: CbmGeaneTrT.h:57
TClonesArray * fTrackParFinal
Definition: CbmGeaneTrT.h:50
TClonesArray * fPointArray
Definition: CbmGeaneTrT.h:42
TClonesArray * fTrackParIni
Definition: CbmGeaneTrT.h:48
FairGeanePro * fPro
Definition: CbmGeaneTrT.h:56
TClonesArray * fHitArray
Definition: CbmGeaneTrT.h:41
TClonesArray * fTrackParGeane
Definition: CbmGeaneTrT.h:49
TClonesArray * fTrackArray
Definition: CbmGeaneTrT.h:43

Member Data Documentation

TFile* FairGeaneTrT::f
private

Definition at line 46 of file CbmGeaneTrT.h.

Int_t FairGeaneTrT::fEvent
private

Definition at line 55 of file CbmGeaneTrT.h.

TClonesArray* FairGeaneTrT::fHitArray
private

Input array of Points

Definition at line 41 of file CbmGeaneTrT.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrT::fPointArray
private

Definition at line 42 of file CbmGeaneTrT.h.

Referenced by Exec(), and Init().

FairGeanePro* FairGeaneTrT::fPro
private

Definition at line 56 of file CbmGeaneTrT.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrT::fTrackArray
private

Definition at line 43 of file CbmGeaneTrT.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrT::fTrackParFinal
private

Definition at line 50 of file CbmGeaneTrT.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrT::fTrackParGeane
private

Definition at line 49 of file CbmGeaneTrT.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrT::fTrackParIni
private

Output array of Hits

Definition at line 48 of file CbmGeaneTrT.h.

Referenced by Exec(), and Init().

TClonesArray* FairGeaneTrT::fTrackParMC
private

Definition at line 51 of file CbmGeaneTrT.h.

FairGeaneUtil* FairGeaneTrT::fUtil
private

Definition at line 57 of file CbmGeaneTrT.h.

Referenced by Exec(), and Init().

TGeant3* FairGeaneTrT::gMC3
private

Definition at line 53 of file CbmGeaneTrT.h.

TTree* FairGeaneTrT::t
private

Definition at line 45 of file CbmGeaneTrT.h.


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