FairRoot/PandaRoot
PndMvdMSAnaTask.cxx
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- PndMvdMSAnaTask source file -----
3 // ----- Created 18/07/08 by T.Stockmanns -----
4 // -------------------------------------------------------------------------
5 // libc includes
6 #include <iostream>
7 
8 // Root includes
9 #include "TROOT.h"
10 #include "TClonesArray.h"
11 #include "TVector3.h"
12 
13 
14 // framework includes
15 #include "FairRootManager.h"
16 #include "PndMvdMSAnaTask.h"
17 #include "FairRun.h"
18 #include "FairRuntimeDb.h"
19 #include "FairHit.h"
20 #include "PndMCTrack.h"
21 // PndMvd includes
22 #include "PndSdsMCPoint.h"
23 #include "FairTrackParH.h"
24 #include "FairTrackParP.h"
25 
26 #include "TDatabasePDG.h"
27 
28 #include <vector>
29 #include <map>
30 
31 
32 // ----- Default constructor -------------------------------------------
33 PndMvdMSAnaTask::PndMvdMSAnaTask() : FairTask("Geane Task for PANDA PndMvd"),
34  fMCHits(NULL),
35  fMCTracks(NULL),
36  fTrackParGeane(NULL),
37  fTrackParIni(NULL),
38  fTrackParFinal(NULL),
39  fDetName(NULL),
40  fPro(NULL),
41  fGeoH(NULL),
42  fEventNr(0),
43  fUseMVDPoint(false),
44  fTrackPixHitIdMap(),
45  fTrackStripHitIdMap()
46 {
48 }
49 // -------------------------------------------------------------------------
50 
51 
52 // ----- Destructor ----------------------------------------------------
54 {
55 }
56 
57 // ----- Public method Init --------------------------------------------
59 {
60  // Get RootManager
61  FairRootManager* ioman = FairRootManager::Instance();
62  if ( !ioman){
63  std::cout << "-E- PndMvdMSAnaTask::Init: "<< "RootManager not instantiated!" << std::endl;
64  return kFATAL;
65  }
66 
67  // Get input arrays
68  fMCTracks = (TClonesArray*) ioman->GetObject("MCTrack");
69  if (!fMCTracks){
70  std::cout << "-W- PndMvdMSAnaTask::Init: "<< "No MCTrack" << " array!" << std::endl;
71  return kERROR;
72  }
73 
74  fMCHits = (TClonesArray*) ioman->GetObject("MVDPoint");
75  if ( !fMCHits) {
76  std::cout << "-W- PndMvdMSAnaTask::Init: "<< "No MVDPoint"<<" array!" << std::endl;
77  return kERROR;
78  }
79 
80  fTrackParGeane = new TClonesArray("FairTrackParP");
81  ioman->Register("GeaneTrackPar","Geane", fTrackParGeane, kTRUE);
82 
83  fTrackParIni = new TClonesArray("FairTrackParP");
84  ioman->Register("GeaneTrackIni","Geane", fTrackParIni, kTRUE);
85 
86  fTrackParFinal = new TClonesArray("FairTrackParP");
87  ioman->Register("GeaneTrackFinal","Geane", fTrackParFinal, kTRUE);
88 
89  fDetName = new TClonesArray("TObjString");
90  ioman->Register("DetName", "Geane", fDetName, kTRUE);
91 
92  fPro = new FairGeanePro();
93  return kSUCCESS;
94 }
95 // -------------------------------------------------------------------------
97 {
98 }
99 
100 
101 // ----- Public method Exec --------------------------------------------
102 void PndMvdMSAnaTask::Exec(Option_t*)
103 {
104  std::map<int, std::vector<int> > mcHitMap; //Track -> MCHits
105  fTrackParGeane->Delete();
106  fTrackParIni->Delete();
107  fTrackParFinal->Delete();
108  fDetName->Delete();
109 
110  mcHitMap = AssignHitsToTracks();
111 
112  std::cout << "------------Event " << fEventNr << "-------------" << std::endl;
113  fEventNr++;
114 
115  for (std::map<int, std::vector<int> >::const_iterator it = mcHitMap.begin(); it!= mcHitMap.end(); it++){ //go through all tracks
116  PndMCTrack* myTrack = (PndMCTrack*)(fMCTracks->At(it->first));
117  std::vector<int> MChits = it->second;
118  if (MChits.size() > 1 && myTrack->GetMotherID() < 0){
119  TVector3 StartPos, StartPosErr, StartMom, StartMomErr, StartO, StartU, StartV;
120  int p = 0;
121  if (fUseMVDPoint){
122  PndSdsMCPoint* startPoint = (PndSdsMCPoint*)(fMCHits->At(0));
123  std::cout << "StartPoint: " << *startPoint << std::endl;
124  StartPos = (startPoint->GetPosition());
125  StartPosErr = TVector3(0.0, 0.0, 0.0);
126  StartMom = TVector3(startPoint->GetPx(), startPoint->GetPy(), startPoint->GetPz());
127  StartMomErr = TVector3(0.0, 0.0, 0.0);
128  fGeoH->GetOUVShortId(startPoint->GetSensorID(), StartO, StartU, StartV);
129  p = 1;
130  }
131  else{
132  StartPos = myTrack->GetStartVertex();
133  StartPosErr = TVector3(0.0, 0.0, 0.0);
134  StartMom = myTrack->GetMomentum();
135  StartMomErr = TVector3(0.0, 0.0, 0.0);
136  StartO = myTrack->GetStartVertex();
137  StartU = TVector3(1.0, 0.0, 0.0);
138  StartV = TVector3(0.0, 1.0, 0.0);
139  p = 0;
140  }
141 
142  fPro->PropagateFromPlane(StartU, StartV);
143 
144  Int_t PDGCode = myTrack->GetPdgCode();
145  TDatabasePDG *fdbPDG= TDatabasePDG::Instance();
146  TParticlePDG *fParticle= fdbPDG->GetParticle(PDGCode);
147  Int_t fCharge= (Int_t)fParticle->Charge();
148 
149  TClonesArray& clref1 = *fTrackParIni;
150  Int_t size1 = clref1.GetEntriesFast();
151  FairTrackParP *fStart= new (clref1[size1]) FairTrackParP(StartPos, StartMom, StartPosErr, StartMomErr, fCharge, StartO, StartU, StartV);
152 
153  for (; p < (int)MChits.size(); p++){ //go through all hits in track
154  PndSdsMCPoint* myPoint = (PndSdsMCPoint*)(fMCHits->At(MChits[p]));
155  std::cout << "PropagationPoint: " << *myPoint << std::endl;
156 
157  TVector3 StopPos(myPoint->GetPosition());
158  TVector3 StopPosErr(0.0, 0.0, 0.0);
159  TVector3 StopMom(myPoint->GetPx(), myPoint->GetPy(), myPoint->GetPz());
160  TVector3 StopMomErr(0.0, 0.0, 0.0);
161 
162  TVector3 o, u, v;
163  fGeoH->GetOUVShortId(myPoint->GetSensorID(), o, u, v);
164 
165  TClonesArray& clref2 = *fTrackParFinal;
166  Int_t size2 = clref2.GetEntriesFast();
167  //FairTrackParP *fStop=
168  new (clref2[size2]) FairTrackParP(StopPos, StopMom, StopPosErr, StopMomErr, fCharge, o, u, v);
169 
170  TClonesArray& clref = *fTrackParGeane;
171  Int_t size = clref.GetEntriesFast();
172  FairTrackParP *fRes= new(clref[size]) FairTrackParP();
173 
174  //std::cout << "DetName: " << fGeoH->GetPath(myPoint->GetDetName()) << std::endl;
175  TClonesArray& cDetRef = *fDetName;
176  Int_t size3 = cDetRef.GetEntriesFast();
177  new (cDetRef[size3]) TObjString(fGeoH->GetPath(myPoint->GetSensorID()));
178 
179  std::cout << "Propagation Plane:" << std::endl;
180  std::cout << "o: " << o[0] << " " << o[1] << " " << o[2] << std::endl;
181  std::cout << "u: " << u[0] << " " << u[1] << " " << u[2] << std::endl;
182  std::cout << "v: " << v[0] << " " << v[1] << " " << v[2] << std::endl;
183 
184  TVector3 endPoint(myPoint->GetPosition());
185  //fPro->SetPoint(endPoint);
186  //fPro->PropagateToPCA(1);
187  fPro->PropagateToPlane(o, u, v);
188  fPro->Propagate(fStart, fRes, PDGCode);
189 
190  std::cout << std::endl;
191  std::cout << "Propagation Points: " << std::endl;
192  std::cout << fRes->GetX() << ", ";//" +/- " << fRes->GetDX() << std::endl;
193  std::cout << fRes->GetY() << ", ";//" +/- " << fRes->GetDY() << std::endl;
194  std::cout << fRes->GetZ(); //" +/- " << fRes->GetDZ() << std::endl;
195  std::cout << std::endl;
196 
197  std::cout << "Propagation Momentum: " << std::endl;
198 
199  std::cout << fRes->GetPx() << ", "; //" +/- " << fRes->GetDPx() << std::endl;
200  std::cout << fRes->GetPy() << ", "; //" +/- " << fRes->GetDPy() << std::endl;
201  std::cout << fRes->GetPz() << std::endl;//<< " +/- " << fRes->GetDPz() << std::endl;
202 
203  std::cout << "Charge: " << fRes->GetQ() << std::endl;
204 
205  std::cout << "Porpagation Covariance Matrix: " << std::endl;
206 
207  Double_t CovMatrix[15];
208  Double_t CovMatrix66[6][6];
209  fRes->GetCov(CovMatrix);
210  fRes->GetMARSCov(CovMatrix66);
211 
212  /* for (int i = 0; i < 5; i++){
213  for (int j = 0; j < 5-i; j++){
214  std::cout << CovMatrix[i*5+j] << " ";
215  }
216  std::cout << std::endl;
217  }
218  */
219  std::cout << CovMatrix66[3][3] << ", " << CovMatrix66[3][4] << ", " << CovMatrix66[3][5] << std::endl;
220  std::cout << CovMatrix66[4][4] << ", " << CovMatrix66[4][5] << ", " << CovMatrix66[5][5] << std::endl;
221  std::cout << CovMatrix66[0][3] << ", " << CovMatrix66[1][3] << ", " << CovMatrix66[2][3] << std::endl;
222  std::cout << CovMatrix66[0][4] << ", " << CovMatrix66[1][4] << ", " << CovMatrix66[2][4] << std::endl;
223  std::cout << CovMatrix66[0][5] << ", " << CovMatrix66[1][5] << ", " << CovMatrix66[2][5] << std::endl;
224  std::cout << CovMatrix66[0][0] << ", " << CovMatrix66[0][1] << ", " << CovMatrix66[0][2] << std::endl;
225  std::cout << CovMatrix66[1][1] << ", " << CovMatrix66[1][2] << ", " << CovMatrix66[2][2] << std::endl;
226  std::cout << std::endl;
227 
228  /* for (int i = 0; i < 6; i++){
229  for (int j = 0; j < 6; j++){
230  std::cout << CovMatrix66[i][j] << " ";
231  }
232  std::cout << std::endl;
233  }
234  */
235  TVector3 global(fRes->GetX(), fRes->GetY(), fRes->GetZ());
236  TVector3 local = fGeoH->MasterToLocalShortId(global, myPoint->GetSensorID());
237 
238  std::cout << "Propagation Point local: " << std::endl;
239  std::cout << local[0] << " " << local[1] << " " << local[2] << std::endl;
240  std::cout << std::endl;
241  }
242  }
243  }
244  fMCTracks->Delete();
245  fMCHits->Delete();
246 }
247 
249 {
250 }
251 
252 std::map<int, std::vector<int> > PndMvdMSAnaTask::AssignHitsToTracks()
253 {
254  std::map<int, std::vector<int> > result;
255  for (int i = 0; i < fMCHits->GetEntriesFast(); i++){ //get all MC Hits
256  PndSdsMCPoint* myPoint = (PndSdsMCPoint*)(fMCHits->At(i)); //sort MCHits with Tracks
257  //PndMCTrack* myTrack = (PndMCTrack*)(fMCTracks->At(myPoint->GetTrackID()));
258  result[myPoint->GetTrackID()].push_back(i);
259 
260  }
261  return result;
262 }
263 
Double_t p
Definition: anasim.C:58
void GetOUVShortId(Int_t shortId, TVector3 &o, TVector3 &u, TVector3 &v)
Int_t i
Definition: run_full.C:25
PndGeoHandling * fGeoH
Definition: anasim.C:34
PndTransMap * map
Definition: sim_emc_apd.C:99
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
TClonesArray * fTrackParGeane
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
Int_t GetSensorID() const
Definition: PndSdsMCPoint.h:89
TClonesArray * fMCHits
TString GetPath(Int_t shortID)
for a given shortID the path is returned
__m128 v
Definition: P4_F32vec4.h:4
FairGeanePro * fPro
TClonesArray * fTrackParIni
PndGeoHandling * fGeoH
Int_t * fParticle
Definition: run_full.C:24
fTrackParFinal
Definition: runPULL1.C:25
TClonesArray * fMCTracks
Double_t
virtual void SetParContainers()
static PndGeoHandling * Instance()
TClonesArray * fTrackParFinal
TVector3 MasterToLocalShortId(const TVector3 &master, const Int_t &shortId)
TVector3 GetPosition() const
Definition: PndSdsMCPoint.h:90
TGeoHMatrix * global
fTrackParGeane
Definition: runPULL1.C:23
virtual void Finish()
ClassImp(PndAnaContFact)
TClonesArray * fDetName
std::map< int, std::vector< int > > AssignHitsToTracks()
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
virtual void Exec(Option_t *opt)
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
virtual InitStatus Init()
fTrackParIni
Definition: runPULL1.C:24