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

#include <PndLmdTrkQTask.h>

Inheritance diagram for PndLmdTrkQTask:

Public Member Functions

 PndLmdTrkQTask (const PndLmdTrkQTask &)=delete
 
 PndLmdTrkQTask (Double_t pBeam=0, TString geaneBranch="GeaneTrackFinal", TString trackBranch="LMDPndTrackFilt")
 Set up beam momuntum value. More...
 
 ~PndLmdTrkQTask ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 
void SetWriteMC (bool wr)
 

Private Member Functions

 ClassDef (PndLmdTrkQTask, 2)
 

Private Attributes

TString fGeaneName
 
TString fTrackBranchName
 
TClonesArray * fMCHits
 
TClonesArray * fMCTracks
 
TClonesArray * fRecHits
 
TClonesArray * fRecCandTracks
 
TClonesArray * fRecTracks
 
TClonesArray * fRecBPTracks
 
TClonesArray * fTrackQ
 
int fEventNr
 
Double_t fPbeam
 
bool fWriteAllMC
 

Detailed Description

Definition at line 31 of file PndLmdTrkQTask.h.

Constructor & Destructor Documentation

PndLmdTrkQTask::PndLmdTrkQTask ( const PndLmdTrkQTask )
delete

Default constructor

PndLmdTrkQTask::PndLmdTrkQTask ( Double_t  pBeam = 0,
TString  geaneBranch = "GeaneTrackFinal",
TString  trackBranch = "LMDPndTrackFilt" 
)

Set up beam momuntum value.

Definition at line 47 of file PndLmdTrkQTask.cxx.

References fGeaneName, fPbeam, fTrackBranchName, and fWriteAllMC.

48  :
49  FairTask("Track Quality Task for PANDA Lmd"), fEventNr(0) {
50  fWriteAllMC = false;
51  fPbeam = pBeam;
52  fGeaneName = geaneBranch;
53  fTrackBranchName = trackBranch;
54  //cout<<"Beam Momentum for particle with PDGid#"<<fPDGid<<" this run is "<<fPbeam<<endl;
55  // vtx = IP;
56  // cout<<"Interaction Point:"<<endl;
57  // vtx.Print();
58 }
TString fTrackBranchName
TString fGeaneName
PndLmdTrkQTask::~PndLmdTrkQTask ( )

Destructor

Definition at line 61 of file PndLmdTrkQTask.cxx.

61  {
62 }

Member Function Documentation

PndLmdTrkQTask::ClassDef ( PndLmdTrkQTask  ,
 
)
private
void PndLmdTrkQTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Set signal/bkg flag -—————————————————


Chech how many MC hits has MC trk -—————————


Chech how many REC hits has MC trk -—————————


Compare trk-cand vs. trk-fit –—————————

Obtain first approximation


Set MC ID for each track ----------------————————————

Matching between MC & Rec on hits level--———————————


Counting number of diff MC ids -------------—————————


Set track quality: good or ghost? and assign MC id to rec trk ———


(II) missed tracks are defined on hits information GHOST -------------------------------------------------------------------------------------------------------------------------------------------------------————— check repeated assignment to one MC trk (maybe one MC trk rec. trk was assigned to several REC trk) ------------------—————

END check repeated assignment to one MC trk ------------------------------------------------------------------------------------------------------————

END GHOST-----------------------------------------------------------------------------------------------------------------------------------------------------------———

MISSED ------------------------------------------------------------------------------------------------------------------------------------------------------------------——

if MC trk was missed, justify why

Definition at line 149 of file PndLmdTrkQTask.cxx.

References Double_t, fabs(), fEventNr, fMCHits, fMCTracks, fPbeam, fRecBPTracks, fRecCandTracks, fRecHits, fRecTracks, fTrackQ, fVerbose, fWriteAllMC, PndTrack::GetChi2(), PndTrackCandHit::GetHitId(), PndLmdGeometryHelper::getInstance(), PndMCTrack::GetMomentum(), PndMCTrack::GetMotherID(), PndTrackCand::GetNHits(), PndTrack::GetParamFirst(), PndMCTrack::GetPdgCode(), PndSdsMCPoint::GetPosition(), PndSdsHit::GetPosition(), PndTrack::GetRefIndex(), PndSdsMergedHit::GetSecondMCHit(), PndSdsHit::GetSensorID(), PndTrackCand::GetSortedHit(), PndMCTrack::GetStartVertex(), m, MCPoint, n, nk, Pi, PndLmdTrackQ::SetEvMCMulti(), PndLmdTrackQ::SetEvRECMulti(), PndLmdTrackQ::SetEvTime(), PndLmdTrackQ::SetHalf(), PndLmdTrackQ::SetIPerrmom(), PndLmdTrackQ::SetIPerrpoint(), PndLmdTrackQ::SetIPmom(), PndLmdTrackQ::SetIPpoint(), PndLmdTrackQ::SetLMDchi2(), PndLmdTrackQ::SetLMDdir(), PndLmdTrackQ::SetLMDpoint(), PndLmdTrackQ::SetMCmom(), PndLmdTrackQ::SetMCmomLMD(), PndLmdTrackQ::SetMCpoint(), PndLmdTrackQ::SetMCpointLMD(), PndLmdTrackQ::SetModule(), PndLmdTrackQ::SetNumDoubleMChits(), PndLmdTrackQ::SetNumMChits(), PndLmdTrackQ::SetPDGcode(), PndLmdTrackQ::SetSecondary(), PndLmdTrackQ::SetSumEvPDG(), PndLmdTrackQ::SetTrkRecStatus(), PndLmdTrackQ::SetTrkTime(), and x.

149  {
150  fTrackQ->Delete();
151 
152  Double_t glXrecLMD, glYrecLMD, glZrecLMD, glThetarecLMD, glPhirecLMD;
153  Double_t glXrec, glYrec, glZrec, glThetarec, glPhirec, glMomrec;
154  Double_t glXmc, glYmc, glZmc, glThetamc, glPhimc, glMommc;
155  Double_t glXmcLMD, glYmcLMD, glZmcLMD, glThetamcLMD, glPhimcLMD, glMommcLMD;
156  Int_t trkRECStatus;
157  Int_t trkMCStatus;
158  Double_t glchi2;
159  int glPDG;
160  int glNumMChits;
161  int glNumDoubleMChits;
162 
163  int glModule, glHalf;
164  FairRootManager* ioman = FairRootManager::Instance();
165  // double glEvTime = FairRootManager::Instance()->GetEventTime();
166  double glEvTime = ioman->GetEventTime();
167  // FairMCEventHeader* iomchead = (FairMCEventHeader*) ioman->GetObject("MCEventHeader");
168  // FairMCEventHeader* mcevhead = (FairMCEventHeader*) fMCHeader->At(0);
169  // int glEvID = mcevhead->GetEventID();
170  // cout<<"this is event # "<<glEvID<<endl;
171  double glTrkTime;
172 
173  const int nGeaneTrks = fRecBPTracks->GetEntriesFast();
174  const int nParticles = fMCTracks->GetEntriesFast();
175  const int numTrk = fRecTracks->GetEntriesFast();
176  const int nRecHits = fRecHits->GetEntriesFast();
177  const int nMCHits = fMCHits->GetEntriesFast();
178 
179  const int nTrkCandidates = fRecCandTracks->GetEntriesFast();
180  //const int nRecTrks = fRecTracks->GetEntriesFast(); //[R.K. 01/2017] unused variable
181 
182  if (fVerbose > 0)
183  cout << "%%%%%! Event #" << fEventNr << " has " << nParticles
184  << " true particles, " << " out of it " << nRecHits << " hits, "
185  << nTrkCandidates << " trk-cands, " << numTrk << " tracks and "
186  << nGeaneTrks << " geane Trks!" << endl;
187 
189  int sumID = 0;
190  //int TotCharge=0; //[R.K. 01/2017] unused variable
191  for (Int_t iN = 0; iN < nParticles; iN++) {
192  PndMCTrack *mctrk = (PndMCTrack*) fMCTracks->At(iN);
193  glPDG = mctrk->GetPdgCode();
194  Int_t mcID = mctrk->GetPdgCode();
195  int motherid = mctrk->GetMotherID();
196  if (motherid < 0 && fabs(mcID) < 1e5) {
197  sumID += fabs(mcID);
198  }
199  // TParticlePDG *fParticle = fdbPDG->GetParticle(mcID);
200  // Double_t fCharge = fParticle->Charge();
201  // TotCharge += fCharge;
202  }
204 
206  int MCtksSIMhits[nParticles];
207  int MCDoubleHits[nParticles];
208  for (int imctrk = 0; imctrk < nParticles; imctrk++) {
209  MCtksSIMhits[imctrk] = 0;
210  MCDoubleHits[imctrk] = 0;
211  }
212  for (int imc = 0; imc < nMCHits; imc++) {
213  PndSdsMCPoint* MCPoint = (PndSdsMCPoint*) (fMCHits->At(imc));
214  int MCtrk = MCPoint->GetTrackID();
215  MCtksSIMhits[MCtrk]++;
216  }
218 
220 
221  int MCtksREChits[nParticles];
222  // const int nSize = 2e4;
223  // int MCtksREChits[nSize];//add arbitrary number of fake trks due to noise hits
224  for (int imctrk = 0; imctrk < nParticles; imctrk++) {
225  MCtksREChits[imctrk] = 0;
226  }
227  if (fVerbose > 7)
228  cout << " ** ALL REChits are made from MChits: " << endl; //start MC hit content
229  //bool fMCnegative = false; //[R.K. 01/2017] unused variable
230  for (int irec = 0; irec < nRecHits; irec++) {
231  PndSdsMergedHit* myHit = (PndSdsMergedHit*) (fRecHits->At(irec));
232  int mcrefbot = myHit->GetSecondMCHit();
233  int MCtrkidbot, MCtrkidtop;
234  if (mcrefbot > 0) {
235  PndSdsMCPoint* MCPointBot = (PndSdsMCPoint*) (fMCHits->At(mcrefbot));
236  int MCtrkid = MCPointBot->GetTrackID();
237 
238  MCtksREChits[MCtrkid]++;
239  MCtrkidbot = MCtrkid;
240  if (fVerbose > 7)
241  cout << " " << MCtrkid << "(MChitID=" << mcrefbot << ",z="
242  << MCPointBot->GetZ() << ")";
243  }
244  int mcreftop = myHit->GetRefIndex();
245  if (mcreftop > 0) {
246  PndSdsMCPoint* MCPointTop = (PndSdsMCPoint*) (fMCHits->At(mcreftop));
247  int MCtrkid = MCPointTop->GetTrackID();
248 
249  MCtksREChits[MCtrkid]++;
250  MCtrkidtop = MCtrkid;
251  if (fVerbose > 7)
252  cout << " " << MCtrkid << "(MChitID=" << mcreftop << ",z="
253  << MCPointTop->GetZ() << ")";
254  }
255  if (mcreftop > 0 && mcrefbot > 0) {
256  if (MCtrkidtop == MCtrkidbot) {
257  MCDoubleHits[MCtrkidbot]++;
258  } else {
259  if (fVerbose > 7)
260  cout << " REChit No." << irec << "contain MCid: " << MCtrkidtop
261  << ", " << MCtrkidbot << " !" << endl;
262  }
263  }
264  // if(fVerbose>7)
265  // cout<<" "<<endl;//next hit content
266  }
267  // if(fMCnegative){
268  // if(fVerbose>7) cout<<"bad event, skip it!"<<endl;
269  // }
270  // if(fMCnegative) continue;
272 
274  // TNtuple *ntupTrkFit = new TNtuple("ntupTrkFit","Info about reconstructed trks: trk-cand vs. trk after fit","xtc:ytc:ztc:thtc:phitc:xtf:ytf:ztf:thtf:phitf:chi2");
275  for (Int_t iN = 0; iN < nGeaneTrks; iN++) { // loop over all reconstructed trks
276  //dummy values
277  glXrecLMD = -9999;
278  glYrecLMD = -9999;
279  glZrecLMD = -9999;
280  glThetarecLMD = -9999;
281  glPhirecLMD = -9999;
282  glXrec = -9999;
283  glYrec = -9999;
284  glZrec = -9999;
285  glThetarec = -9999;
286  glPhirec = -9999;
287  glMomrec = -9999;
288  glXmc = -9999;
289  glYmc = -9999;
290  glZmc = -9999;
291  glThetamc = -9999;
292  glPhimc = -9999;
293  glMommc = -9999;
294  glXmcLMD = -9999;
295  glYmcLMD = -9999;
296  glZmcLMD = -9999;
297  glThetamcLMD = -9999;
298  glPhimcLMD = -9999;
299  glMommcLMD = -9999;
300  trkRECStatus = -9999;
301  trkMCStatus = -9999;
302  glchi2 = -9999;
303  glPDG = -9999;
304  glNumMChits = -9999;
305  glNumDoubleMChits = -9999;
306 
307  PndTrack *trkpnd = (PndTrack*) fRecTracks->At(iN);
308  //double chi2 = trkpnd->GetChi2(); //[R.K. 01/2017] unused variable
309  FairTrackParP fFittedTrkP = trkpnd->GetParamFirst();
310  TVector3 PosRecLMD(fFittedTrkP.GetX(), fFittedTrkP.GetY(),
311  fFittedTrkP.GetZ());
312  TVector3 MomRecLMD(fFittedTrkP.GetPx(), fFittedTrkP.GetPy(),
313  fFittedTrkP.GetPz());
314  MomRecLMD *= fPbeam / MomRecLMD.Mag();
315 
316  int candID = trkpnd->GetRefIndex();
317  // cout<<"candID = "<<candID<<endl;
318  PndTrackCand *trkcand = (PndTrackCand*) fRecCandTracks->At(candID);
320  const int numPts = trkcand->GetNHits(); //read how many points in this track
321  PndTrackCandHit theHit1 = trkcand->GetSortedHit(0); //get 1st hit
322  Int_t index1 = theHit1.GetHitId();
323  PndSdsHit* Hit1sds = (PndSdsHit*) fRecHits->At(index1);
324  TVector3 posSeed = Hit1sds->GetPosition();
325  PndTrackCandHit theHit2 = trkcand->GetSortedHit(numPts - 1); //get last hit
326  Int_t index2 = theHit2.GetHitId();
327  PndSdsHit* Hit2sds = (PndSdsHit*) fRecHits->At(index2);
328  TVector3 pos2 = Hit2sds->GetPosition();
329  TVector3 MomRecCandLMD = pos2 - posSeed;
330  MomRecCandLMD *= fPbeam / MomRecCandLMD.Mag();
331  // TVector3 PosRecCandLMD = trkcand->getPosSeed();
332  // TVector3 MomRecCandLMD = trkcand->getDirSeed();
333  // MomRecCandLMD *=fPbeam/MomRecCandLMD.Mag();
334  }
336 
338  //Int_t nRecGEANEtrk = 0; //[R.K. 01/2017] unused variable
339  int MCtrk[nParticles]; //Number of participation this MCid in rec.tracks
340  int RECtrkMCid[nGeaneTrks]; //Assignment MC id to REC trk;
341  for (int nk = 0; nk < nParticles; nk++)
342  MCtrk[nk] = 0;
343  bool goodTrk[nGeaneTrks];
344  bool ghostTrk[nGeaneTrks];
345  for (Int_t iN = 0; iN < nGeaneTrks; iN++) {
346  goodTrk[iN] = false;
347  ghostTrk[iN] = false;
348  RECtrkMCid[iN] = -1;
349  }
350 
351  //int goodRectrk=0;//for missed trk-search //[R.K. 01/2017] unused variable
352  for (Int_t iN = 0; iN < nGeaneTrks; iN++) { // loop over all reconstructed trks
353  FairTrackParH *fRes = (FairTrackParH*) fRecBPTracks->At(iN);
354  TVector3 PosRec = fRes->GetPosition();
355  Double_t lyambda = fRes->GetLambda();
356  if (lyambda == 0) {
357  cout << "GEANE didn't propagate trk No." << iN << "!" << endl;
358  }
359 
360  PndTrack *trkpnd = (PndTrack*) fRecTracks->At(iN);
361  //double chi2 = trkpnd->GetChi2(); //[R.K. 01/2017] unused variable
362  int candID = trkpnd->GetRefIndex();
363  //if(fVerbose>5) cout<<"candID = "<<candID<<endl;
364  PndTrackCand *trkcand = (PndTrackCand*) fRecCandTracks->At(candID);
365  const int Ntrkcandhits = trkcand->GetNHits();
366  // PndSdsMCPoint* MCPointHit;
367  int MCtrkID[Ntrkcandhits]; //arrray of hits MCid
368  for (Int_t iHit = 0; iHit < Ntrkcandhits; iHit++) {
369  MCtrkID[iHit] = 9999;
370  }
371  Int_t diffIDs = 1;
372 
374  //bool emergExit=false; //[R.K. 01/2017] unused variable
375  if (fVerbose > 7)
376  cout << " *** REChits in trk (with " << Ntrkcandhits
377  << " hits) are made from MChits: " << endl; //start MC hit content
378  for (Int_t iHit = 0; iHit < Ntrkcandhits; iHit++) { // loop over rec.hits
379  if (fVerbose > 7)
380  cout << "No." << iHit << ": ";
381  PndTrackCandHit candhit = (PndTrackCandHit) (trkcand->GetSortedHit(iHit));
382  Int_t hitID = candhit.GetHitId();
383  PndSdsMergedHit* myHit = (PndSdsMergedHit*) (fRecHits->At(hitID));
384  int mcrefbot = myHit->GetSecondMCHit();
385  //bool badpxbot=false; //[R.K. 01/2017] unused variable
386  //bool badpxtop=false; //[R.K. 01/2017] unused variable
387  if (mcrefbot >= 0) {
388  PndSdsMCPoint* MCPointBot = (PndSdsMCPoint*) (fMCHits->At(mcrefbot));
389  int MCtrkid = MCPointBot->GetTrackID();
390  if (MCtrkid < 0) {
391  //emergExit=true; //[R.K. 01/2017] unused variable
392  if (fVerbose > 7)
393  cout << " " << MCtrkid << "!!!";
394  break; //TODO: how it is possible???
395  }
396  // if(MCtrkid>-1)
397  MCtrkID[iHit] = MCtrkid;
398  if (fVerbose > 7)
399  cout << " " << MCtrkid;
400  } else {
401  //badpxbot=true; //[R.K. 01/2017] unused variable
402  if (fVerbose > 7)
403  cout << " Ooops, mcrefbot = " << mcrefbot;
404  }
405  int mcreftop = myHit->GetRefIndex();
406 
407  if (mcreftop >= 0) {
408  // if(fVerbose>5)
409  // cout<<", "<<mcreftop<<")";
410  PndSdsMCPoint* MCPointTop = (PndSdsMCPoint*) (fMCHits->At(mcreftop));
411  int MCtrkid = MCPointTop->GetTrackID();
412  if (MCtrkid < 0) {
413  //emergExit=true; //[R.K. 01/2017] unused variable
414  if (fVerbose > 7)
415  cout << " " << MCtrkid << "!!!";
416  break; //TODO: how it is possible???
417  }
418  // if(MCtrkid>-1)
419  MCtrkID[iHit] = MCtrkid;
420  if (fVerbose > 7)
421  cout << " " << MCtrkid;
422  } else {
423  //badpxtop=true; //[R.K. 01/2017] unused variable
424  if (fVerbose > 7)
425  cout << " Ooops, mcreftop = " << mcreftop;
426  }
427 
428  if (fVerbose > 7)
429  cout << " " << endl; //next hit content
430  }
431  // if(emergExit) continue;
432 
433  // Sorting MC IDs ----------------------------------------
434  //TODO: sort it by c++ function
435  Int_t k, x;
436  //bool ch=false; //Was element changed? //[R.K. 01/2017] unused variable
437  Int_t nch = 0; //How many times?
438  for (Int_t n = 0; n < Ntrkcandhits; n++) { // n - current position
439  k = n;
440  x = MCtrkID[n];
441  for (Int_t m = n + 1; m < Ntrkcandhits; m++) // find the least element
442  if (MCtrkID[m] < x) {
443  k = m;
444  x = MCtrkID[m]; // k - index for the least element
445  //ch=true; //[R.K. 01/2017] unused variable
446  nch++;
447  }
448  MCtrkID[k] = MCtrkID[n];
449  MCtrkID[n] = x; // change position between the least and current elements
450  }
452 
454  Int_t prevID = MCtrkID[0];
455  for (int nk = 0; nk < nParticles; nk++) {
456  if (MCtrkID[0] == nk)
457  MCtrk[nk] = MCtrk[nk] + 1;
458  }
459  diffIDs = 1;
460  for (Int_t n = 1; n < Ntrkcandhits; n++) {
461  if (MCtrkID[n] > 9998)
462  break;
463  if (prevID < MCtrkID[n]) {
464  diffIDs++;
465  prevID = MCtrkID[n];
466  for (int nk = 0; nk < nParticles; nk++) {
467  if (MCtrkID[n] == nk)
468  MCtrk[nk] = MCtrk[nk] + 1;
469  }
470  }
471  }
472 
474 
476  if (diffIDs < 2) {
477  RECtrkMCid[iN] = MCtrkID[0];
478  // cout<<"MCtrkID[0] = "<<MCtrkID[0]<<endl;
479  if (MCtrkID[0] < 9999)
480  goodTrk[iN] = true;
481  else {
482  goodTrk[iN] = false;
483  ghostTrk[iN] = true;
484  }
485  } else {
486  vector<int> countMC_IDs(diffIDs);
487  for (int kn = 0; kn < diffIDs; kn++) {
488  countMC_IDs[kn] = 0;
489  }
490 
491  prevID = MCtrkID[0];
492  int diffCount = 0;
493  for (Int_t n = 0; n < Ntrkcandhits; n++) {
494  // cout<<"MCtrkID[n] = "<<MCtrkID[n]<<endl;
495  if (MCtrkID[n] > 9998) { //Trk from noise hits
496  goodTrk[iN] = false;
497  ghostTrk[iN] = true;
498  // break;
499  }
500 
501  if (prevID < MCtrkID[n]) {
502  diffCount++;
503  prevID = MCtrkID[n];
504  }
505  countMC_IDs[diffCount]++;
506  }
507  // if(fVerbose>1) cout<<""<<endl;
508  int maxID = countMC_IDs[0];
509  int posIDmax = 0;
510  for (int kn = 0; kn < diffIDs; kn++) {
511  if (fVerbose > 3)
512  cout << "countMC_IDs[" << kn << "]=" << countMC_IDs[kn]
513  << " 0.65*Ntrkcandhits = " << 0.65 * Ntrkcandhits << endl;
514  if (countMC_IDs[kn] > 0.65 * Ntrkcandhits) { //more then 65% of hits come from the same MC id
515  goodTrk[iN] = true;
516  ghostTrk[iN] = false;
517  } else {
518  if (!goodTrk[iN]) {
519  ghostTrk[iN] = true;
520  // cout<<" aaaaa!"<<iN<< " is ghost trk !aaaa"<<endl;
521  }
522  }
523 
524  if (countMC_IDs[kn] > maxID) {
525  maxID = countMC_IDs[kn];
526  posIDmax = kn;
527  }
528  }
529  prevID = MCtrkID[0];
530  diffCount = 0;
531  for (Int_t n = 0; n < Ntrkcandhits; n++) {
532  if (diffCount == posIDmax)
533  RECtrkMCid[iN] = prevID;
534  if (prevID < MCtrkID[n]) {
535  diffCount++;
536  prevID = MCtrkID[n];
537  }
538  }
539  }
541 
542  // /// Comporision MC, trk-candidates and REConstructed tracks --------------------------------------------
543 
544  // if(diffIDs>-1){ // All tracks
545 
546  // /// Read MC track parameters ----------------------------------------------------------
547  // int MCidforREC = RECtrkMCid[iN];
548  // PndMCTrack *mctrk =(PndMCTrack*) fMCTracks->At(MCidforREC);
549  // Int_t mcID = mctrk->GetPdgCode();
550  // // glPDG = mctrk->GetPdgCode();
551  // TVector3 MomMC = mctrk->GetMomentum();
552  // Double_t thetaMC = MomMC.Theta();
553  // Double_t phiMC = MomMC.Phi();
554  // ///------------------------------------------------------------------------------------
555 
556  // /// Read track-parameters after back-propagation ---------------------------------------
557  // //TODO: problem with covarance matrix in FairTrackParH???
558  // Double_t CovGEANELAB[6][6];
559  // fRes->GetMARSCov(CovGEANELAB);
560  // Double_t errX = fRes->GetDX();
561  // Double_t errY = fRes->GetDY();
562  // Double_t errZ = fRes->GetDZ();
563 
564  // Double_t errPx = fRes->GetDPx();
565  // Double_t errPy = fRes->GetDPy();
566  // Double_t errPz = fRes->GetDPz();
567  // TVector3 errMomBP(errPx,errPy,errPz);
568 
569  // Double_t thetaBP = TMath::Pi()/2. - lyambda;
570  // Double_t err_lyambda = fRes->GetDLambda();
571  // if(err_lyambda==0) err_lyambda = errMomBP.Theta();
572  // // Double_t err_lyambda = errMom.Theta();
573  // Double_t phiBP = fRes->GetPhi();
574  // Double_t err_phi = fRes->GetDPhi();
575  // if(err_phi==0) err_phi=errMomBP.Phi();
576  // Double_t errMomRecBP = fRes->GetDQp();
577  // TVector3 MomRecBP = fRes->GetMomentum();
578 
579  // ///==================================
580  // }
581  // ///-------------------------------------------------------------------------------------
582  }
583 
587  //NB: it's only quantity and _ NOT_ quality check. 1st REC trk will be called GOOD and 2nd and others GHOST
588  int RECassigMC[nGeaneTrks]; //how many time MCid from this RECtrk met in other RECtrks
589  for (Int_t iN = 0; iN < nGeaneTrks; iN++)
590  RECassigMC[iN] = 1;
591 
592  for (Int_t iN = 0; iN < nGeaneTrks; iN++) {
593  int mc_i = RECtrkMCid[iN];
594  for (Int_t jN = (iN + 1); jN < nGeaneTrks; jN++) {
595  int mc_j = RECtrkMCid[jN];
596  if (mc_i == mc_j) {
597  RECassigMC[iN]++;
598  if (fVerbose > 3)
599  cout << "GHOST?!: rec.trks #" << iN << " and #" << jN
600  << "were assigned to one MC trk =|" << endl;
601  }
602  }
603  }
605 
606  int goodRecII = 0, ghostRecI = 0, ghostRecII = 0;
607  for (Int_t iN = 0; iN < nGeaneTrks; iN++) { //ghost/good trk determination
608 
609  int trkType = -1;
610  if (RECassigMC[iN] > 1) { //one RECtrk was already called GOOD for this MCid -> all next are GHOSTs
611  goodTrk[iN] = false;
612  ghostTrk[iN] = true;
613  trkType = +2;
614  ghostRecI++;
615  }
616 
617  if (goodTrk[iN]) {
618  if (fVerbose > 7)
619  cout << "... RECtrk#" << iN << " was defined as GOOD ..." << endl;
620  goodRecII++;
621  trkType = 0;
622  }
623  if (ghostTrk[iN]) {
624  if (fVerbose > 7)
625  cout << "... RECtrk#" << iN << " was defined as GHOST ..." << endl;
626  if (trkType < 0) {
627  trkType = +1;
628  ghostRecII++;
629  }
630  }
631  if (trkType < 0) {
632  cout << "Ooops, RECtrk isn't GOOD and isn't GHOST!!!" << endl;
633  trkType = -10;
634  }
635  // if(trkType<0) continue; //GEANE didn't propagate this track
636  FairTrackParH *fRes = (FairTrackParH*) fRecBPTracks->At(iN);
637  Double_t lyambda = fRes->GetLambda();
638  if (lyambda == 0)
639  trkType = -100; //Trk was not back-propagated!
640  // cout<<"GEANE didn't propagate "<<iN<<" trk!"<<endl;
641  Double_t thetaBP = TMath::Pi() / 2. - lyambda;
642  Double_t phiBP = fRes->GetPhi();
643  TVector3 MomRecBP = fRes->GetMomentum();
644  TVector3 PosBP = fRes->GetPosition();
645 
646  //Fill tree with rec vs. mc trk info ------------------------------------------------
647  glXrec = PosBP.X();
648  glYrec = PosBP.Y();
649  glZrec = PosBP.Z();
650  glThetarec = thetaBP;
651  glPhirec = phiBP;
652  glMomrec = MomRecBP.Mag();
653  Double_t errX = fRes->GetDX();
654  Double_t errY = fRes->GetDY();
655  Double_t errZ = fRes->GetDZ();
656 
657  Double_t errPx = fRes->GetDPx();
658  Double_t errPy = fRes->GetDPy();
659  Double_t errPz = fRes->GetDPz();
660  TVector3 errMomBP(errPx, errPy, errPz);
661 
662  Double_t err_lyambda = fRes->GetDLambda();
663  // if(err_lyambda==0) err_lyambda = errMomBP.Theta();
664  Double_t err_phi = fRes->GetDPhi();
665 
666  PndTrack *trkpnd = (PndTrack*) fRecTracks->At(iN);
667  glchi2 = trkpnd->GetChi2();
668  FairTrackParP fFittedTrkP = trkpnd->GetParamFirst();
669  TVector3 PosRecLMD(fFittedTrkP.GetX(), fFittedTrkP.GetY(),
670  fFittedTrkP.GetZ());
671  TVector3 MomRecLMD(fFittedTrkP.GetPx(), fFittedTrkP.GetPy(),
672  fFittedTrkP.GetPz());
673  MomRecLMD *= fPbeam / MomRecLMD.Mag();
674  glXrecLMD = PosRecLMD.X();
675  glYrecLMD = PosRecLMD.Y();
676  glZrecLMD = PosRecLMD.Z();
677  glThetarecLMD = MomRecLMD.Theta();
678  glPhirecLMD = MomRecLMD.Phi();
679  trkRECStatus = trkType;
680  int candID = trkpnd->GetRefIndex();
681  PndTrackCand *trkcand = (PndTrackCand*) fRecCandTracks->At(candID);
682  PndTrackCandHit candhit = (PndTrackCandHit) (trkcand->GetSortedHit(0)); //1st hit info
683  Int_t hitID = candhit.GetHitId();
684  PndSdsMergedHit* myHit = (PndSdsMergedHit*) (fRecHits->At(hitID));
685  int sensorID = myHit->GetSensorID();
686 
687  auto& lmd_geo_helper(PndLmdGeometryHelper::getInstance());
688  auto const& digi_info(lmd_geo_helper.getHitLocationInfo(sensorID));
689  glModule = digi_info.module;
690  glHalf = digi_info.detector_half;
691  glTrkTime = myHit->GetTimeStamp();
692  //TODO: include info from GetTimeStampError()
693 
694  // if(trkType>0){ //TODO: ghost-doubled trks has MC trk!!!
695  // // int MCidforREC = RECtrkMCid[iN];
696 
697  // // glXmc= -9999; glYmc =-9999; glZmc = -9999; glThetamc =-9999; glPhimc = -9999; glMommc = -9999;
698  // // glXmcLMD = -9999; glYmcLMD =-9999; glZmcLMD = -9999; glThetamcLMD =-9999; glPhimcLMD = -9999; glMommcLMD = -9999;
699  // // glNumMChits = -9999; glNumDoubleMChits = -9999; // glEvTime = -9999;
700  // // trkMCStatus = -9999;
701  // // glPDG = -9999;
702  // }
703  // else{
704 
705  // }
706  int MCidforREC = RECtrkMCid[iN];
707  // if(!goodTrk[iN]) continue;
708  //cout<<" MCidforREC = "<<MCidforREC<<endl;
709  if (MCidforREC > 9998) {
710  glXmc = -9999;
711  glYmc = -9999;
712  glZmc = -9999;
713  glThetamc = -9999;
714  glPhimc = -9999;
715  glMommc = -9999;
716  glXmcLMD = -9999;
717  glYmcLMD = -9999;
718  glZmcLMD = -9999;
719  glThetamcLMD = -9999;
720  glPhimcLMD = -9999;
721  glMommcLMD = -9999;
722  glNumMChits = -9999;
723  glNumDoubleMChits = -9999; // glEvTime = -9999;
724  trkMCStatus = -9999;
725  glPDG = -9999;
726  } else {
727  PndMCTrack *mctrk = (PndMCTrack*) fMCTracks->At(MCidforREC);
728  glPDG = mctrk->GetPdgCode();
729  // glEvTime = mctrk->GetStartTime(); //this gives start time of trk, doesn't matter since all el. pbar start in IP
730  TVector3 MomMC = mctrk->GetMomentum();
731  glThetamc = MomMC.Theta();
732  glPhimc = MomMC.Phi();
733  glMommc = MomMC.Mag();
734  TVector3 StartMC = mctrk->GetStartVertex();
735  glXmc = StartMC.X();
736  glYmc = StartMC.Y();
737  glZmc = StartMC.Z();
738  int movID = mctrk->GetMotherID();
739  trkMCStatus = movID;
740  // if(movID<0) trkMCStatus=1;
741  // else
742  // trkMCStatus=0;
743  glNumMChits = MCtksSIMhits[MCidforREC];
744  glNumDoubleMChits = MCDoubleHits[MCidforREC];
745 
746  // //Get MC info in LMD
747 
748  int mcrefbot = myHit->GetSecondMCHit();
749  int mcreftop = myHit->GetRefIndex();
750  if (fVerbose > 0) {
751  cout << "mcrefbot = " << mcrefbot << " mcreftop = " << mcreftop << endl;
752  }
753  PndSdsMCPoint* MCPointHit;
754  if (mcreftop >= 0) {
755  MCPointHit = (PndSdsMCPoint*) (fMCHits->At(mcreftop));
756  } else {
757  if (mcrefbot >= 0)
758  MCPointHit = (PndSdsMCPoint*) (fMCHits->At(mcrefbot));
759  else
760  MCPointHit = NULL;
761  }
762  if (MCPointHit != NULL) {
763  TVector3 PosMClmd = MCPointHit->GetPosition();
764 
765  double pxTrue = MCPointHit->GetPx();
766  double pyTrue = MCPointHit->GetPy();
767  double pzTrue = MCPointHit->GetPz();
768  TVector3 MomMClmd(pxTrue, pyTrue, pzTrue);
769  glXmcLMD = PosMClmd.X();
770  glYmcLMD = PosMClmd.Y();
771  glZmcLMD = PosMClmd.Z();
772  glThetamcLMD = MomMClmd.Theta();
773  glPhimcLMD = MomMClmd.Phi();
774  glMommcLMD = MomMClmd.Mag();
775  } else { // Noise MC hit
776  // glXmcLMD = -9999; glYmcLMD = -9999; glZmcLMD = -9999;
777  // glThetamcLMD = -9999; glPhimcLMD = -9999; glMommcLMD = -9999;
778  glXmc = -9999;
779  glYmc = -9999;
780  glZmc = -9999;
781  glThetamc = -9999;
782  glPhimc = -9999;
783  glMommc = -9999;
784  glXmcLMD = -9999;
785  glYmcLMD = -9999;
786  glZmcLMD = -9999;
787  glThetamcLMD = -9999;
788  glPhimcLMD = -9999;
789  glMommcLMD = -9999;
790  glNumMChits = -9999;
791  glNumDoubleMChits = -9999; // glEvTime = -9999;
792  trkMCStatus = -99; //1st rec hit = Noise MC hit TODO: check if all hits of this rec.trk are noise ones
793  glPDG = -99; //1st rec hit = Noise MC hit
794  }
795  }
796  // tRECMCtrks->Fill();
797  PndLmdTrackQ *trkqlmd =
798  new ((*fTrackQ)[fTrackQ->GetEntriesFast()]) PndLmdTrackQ(glMommc);
799 
800  trkqlmd->SetTrkRecStatus(trkRECStatus);
801  trkqlmd->SetSecondary(trkMCStatus);
802  trkqlmd->SetPDGcode(glPDG);
803  trkqlmd->SetLMDpoint(glXrecLMD, glYrecLMD, glZrecLMD);
804  trkqlmd->SetLMDdir(glThetarecLMD, glPhirecLMD);
805  trkqlmd->SetLMDchi2(glchi2);
806  trkqlmd->SetIPpoint(glXrec, glYrec, glZrec);
807  trkqlmd->SetIPmom(glThetarec, glPhirec, glMomrec);
808  trkqlmd->SetIPerrpoint(errX, errY, errZ);
809  trkqlmd->SetIPerrmom(err_lyambda, err_phi, errMomBP.Mag());
810 
811  trkqlmd->SetMCpoint(glXmc, glYmc, glZmc);
812  trkqlmd->SetMCmom(glThetamc, glPhimc, glMommc);
813  // trkqlmd->SetSecondary(trkMCStatus);
814  trkqlmd->SetNumMChits(glNumMChits);
815  trkqlmd->SetNumDoubleMChits(glNumDoubleMChits);
816  trkqlmd->SetMCpointLMD(glXmcLMD, glYmcLMD, glZmcLMD);
817  trkqlmd->SetMCmomLMD(glThetamcLMD, glPhimcLMD, glMommcLMD);
818  // trkqlmd->SetTotEvCharge(TotCharge);
819  trkqlmd->SetSumEvPDG(sumID);
820  trkqlmd->SetEvMCMulti(nParticles);
821  trkqlmd->SetEvRECMulti(nGeaneTrks);
822  trkqlmd->SetEvTime(glEvTime);
823  trkqlmd->SetTrkTime(glTrkTime);
824  trkqlmd->SetModule(glModule);
825  trkqlmd->SetHalf(glHalf);
826  //(end) Fill tree with rec vs. mc trk info ---------------------------------------
827 
828  }
830 
832  if ((nParticles - goodRecII) > 0) { //missed trks
833 
834  int nMCmissedTrkSearch = 0;
835  int nMCmissedLossHits = 0;
836  for (int imc = 0; imc < nParticles; imc++) { //MC trks
837  bool missTrk = true;
838  for (int irec = 0; irec < nGeaneTrks; irec++) { //RECids assigment
839  int mc_comp = RECtrkMCid[irec];
840  if (mc_comp == imc && goodTrk[irec])
841  missTrk = false;
842  }
843  PndMCTrack *mctrk = (PndMCTrack*) fMCTracks->At(imc);
844  int movID = mctrk->GetMotherID();
845  TVector3 MomMC = mctrk->GetMomentum();
846  TVector3 PosMC = mctrk->GetStartVertex();
847  glXmc = PosMC.X();
848  glYmc = PosMC.Y();
849  glZmc = PosMC.Z();
850  glThetamc = MomMC.Theta();
851  glPhimc = MomMC.Phi();
852  glMommc = MomMC.Mag();
853 
854  int trkQ = 0;
856  int minHits = MCDoubleHits[imc] + 2;
857  if (MCDoubleHits[imc] == 4)
858  minHits -= 1;
859 
860  if (missTrk) {
861  if (MCtksREChits[imc] > minHits) { //missed during trk-search
862  if (fVerbose > 7) {
863  cout << " --- MCtrk#" << imc
864  << " was defined as MISSED during trk-search (#MChits="
865  << MCtksREChits[imc] << " with limit>" << minHits << ")"
866  << endl;
867  if (MCDoubleHits[imc] > 0)
868  cout << " NB:this MCtrk contains " << MCDoubleHits[imc]
869  << " double hits!" << endl;
870  }
871  trkQ = -1;
872  nMCmissedTrkSearch++;
873  } else { //missed due to small amount of hits
874  if (MCtksREChits[imc] > 0) { // if MCtksREChits[imc]==0 trk is out of measurment range
875  // if(MCtksREChits[imc]>0){// if MCtksREChits[imc]==0 trk is out of measurment range
876  nMCmissedLossHits++;
877  trkQ = -2;
878  if (fVerbose > 7)
879  cout << " --- MCtrk#" << imc
880  << " was defined as MISSED due to little amount of hits (#MChits="
881  << MCtksREChits[imc] << " with limit>" << minHits << ")"
882  << endl;
883  int hitArr = -1;
884  // cout<<"nMCHits = "<<nMCHits<<endl;
885  for (int imhc = 0; imhc < nMCHits; imhc++) {//find corresponding MC hit for this MC trk
886  PndSdsMCPoint* MCPoint = (PndSdsMCPoint*) (fMCHits->At(imhc));
887  int MCtrkID = MCPoint->GetTrackID();
888  // cout<<"MCtrkID = "<<MCtrkID<<" imc = "<<imc<<endl;
889  if (MCtrkID == imc) {
890  hitArr = imhc;
891  break;
892  }
893  }
894  if (hitArr > 0) {
895  PndSdsMCPoint* MCPointHit = (PndSdsMCPoint*) (fMCHits->At(hitArr));
896  TVector3 PosMClmd = MCPointHit->GetPosition();
897  double pxTrue = MCPointHit->GetPx();
898  double pyTrue = MCPointHit->GetPy();
899  double pzTrue = MCPointHit->GetPz();
900  TVector3 MomMClmd(pxTrue, pyTrue, pzTrue);
901  glXmcLMD = PosMClmd.X();
902  glYmcLMD = PosMClmd.Y();
903  glZmcLMD = PosMClmd.Z();
904  glThetamcLMD = MomMClmd.Theta();
905  glPhimcLMD = MomMClmd.Phi();
906  glMommcLMD = MomMClmd.Mag();
907  } else {
908  trkQ = -3;
909  glXmcLMD = -9999;
910  glYmcLMD = -9999;
911  glZmcLMD = -9999;
912  glThetamcLMD = -9999;
913  glPhimcLMD = -9999;
914  glMommcLMD = -9999;
915  }
916  } else {
917  trkQ = -3;
918  glXmcLMD = -9999;
919  glYmcLMD = -9999;
920  glZmcLMD = -9999;
921  glThetamcLMD = -9999;
922  glPhimcLMD = -9999;
923  glMommcLMD = -9999;
924  }
925  }
926  // cout<<"trk was marked ad missed, reason trkQ="<<trkQ<<endl;
927 
928  //Fill tree with rec vs. mc trk info for missed trks ------------------------------------------------
929  glXrec = -9999;
930  glYrec = -9999;
931  glZrec = -9999;
932  glThetarec = -9999;
933  glPhirec = -9999;
934  glXrecLMD = -9999;
935  glYrecLMD = -9999;
936  glZrecLMD = -9999;
937  glThetarecLMD = -9999;
938  glPhirecLMD = -9999;
939  glMomrec = -9999;
940  glchi2 = -9999;
941  glTrkTime = -9999;
942  trkRECStatus = trkQ;
943 
944  trkMCStatus = movID;
945  glPDG = mctrk->GetPdgCode();
946  // glEvTime = mctrk->GetStartTime();
947  // if(movID<0) trkMCStatus=0;
948  // if(movID>=0) trkMCStatus=+1;
949  // tRECMCtrks->Fill();
950  if (!fWriteAllMC && trkRECStatus == -3) {
951  } else {
952  TClonesArray& clref = *fTrackQ;
953  Int_t size = clref.GetEntriesFast();
954  PndLmdTrackQ *trkqlmd = new (clref[size]) PndLmdTrackQ(glMommc);
955  trkqlmd->SetTrkRecStatus(trkRECStatus);
956  trkqlmd->SetSecondary(trkMCStatus);
957  trkqlmd->SetPDGcode(glPDG);
958  trkqlmd->SetLMDpoint(glXrecLMD, glYrecLMD, glZrecLMD);
959  trkqlmd->SetLMDdir(glThetarecLMD, glPhirecLMD);
960  trkqlmd->SetLMDchi2(glchi2);
961  trkqlmd->SetIPpoint(glXrec, glYrec, glZrec);
962  trkqlmd->SetIPmom(glThetarec, glPhirec, glMomrec);
963  trkqlmd->SetMCpoint(glXmc, glYmc, glZmc);
964  trkqlmd->SetMCmom(glThetamc, glPhimc, glMommc);
965  // trkqlmd->SetSecondary(trkMCStatus);
966  trkqlmd->SetIPerrpoint(-9999, -9999, -9999);
967  trkqlmd->SetIPerrmom(-9999, -9999, -9999);
968  trkqlmd->SetNumMChits(MCtksREChits[imc]);
969  trkqlmd->SetNumDoubleMChits(MCDoubleHits[imc]);
970  // cout<<"glXmcLMD = "<<glXmcLMD<<endl;
971  trkqlmd->SetMCpointLMD(glXmcLMD, glYmcLMD, glZmcLMD);
972  trkqlmd->SetMCmomLMD(glThetamcLMD, glPhimcLMD, glMommcLMD);
973  // cout<<"glPhimcLMD = "<<glPhimcLMD<<endl;
974  // trkqlmd->SetTotEvCharge(TotCharge);
975  trkqlmd->SetSumEvPDG(sumID);
976  trkqlmd->SetEvMCMulti(nParticles);
977  trkqlmd->SetEvRECMulti(nGeaneTrks);
978  trkqlmd->SetEvTime(glEvTime);
979  trkqlmd->SetTrkTime(glTrkTime);
980  trkqlmd->SetModule(glModule);
981  trkqlmd->SetHalf(glHalf);
982  }
983  //(end) Fill tree with rec vs. mc trk info for missed trks ---------------------------------------
984  }
985  }
986  }
987 
988  fMCHits->Delete();
989  fMCTracks->Delete();
990  fRecHits->Delete();
991  fRecCandTracks->Delete();
992  fRecTracks->Delete();
993  fRecBPTracks->Delete();
994 
995  if (fVerbose > 2)
996  cout << "PndLmdTrkQTask::Exec END!" << endl;
997  fEventNr++;
998 }
void SetModule(int mod)
Definition: PndLmdTrackQ.h:197
TClonesArray * fRecBPTracks
TClonesArray * fMCTracks
void SetEvRECMulti(int tot)
Definition: PndLmdTrackQ.h:185
int fVerbose
Definition: poormantracks.C:24
__m128 m
Definition: P4_F32vec4.h:28
TVector3 GetPosition() const
Definition: PndSdsHit.h:93
int nk
Definition: toy_core.C:124
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
void SetMCmom(double theta, double phi, double mom)
Definition: PndLmdTrackQ.h:120
void SetMCpoint(double x, double y, double z)
Definition: PndLmdTrackQ.h:114
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
void SetTrkRecStatus(int st)
Definition: PndLmdTrackQ.h:47
int n
static PndLmdGeometryHelper & getInstance()
void SetLMDdir(double theta, double phi)
Definition: PndLmdTrackQ.h:69
TFile * MCPoint
Definition: anaLmdDigi.C:25
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
void SetIPpoint(double x, double y, double z)
Definition: PndLmdTrackQ.h:84
void SetMCpointLMD(double x, double y, double z)
Definition: PndLmdTrackQ.h:135
TClonesArray * fTrackQ
void SetIPerrpoint(double errx, double erry, double errz)
Definition: PndLmdTrackQ.h:102
void SetMCmomLMD(double theta, double phi, double mom)
Definition: PndLmdTrackQ.h:138
void SetPDGcode(int pdg)
Definition: PndLmdTrackQ.h:56
Double_t
void SetSecondary(int sec)
Definition: PndLmdTrackQ.h:151
void SetIPerrmom(double errtheta, double errphi, double errmom)
Definition: PndLmdTrackQ.h:108
TClonesArray * fRecCandTracks
TClonesArray * fMCHits
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
void SetEvTime(double evtm)
Definition: PndLmdTrackQ.h:191
void SetTrkTime(double trktm)
Definition: PndLmdTrackQ.h:192
TClonesArray * fRecTracks
void SetSumEvPDG(int sumid)
Definition: PndLmdTrackQ.h:170
void SetLMDchi2(double chi2)
Definition: PndLmdTrackQ.h:78
TVector3 GetPosition() const
Definition: PndSdsMCPoint.h:90
void SetLMDpoint(double x, double y, double z)
Definition: PndLmdTrackQ.h:63
Int_t GetSecondMCHit() const
Double_t x
void SetNumMChits(int num)
Definition: PndLmdTrackQ.h:158
Double_t GetChi2() const
Definition: PndTrack.h:34
void SetNumDoubleMChits(int num)
Definition: PndLmdTrackQ.h:164
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
void SetEvMCMulti(int tot)
Definition: PndLmdTrackQ.h:179
Int_t GetRefIndex() const
Definition: PndTrack.h:36
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
Int_t GetSensorID() const
Definition: PndSdsHit.h:90
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
Int_t GetHitId() const
Double_t Pi
TClonesArray * fRecHits
void SetIPmom(double theta, double phi, double mom)
Definition: PndLmdTrackQ.h:90
FairTrackParP GetParamFirst()
Definition: PndTrack.h:49
void SetHalf(int hf)
Definition: PndLmdTrackQ.h:198
void PndLmdTrkQTask::Finish ( )
virtual

Definition at line 1000 of file PndLmdTrkQTask.cxx.

1000  {
1001 }
InitStatus PndLmdTrkQTask::Init ( )
virtual

Virtual method Init

Definition at line 65 of file PndLmdTrkQTask.cxx.

References fGeaneName, fMCHits, fMCTracks, fRecBPTracks, fRecCandTracks, fRecHits, fRecTracks, fTrackBranchName, and fTrackQ.

65  {
66 
67  // Get RootManager
68  FairRootManager* ioman = FairRootManager::Instance();
69  if (!ioman) {
70  std::cout << "-E- PndLmdTrkQTask::Init: " << "RootManager not instantiated!"
71  << std::endl;
72  return kFATAL;
73  }
74 
75  // fMCHeader = (TClonesArray*) ioman->GetObject("MCEventHeader.");
76  // if ( !fMCHeader) {
77  // std::cout << "-W- PndLmdTrkQTask::Init: "<< "No MCEventHeader "<<" object!" << std::endl;
78  // return kERROR;
79  // }
80 
81  //Get MC points
82  fMCHits = (TClonesArray*) ioman->GetObject("LMDPoint");
83  if (!fMCHits) {
84  std::cout << "-W- PndLmdTrkQTask::Init: " << "No LMDPoint" << " array!"
85  << std::endl;
86  return kERROR;
87  }
88 
89  //Get MC tracks
90  fMCTracks = (TClonesArray*) ioman->GetObject("MCTrack");
91  if (!fMCTracks) {
92  std::cout << "-W- PndLmdTrkQTask::Init: " << "No MCTrack" << " array!"
93  << std::endl;
94  return kERROR;
95  }
96 
97  //Get rec. hits
98  fRecHits = (TClonesArray*) ioman->GetObject("LMDHitsMerged");
99  if (!fRecHits) {
100  std::cout << "-W- PndLmdTrkQTask::Init: " << "No LMDHitsMerged" << " array!"
101  << std::endl;
102  return kERROR;
103  }
104  //Get trk cand
105  fRecCandTracks = (TClonesArray*) ioman->GetObject("LMDTrackCand");
106  if (!fRecCandTracks) {
107  std::cout << "-W- PndLmdTrkQTask::Init: " << "No LMDTrackCand" << " array!"
108  << std::endl;
109  return kERROR;
110  }
111  //Get rec.tracks before back propagation
112  fRecTracks = (TClonesArray*) ioman->GetObject(fTrackBranchName);
113  if (!fRecTracks) {
114  std::cout << "-W- PndLmdTrkQTask::Init: " << "No LMDPndTrackFilt"
115  << " array!" << std::endl;
116  return kERROR;
117  }
118 
119  //Get rec.tracks after back propagation
120  // fRecBPTracks = (TClonesArray*) ioman->GetObject("GeaneTrackFinal");//use raw reconstructed data
121  // fRecBPTracks = (TClonesArray*) ioman->GetObject("LMDCleanTrack");//use cleaned data
122  fRecBPTracks = (TClonesArray*) ioman->GetObject(fGeaneName);
123  if (!fRecBPTracks) {
124  std::cout << "-W- PndLmdTrkQTask::Init: " << "No GeaneTrackFinal"
125  << " array!" << std::endl;
126  return kERROR;
127  }
128 
129  fTrackQ = new TClonesArray("PndLmdTrackQ");
130  ioman->Register("LMDTrackQ", "TrkQ", fTrackQ, kTRUE);
131 
132  // fDetName = new TClonesArray("TObjString");
133  // ioman->Register("DetName", "TrkQ", fDetName, kTRUE);
134 
135  // fGeoH = PndGeoHandling::Instance();
136  //FairRun* fRun = FairRun::Instance(); //[R.K. 01/2017] unused variable?
137  //FairRuntimeDb* rtdb = fRun->GetRuntimeDb(); //[R.K. 01/2017] unused variable
138  // TDatabasePDG *fdbPDG = TDatabasePDG::Instance();
139  // fdbPDG = TDatabasePDG::Instance();
140 
141  return kSUCCESS;
142 }
TClonesArray * fRecBPTracks
TClonesArray * fMCTracks
TClonesArray * fTrackQ
TString fTrackBranchName
TClonesArray * fRecCandTracks
TClonesArray * fMCHits
TClonesArray * fRecTracks
TClonesArray * fRecHits
TString fGeaneName
void PndLmdTrkQTask::SetWriteMC ( bool  wr)
inline

Definition at line 57 of file PndLmdTrkQTask.h.

References fWriteAllMC.

Referenced by runLumiPixel7TrksQA().

57 {fWriteAllMC = wr;}

Member Data Documentation

int PndLmdTrkQTask::fEventNr
private

Definition at line 72 of file PndLmdTrkQTask.h.

Referenced by Exec().

TString PndLmdTrkQTask::fGeaneName
private

Definition at line 61 of file PndLmdTrkQTask.h.

Referenced by Init(), and PndLmdTrkQTask().

TClonesArray* PndLmdTrkQTask::fMCHits
private

Definition at line 63 of file PndLmdTrkQTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndLmdTrkQTask::fMCTracks
private

Definition at line 64 of file PndLmdTrkQTask.h.

Referenced by Exec(), and Init().

Double_t PndLmdTrkQTask::fPbeam
private

Definition at line 74 of file PndLmdTrkQTask.h.

Referenced by Exec(), and PndLmdTrkQTask().

TClonesArray* PndLmdTrkQTask::fRecBPTracks
private

Definition at line 68 of file PndLmdTrkQTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndLmdTrkQTask::fRecCandTracks
private

Definition at line 66 of file PndLmdTrkQTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndLmdTrkQTask::fRecHits
private

Definition at line 65 of file PndLmdTrkQTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndLmdTrkQTask::fRecTracks
private

Definition at line 67 of file PndLmdTrkQTask.h.

Referenced by Exec(), and Init().

TString PndLmdTrkQTask::fTrackBranchName
private

Definition at line 62 of file PndLmdTrkQTask.h.

Referenced by Init(), and PndLmdTrkQTask().

TClonesArray* PndLmdTrkQTask::fTrackQ
private

Definition at line 70 of file PndLmdTrkQTask.h.

Referenced by Exec(), and Init().

bool PndLmdTrkQTask::fWriteAllMC
private

Definition at line 77 of file PndLmdTrkQTask.h.

Referenced by Exec(), PndLmdTrkQTask(), and SetWriteMC().


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