]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliCFVertexingHFLctoV0bachelor.cxx
c7c63b07a0f7f73e62ec872183860faad3c4b5cc
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / AliCFVertexingHFLctoV0bachelor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2011, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */ 
17
18 //-----------------------------------------------------------------------
19 // Class for HF corrections as a function of many variables and steps
20 // For Lc->V0+bachelor
21 // 
22 // - Based on AliCFVertexingHFCascade -
23 //
24 // Contact : A.De Caro - decaro@sa.infn.it
25 //           Centro 'E.Fermi' - Rome (Italy)
26 //           INFN and University of Salerno (Italy)
27 //
28 //-----------------------------------------------------------------------
29
30 #include "TDatabasePDG.h"
31 #include "TClonesArray.h"
32 #include "AliAODv0.h"
33 #include "AliAODMCParticle.h"
34 #include "AliAODRecoDecayHF.h"
35 #include "AliAODRecoCascadeHF.h"
36 #include "AliCFTaskVertexingHF.h"
37 #include "AliCFContainer.h"
38 #include "AliCFVertexingHF.h"
39 #include "AliCFVertexingHFLctoV0bachelor.h"
40
41 #include <Riostream.h>
42
43 using std::cout;
44 using std::endl;
45
46 ClassImp(AliCFVertexingHFLctoV0bachelor)
47
48 //_________________________________________
49 AliCFVertexingHFLctoV0bachelor::AliCFVertexingHFLctoV0bachelor():
50 fGenLcOption(0)
51 {
52   // standard constructor
53 }
54
55 //_____________________________________
56 AliCFVertexingHFLctoV0bachelor::AliCFVertexingHFLctoV0bachelor(TClonesArray *mcArray, UShort_t originDselection, Int_t lcDecay):
57 AliCFVertexingHF(mcArray, originDselection),
58   fGenLcOption(lcDecay)
59 {
60   // standard constructor
61
62   SetNProngs(3);
63   fPtAccCut=new Float_t[fProngs];
64   fEtaAccCut=new Float_t[fProngs];
65   for(Int_t iP=0; iP<fProngs; iP++){
66     fPtAccCut[iP]=0.1;
67     fEtaAccCut[iP]=0.9;
68   }
69
70 }
71
72
73 //_____________________________________
74 AliCFVertexingHFLctoV0bachelor& AliCFVertexingHFLctoV0bachelor::operator=(const AliCFVertexingHFLctoV0bachelor& c)
75 {
76   // operator =
77  
78   if  (this != &c) {
79     AliCFVertexingHF::operator=(c);
80   }
81
82   return *this;
83
84 }
85
86 //__________________________________________
87 Bool_t AliCFVertexingHFLctoV0bachelor::SetRecoCandidateParam(AliAODRecoDecayHF *recoCand)
88 {
89   // set the AliAODRecoDecay candidate
90   
91   Bool_t bSignAssoc = kFALSE;
92  
93   fRecoCandidate = recoCand;
94   if (!fRecoCandidate) {
95     AliError("fRecoCandidate not found, problem in assignement\n");
96     return bSignAssoc;
97   }
98   
99   if (fRecoCandidate->GetPrimaryVtx()) AliDebug(4,"fReco Candidate has a pointer to PrimVtx\n");
100   
101   AliAODRecoCascadeHF* lcV0bachelor = dynamic_cast<AliAODRecoCascadeHF*>(fRecoCandidate);
102   if (!lcV0bachelor) {
103     AliDebug(2,"cascade doesn't exist, skipped!");
104     return bSignAssoc;
105   }
106
107   if (!lcV0bachelor->GetSecondaryVtx()) {
108     AliDebug(2,"No secondary vertex for cascade");
109     return bSignAssoc;
110   }
111
112   if (lcV0bachelor->GetNDaughters()!=2) {
113     AliDebug(2,Form("No 2 daughters for current cascade (nDaughters=%d)",lcV0bachelor->GetNDaughters()));
114     return bSignAssoc;
115   }
116
117   AliVTrack *cascTrk0 = dynamic_cast<AliVTrack*>(lcV0bachelor->GetDaughter(0));
118   AliVTrack *cascTrk1 = dynamic_cast<AliVTrack*>(lcV0bachelor->GetDaughter(1));
119   if (!cascTrk0 || !cascTrk1) {
120     AliDebug(2,"At least one of V0daughters doesn't exist");
121     return bSignAssoc;
122   }
123
124   AliAODv0 * v0part = dynamic_cast<AliAODv0*>(lcV0bachelor->Getv0());
125   AliAODTrack * bachPart = dynamic_cast<AliAODTrack*>(lcV0bachelor->GetBachelor());
126   if (!v0part || !bachPart) {
127     AliDebug(2,"No V0 or no bachelor for current cascade");
128     return bSignAssoc;
129   }
130
131   if (bachPart->GetID()<0) {
132     AliDebug(2,Form("Bachelor has negative ID %d",bachPart->GetID()));
133     return bSignAssoc;
134   }
135
136   if (!v0part->GetSecondaryVtx()) {
137     AliDebug(2,"No secondary vertex for V0 by cascade");
138     return bSignAssoc;
139   }
140
141   if (v0part->GetNDaughters()!=2) {
142     AliDebug(2,Form("No 2 daughters for V0 of current cascade (onTheFly=%d, nDaughters=%d)",v0part->GetOnFlyStatus(),v0part->GetNDaughters()));
143     return bSignAssoc;
144   }
145
146   AliVTrack *trk0 = dynamic_cast<AliVTrack*>(v0part->GetDaughter(0));
147   AliVTrack *trk1 = dynamic_cast<AliVTrack*>(v0part->GetDaughter(1));
148   if (!trk0 || !trk1) {
149     AliDebug(2,"At least one of V0daughters doesn't exist");
150     return bSignAssoc;
151   }
152
153   //if (trk0->GetLabel()<0 || trk1->GetLabel()<0) {
154   //AliDebug(2,Form("At least one of V0daughters has label negative (%d %d)",trk0->GetLabel(),trk1->GetLabel()));
155   //return bSignAssoc;
156   //}
157
158  if (trk0->GetID()<0 || trk1->GetID()<0) {
159     AliDebug(2,Form("At least one of V0 daughters has negative ID %d %d",trk0->GetID(),trk1->GetID()));
160     return bSignAssoc;
161   }
162
163   Int_t pdgCand = 4122;
164   Int_t mcLabel = -1;
165   Int_t mcLabelK0S = -1;
166   Int_t mcLabelLambda = -1;
167
168   // Lc->K0S+p and cc
169   Int_t pdgDgLctoV0bachelor[2]={2212,310}; // first bachelor, second V0
170   Int_t pdgDgV0toDaughters[2]={211,211};
171   mcLabelK0S = lcV0bachelor->MatchToMC(pdgCand,pdgDgLctoV0bachelor[1],pdgDgLctoV0bachelor,pdgDgV0toDaughters,fmcArray,kTRUE);
172
173   // Lc->Lambda+pi and cc
174   pdgDgLctoV0bachelor[0]=211, pdgDgLctoV0bachelor[1]=3122; // first bachelor, second V0
175   pdgDgV0toDaughters[0]=2212,  pdgDgV0toDaughters[1]=211;
176   mcLabelLambda = lcV0bachelor->MatchToMC(pdgCand,pdgDgLctoV0bachelor[1],pdgDgLctoV0bachelor,pdgDgV0toDaughters,fmcArray,kTRUE);
177
178   if (mcLabelK0S!=-1 && mcLabelLambda!=-1)
179     AliDebug(2,"Strange: current Lc->V0+bachelor candidate has two MC different labels!");
180
181   if (fGenLcOption==kCountK0Sp) {
182     if (mcLabelK0S!=-1) mcLabel=mcLabelK0S;
183     if (mcLabelLambda!=-1) {
184       mcLabel=-1;
185       fFake = 0;    // fake candidate
186       if (fFakeSelection==1) return bSignAssoc;
187     }
188   }
189   else if (fGenLcOption==kCountLambdapi) {
190     if (mcLabelLambda!=-1) mcLabel=mcLabelLambda;
191     if (mcLabelK0S!=-1) {
192       mcLabel=-1;
193       fFake = 0;    // fake candidate
194       if (fFakeSelection==1) return bSignAssoc;
195     }
196   }
197
198   if (mcLabel==-1) {
199     AliDebug(4,"No mcLabel found for current candidate");
200     return bSignAssoc;
201   }
202   AliDebug(1,Form("Found mcLabel (%d) for current candidate (onTheFly=%1d)",mcLabel,v0part->GetOnFlyStatus()));
203
204   if (fRecoCandidate->NumberOfFakeDaughters()>0){
205     fFake = 0;    // fake candidate
206     if (fFakeSelection==1) return bSignAssoc;
207   }
208   if (fRecoCandidate->NumberOfFakeDaughters()==0){
209     fFake = 2;    // non-fake candidate
210     if (fFakeSelection==2) return bSignAssoc;
211   }
212   
213   SetMCLabel(mcLabel); // fmcLabel=mcLabel
214   fmcPartCandidate = dynamic_cast<AliAODMCParticle*>(fmcArray->At(fmcLabel)); 
215   if (!fmcPartCandidate){
216     AliDebug(3,"No MC object for current candidate");
217     return bSignAssoc;
218   }
219  
220   bSignAssoc = kTRUE;
221   return bSignAssoc;
222
223 }
224
225 //______________________________________________
226 Bool_t AliCFVertexingHFLctoV0bachelor::GetGeneratedValuesFromMCParticle(Double_t* vectorMC) 
227 {
228   // 
229   // collecting all the necessary info (pt, y, invMassV0, cosPAwrtPVxV0, onTheFlyStatusV0) from MC particle
230   // (additional infos: pTbachelor, pTV0pos, pTV0neg, phi, dcaV0, cTV0, cT, cosPA)
231   //
232
233   Bool_t bGenValues = kFALSE;
234
235   if (fmcPartCandidate->GetNDaughters()!=2) {
236     AliDebug(2,"Lc MC particle doesn't decay in 2 daughters");
237     return bGenValues;
238   }
239
240   Int_t daughter0lc = fmcPartCandidate->GetDaughter(0);
241   Int_t daughter1lc = fmcPartCandidate->GetDaughter(1);
242   if (daughter0lc<=0 || daughter1lc<=0) {
243     AliDebug(2,"Lc daughters are not in MC array");
244     return bGenValues;
245   }
246
247   AliAODMCParticle* mcPartDaughter0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter0lc));
248   AliAODMCParticle* mcPartDaughter1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1lc));
249   if (!mcPartDaughter0 || !mcPartDaughter1) {
250     AliDebug(2,"Problems in the MC Daughters\n");
251     return bGenValues;
252   }
253
254   if ( fGenLcOption==kCountLambdapi &&
255        !(TMath::Abs(mcPartDaughter0->GetPdgCode())==3122 &&
256          TMath::Abs(mcPartDaughter1->GetPdgCode())==211) &&
257        !(TMath::Abs(mcPartDaughter1->GetPdgCode())==3122 &&
258          TMath::Abs(mcPartDaughter0->GetPdgCode())==211) ) return bGenValues;
259   if ( fGenLcOption==kCountK0Sp &&
260        !(TMath::Abs(mcPartDaughter0->GetPdgCode())==2212 &&
261          TMath::Abs(mcPartDaughter1->GetPdgCode())==311) &&
262        !(TMath::Abs(mcPartDaughter1->GetPdgCode())==2212 &&
263          TMath::Abs(mcPartDaughter0->GetPdgCode())==311) ) return bGenValues;
264
265   if ( (TMath::Abs(mcPartDaughter0->GetPdgCode())==311   &&
266         TMath::Abs(mcPartDaughter1->GetPdgCode())==2212) ||
267        (TMath::Abs(mcPartDaughter0->GetPdgCode())==3122  &&
268         TMath::Abs(mcPartDaughter1->GetPdgCode())==211) )
269     bGenValues = FillVectorFromMCarray(mcPartDaughter1,mcPartDaughter0,vectorMC);
270   else if ( (TMath::Abs(mcPartDaughter1->GetPdgCode())==311   &&
271              TMath::Abs(mcPartDaughter0->GetPdgCode())==2212) ||
272             (TMath::Abs(mcPartDaughter1->GetPdgCode())==3122  &&
273              TMath::Abs(mcPartDaughter0->GetPdgCode())==211) )
274     bGenValues = FillVectorFromMCarray(mcPartDaughter0,mcPartDaughter1,vectorMC);
275
276   if (!bGenValues)
277     AliDebug(2,"There is something wrong in filling MC vector");
278
279   return bGenValues;
280
281 }
282
283 //____________________________________________
284 Bool_t AliCFVertexingHFLctoV0bachelor::GetRecoValuesFromCandidate(Double_t *vectorReco) const
285
286   // read the variables for the container
287
288   Bool_t bFillRecoValues = kFALSE;
289
290   //Get the Lc and the V0 from Lc
291   AliAODRecoCascadeHF* lcV0bachelor = dynamic_cast<AliAODRecoCascadeHF*>(fRecoCandidate);
292
293   if (!lcV0bachelor) {
294     AliDebug(2,"Current cascade doesn't exist, skipped");
295     return bFillRecoValues;
296   }
297
298   if (!lcV0bachelor->GetSecondaryVtx()) {
299     AliDebug(2,"No secondary vertex for cascade");
300     return bFillRecoValues;
301   }
302
303   if (lcV0bachelor->GetNDaughters()!=2) {
304     AliDebug(2,Form("No 2 daughters for current cascade (nDaughters=%d)",lcV0bachelor->GetNDaughters()));
305     return bFillRecoValues;
306   }
307
308   AliAODTrack* bachelor = dynamic_cast<AliAODTrack*>(lcV0bachelor->GetBachelor());
309   AliAODv0* v0toDaughters = dynamic_cast<AliAODv0*>(lcV0bachelor->Getv0());
310   if (!bachelor || !v0toDaughters) {
311     AliDebug(2,"No V0 or bachelor in this reco candidate, skipping!");
312     return bFillRecoValues;
313   }
314
315   if (!v0toDaughters->GetSecondaryVtx()) {
316     AliDebug(2,"No secondary vertex for V0 by cascade");
317     return bFillRecoValues;
318   }
319
320   if (v0toDaughters->GetNDaughters()!=2) {
321     AliDebug(2,Form("current V0 has not 2 daughters (onTheFly=%d, nDaughters=%d)",v0toDaughters->GetOnFlyStatus(),v0toDaughters->GetNDaughters()));
322     return bFillRecoValues;
323   }
324
325   Bool_t onTheFlyStatus = v0toDaughters->GetOnFlyStatus();
326   AliAODTrack* v0positiveTrack;
327   AliAODTrack* v0negativeTrack;
328   if (onTheFlyStatus) {
329     v0positiveTrack = dynamic_cast<AliAODTrack*>(lcV0bachelor->Getv0NegativeTrack());
330     v0negativeTrack = dynamic_cast<AliAODTrack*>(lcV0bachelor->Getv0PositiveTrack());
331   } else {
332     v0positiveTrack = dynamic_cast<AliAODTrack*>(lcV0bachelor->Getv0PositiveTrack());
333     v0negativeTrack = dynamic_cast<AliAODTrack*>(lcV0bachelor->Getv0NegativeTrack());
334   }
335
336   if (!v0positiveTrack || !v0negativeTrack) {
337     AliDebug(2,"No V0daughters in this reco candidate, skipping!");
338     return bFillRecoValues;
339   }
340
341   Double_t pt = lcV0bachelor->Pt();
342   Double_t rapidity = lcV0bachelor->Y(4122);
343
344   Double_t cosPAwrtPrimVtxV0 = lcV0bachelor->CosV0PointingAngle();
345
346   //Double_t pTbachelor = bachelor->Pt();
347   Double_t pbachelor = bachelor->P();
348   //Double_t pTV0pos = v0positiveTrack->Pt();
349   //Double_t pTV0neg = v0negativeTrack->Pt();
350   Double_t phi = lcV0bachelor->Phi();
351   Double_t dcaV0 = v0toDaughters->GetDCA();
352   //Double_t cTLc = lcV0bachelor->Ct(4122); // wrt PrimVtx
353   //Double_t dcaLc = lcV0bachelor->GetDCA();
354   Double_t cosPointingAngleLc = lcV0bachelor->CosPointingAngle();
355
356   /*
357   Double_t cTV0 = 0.;
358   AliAODVertex *vtx0 = (AliAODVertex*)lcV0bachelor->GetPrimaryVtx();
359   if (!vtx0) {
360     AliDebug(2,"Candidate has not primary vtx");
361   } else {
362     Double_t primVtxPos[3] = {0.,0.,0.}; vtx0->GetXYZ(primVtxPos);
363     if (fGenLcOption==kCountK0Sp) {
364       cTV0 = v0toDaughters->Ct(310,primVtxPos);
365     } else if (fGenLcOption==kCountLambdapi) {
366       cTV0 = v0toDaughters->Ct(3122,primVtxPos);
367     }
368   }
369   */
370
371   Double_t invMassV0 = 0.;
372   if (fGenLcOption==kCountLambdapi) {
373
374     Short_t bachelorCharge = bachelor->Charge();
375     if (bachelorCharge==1) {
376       invMassV0 = v0toDaughters->MassLambda();
377     } else if (bachelorCharge==-1) {
378       invMassV0 = v0toDaughters->MassAntiLambda();
379     }
380
381   } else if (fGenLcOption==kCountK0Sp) {
382
383     invMassV0 = v0toDaughters->MassK0Short();
384
385   }
386
387   vectorReco[0]  = pt;
388   vectorReco[1]  = rapidity;
389   vectorReco[2]  = phi;
390   vectorReco[3]  = onTheFlyStatus;
391   vectorReco[4]  = fzPrimVertex;
392   vectorReco[5]  = fCentValue;
393   vectorReco[6]  = fFake; // whether the reconstructed candidate was a fake (fFake = 0) or not (fFake = 2) 
394   vectorReco[7]  = fMultiplicity;
395
396   if (fConfiguration==AliCFTaskVertexingHF::kSnail) {
397     //vectorReco[8]  = pTbachelor;
398     vectorReco[8]  = pbachelor;
399     vectorReco[9]  = v0toDaughters->Pt();
400     if (fGenLcOption==kCountLambdapi) {
401       vectorReco[10] = v0toDaughters->Y(3122);
402     } else if (fGenLcOption==kCountK0Sp) {
403       vectorReco[10] = v0toDaughters->Y(310);
404     }
405     vectorReco[11] = v0toDaughters->Phi();
406     vectorReco[12] = invMassV0;
407     vectorReco[13] = dcaV0;
408     vectorReco[14] = cosPAwrtPrimVtxV0;
409     vectorReco[15] = cosPointingAngleLc;
410     //vectorReco[16] = cTV0*1.E4; // in micron
411     //vectorReco[17] = cTLc*1.E4; // in micron
412   }
413
414   bFillRecoValues = kTRUE;
415
416   return bFillRecoValues;
417 }
418
419 //_____________________________________________________________
420 Bool_t AliCFVertexingHFLctoV0bachelor::CheckMCChannelDecay() const
421
422   // check the required decay channel
423
424   Bool_t checkCD = kFALSE;
425   
426   if (fmcPartCandidate->GetNDaughters()!=2) {
427     AliDebug(2, Form("The MC particle doesn't decay in 2 particles, skipping!!"));
428     return checkCD;
429   }
430
431   Int_t daughter0 = fmcPartCandidate->GetDaughter(0);
432   Int_t daughter1 = fmcPartCandidate->GetDaughter(1);
433   if (daughter0<=0 || daughter1<=0){
434     AliDebug(2, Form("The MC particle doesn't have correct daughters, skipping!!"));
435     return checkCD;
436   }
437   AliAODMCParticle* mcPartDaughter0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter0));
438   AliAODMCParticle* mcPartDaughter1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1));
439   if (!mcPartDaughter0 || !mcPartDaughter1) {
440     AliDebug(2,"Problems in the MC Daughters\n");
441     return checkCD;
442   }
443
444   // Lc -> Lambda + pion AND cc
445   if (fGenLcOption==kCountLambdapi) {
446
447     if (!(TMath::Abs(mcPartDaughter0->GetPdgCode())==3122 &&
448           TMath::Abs(mcPartDaughter1->GetPdgCode())==211) && 
449         !(TMath::Abs(mcPartDaughter0->GetPdgCode())==211  &&
450           TMath::Abs(mcPartDaughter1->GetPdgCode())==3122)) {
451       AliDebug(2, "The Lc MC doesn't decay in Lambda+pion (or cc), skipping!!");
452       return checkCD;  
453     }
454
455     if (TMath::Abs(mcPartDaughter0->GetPdgCode())==3122) {
456       mcPartDaughter0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1)); // the bachelor
457       mcPartDaughter1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter0)); // the V0
458     }
459     if (!mcPartDaughter0 || !mcPartDaughter1) {
460       AliDebug(2,"Problems in the MC Daughters\n");
461       return checkCD;
462     }
463
464     if (mcPartDaughter1->GetNDaughters()!=2) {
465       AliDebug(2, "The Lambda MC particle doesn't decay in 2 particles, skipping!!");
466       return checkCD;
467     }
468
469     Int_t daughter1D0 = mcPartDaughter1->GetDaughter(0);
470     Int_t daughter1D1 = mcPartDaughter1->GetDaughter(1);
471     if (daughter1D0<=0 || daughter1D1<=0) {
472       AliDebug(2, Form("The Lambda MC particle doesn't have correct daughters, skipping!!"));
473       return checkCD;
474     }
475
476     AliAODMCParticle* mcPartDaughter1D0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1D0));
477     AliAODMCParticle* mcPartDaughter1D1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1D1));
478     if(!mcPartDaughter1D0 || !mcPartDaughter1D1) {
479       AliError("The Lambda daughter particle not found in MC array");
480       return checkCD;
481     }
482
483     if (!(TMath::Abs(mcPartDaughter1D0->GetPdgCode())==211   &&
484           TMath::Abs(mcPartDaughter1D1->GetPdgCode())==2212) &&
485         !(TMath::Abs(mcPartDaughter1D0->GetPdgCode())==2212  &&
486           TMath::Abs(mcPartDaughter1D1->GetPdgCode())==211)) {
487       AliDebug(2, "The Lambda MC doesn't decay in pi+proton (or cc), skipping!!");
488       return checkCD;
489     }
490
491   } else if (fGenLcOption==kCountK0Sp) { // Lc -> K0bar + proton AND cc
492
493     if (!(TMath::Abs(mcPartDaughter0->GetPdgCode())==311   &&
494           TMath::Abs(mcPartDaughter1->GetPdgCode())==2212) &&
495         !(TMath::Abs(mcPartDaughter0->GetPdgCode())==2212  &&
496           TMath::Abs(mcPartDaughter1->GetPdgCode())==311)) {
497       AliDebug(2, "The Lc MC doesn't decay in K0+proton (or cc), skipping!!");
498       return checkCD;  
499     }
500     
501     if (TMath::Abs(mcPartDaughter0->GetPdgCode())==311) {
502       mcPartDaughter0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1)); // the bachelor
503       mcPartDaughter1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter0)); // the V0
504     }
505     if (!mcPartDaughter0 || !mcPartDaughter1) {
506       AliDebug(2,"Problems in the MC Daughters after swapping V0 and bachelor\n");
507       return checkCD;
508     }
509
510     if (mcPartDaughter1->GetNDaughters()!=1) {
511       AliDebug(2, "The K0/K0bar MC particle doesn't decay in 1 particles, skipping!!");
512       return checkCD;
513     }
514
515     Int_t daughter = mcPartDaughter1->GetDaughter(0);
516     if (daughter<=0) {
517       AliDebug(2, Form("The K0/K0bar MC particle doesn't have correct daughter, skipping!!"));
518       return checkCD;
519     }
520
521     AliAODMCParticle* mcPartDaughter = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter));
522     if(!mcPartDaughter){
523       AliError("The K0/K0bar daughter particle not found in MC array");
524       return checkCD;
525     }
526
527     if (!(TMath::Abs(mcPartDaughter->GetPdgCode())==310)) {
528       AliDebug(2, "The K0/K0bar MC doesn't go in K0S, skipping!!");
529       return checkCD;
530     }
531
532     if (mcPartDaughter->GetNDaughters()!=2) {
533       AliDebug(2, "The K0S MC doesn't decay in 2 particles, skipping!!");
534       return checkCD;
535     }
536
537     Int_t daughterD0 = mcPartDaughter->GetDaughter(0);
538     Int_t daughterD1 = mcPartDaughter->GetDaughter(1);
539     if (daughterD0<=0 || daughterD1<=0) {
540       AliDebug(2, Form("The K0S MC particle doesn't have correct daughters, skipping!!"));
541       return checkCD;
542     }
543
544     AliAODMCParticle* mcPartDaughterD0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughterD0));
545     AliAODMCParticle* mcPartDaughterD1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughterD1));
546     if (!mcPartDaughterD0 || !mcPartDaughterD1) {
547       AliError("Daughter particle not found in MC array");
548       return checkCD;
549     }
550
551     if (! ( TMath::Abs(mcPartDaughterD0->GetPdgCode())==211 &&
552             TMath::Abs(mcPartDaughterD1->GetPdgCode())==211 ) ) {
553       AliDebug(2, "The K0S MC doesn't decay in pi+ pi-, skipping!!");
554       return checkCD;
555     }
556
557   }
558   
559   checkCD = kTRUE;
560   return checkCD;
561   
562 }
563
564 //_____________________________________________________________
565 Double_t AliCFVertexingHFLctoV0bachelor::GetEtaProng(Int_t iProng) const 
566 {
567   //
568   // getting eta of the prong - overload the mother class method
569   //
570
571   Double_t etaProng =-9999;
572
573   if (!fRecoCandidate) {
574     AliDebug(2,"No reco candidate selected");
575     return etaProng;
576   }
577
578   AliAODRecoCascadeHF* lcV0bachelor = (AliAODRecoCascadeHF*)fRecoCandidate;
579   AliAODTrack* bachelor = (AliAODTrack*)lcV0bachelor->GetBachelor();
580   AliAODTrack* v0Pos = (AliAODTrack*)lcV0bachelor->Getv0PositiveTrack();
581   AliAODTrack* v0Neg = (AliAODTrack*)lcV0bachelor->Getv0NegativeTrack();
582   if (!(lcV0bachelor->Getv0()) || !bachelor || !v0Pos || !v0Neg) {
583     AliDebug(2,"No V0 for this reco candidate selected");
584     return etaProng;
585   }
586
587   if (iProng==0) etaProng = bachelor->Eta();
588   else if (iProng==1) etaProng = v0Pos->Eta();
589   else if (iProng==2) etaProng = v0Neg->Eta();
590
591   AliDebug(4,Form("Eta value for prong number %1d = %f",iProng,etaProng));
592
593   return etaProng;
594
595 }
596
597 //_____________________________________________________________
598
599 Double_t AliCFVertexingHFLctoV0bachelor::GetPtProng(Int_t iProng) const 
600 {
601   //
602   // getting pt of the prong
603   //
604
605   Double_t ptProng=-9999.;
606
607   if (!fRecoCandidate) {
608     AliDebug(2,"No reco candidate selected");
609     return ptProng;
610   }
611
612   AliAODRecoCascadeHF* lcV0bachelor = (AliAODRecoCascadeHF*)fRecoCandidate;
613   AliAODTrack* bachelor = (AliAODTrack*)lcV0bachelor->GetBachelor();
614   AliAODTrack* v0Pos = (AliAODTrack*)lcV0bachelor->Getv0PositiveTrack();
615   AliAODTrack* v0Neg = (AliAODTrack*)lcV0bachelor->Getv0NegativeTrack();
616   if (!(lcV0bachelor->Getv0()) || !bachelor || !v0Pos || !v0Neg) {
617     AliDebug(2,"No V0 for this reco candidate selected");
618     return ptProng;
619   }
620
621   if (iProng==0) ptProng = bachelor->Pt();
622   else if (iProng==1) ptProng = v0Pos->Pt();
623   else if (iProng==2) ptProng = v0Neg->Pt();
624     
625   AliDebug(4,Form("Pt value for prong number %1d = %f",iProng,ptProng));
626
627   return ptProng;
628   
629 }
630
631 //_____________________________________________________________
632
633 Double_t AliCFVertexingHFLctoV0bachelor::Ctau(AliAODMCParticle *mcPartCandidate)
634 {
635
636   Double_t cTau = 999999.;
637
638   Int_t daughterD0 = mcPartCandidate->GetDaughter(0);
639   Int_t daughterD1 = mcPartCandidate->GetDaughter(1);
640   if (daughterD0<=0 || daughterD1<=0) {
641     AliDebug(2, Form("The Lc MC particle doesn't have correct daughters, skipping!!"));
642     return cTau;
643   }
644
645   AliAODMCParticle *mcPartDaughter0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughterD0));
646   AliAODMCParticle *mcPartDaughter1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughterD1));
647   if (!mcPartDaughter0 || !mcPartDaughter1) {
648     AliDebug(2,"The candidate daughter particles not found in MC array");
649     return cTau;
650   }
651
652   Double_t vtx1[3] = {0,0,0};   // primary vertex               
653   Bool_t hasProdVertex = mcPartCandidate->XvYvZv(vtx1);  // cm
654
655   Double_t vtx1daughter[3] = {0,0,0};   // secondary vertex
656   Bool_t v0Vertex = mcPartDaughter0->XvYvZv(vtx1daughter);  //cm
657   Double_t vtx2daughter[3] = {0,0,0};   // secondary vertex
658   Bool_t bachVertex = hasProdVertex && mcPartDaughter1->XvYvZv(vtx2daughter);  //cm
659
660   if (!hasProdVertex || !v0Vertex || !bachVertex) {
661     AliDebug(2,"At least one of Prim.vtx, V0vtx, BachelorVtx doesn't exist!");
662     return cTau;
663   }
664
665   if (TMath::Abs(vtx1daughter[0]-vtx2daughter[0])>1E-5 ||
666       TMath::Abs(vtx1daughter[1]-vtx2daughter[1])>1E-5 ||
667       TMath::Abs(vtx1daughter[2]-vtx2daughter[2])>1E-5) {
668     AliDebug(2,"Bachelor and V0 haven't common vtx!");
669     return cTau;
670   }
671
672   Double_t decayLength = 0.;
673   for (Int_t ii=0; ii<3; ii++) decayLength += (vtx1daughter[ii]-vtx1[ii])*(vtx1daughter[ii]-vtx1[ii]);
674   decayLength = TMath::Sqrt(decayLength);
675
676   cTau = decayLength * mcPartCandidate->M()/mcPartCandidate->P();
677
678   AliDebug(2,Form(" cTau(4122)=%f",cTau));
679
680   return cTau;
681
682 }
683
684 //------------
685 Bool_t AliCFVertexingHFLctoV0bachelor::SetLabelArray()
686 {
687   //
688   // setting the label arrays
689   //
690
691   Bool_t checkCD = kFALSE;
692   
693   if (fmcPartCandidate->GetNDaughters()!=2) {
694     AliDebug(2, Form("The MC particle have %d daughters (not 2), skipping!!",fmcPartCandidate->GetNDaughters()));
695     fmcPartCandidate->Print();
696     return checkCD;
697   }
698
699   Int_t daughter0 = fmcPartCandidate->GetDaughter(0);
700   Int_t daughter1 = fmcPartCandidate->GetDaughter(1);
701   if (daughter0<=0 || daughter1<=0){
702     AliDebug(2, Form("The MC particle doesn't have correct daughters, skipping!!"));
703     return checkCD;
704   }
705
706   AliAODMCParticle* mcPartDaughter0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter0));
707   AliAODMCParticle* mcPartDaughter1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1));
708   if (!mcPartDaughter0 || !mcPartDaughter1) {
709     AliDebug(2,"Problems in the MC Daughters\n");
710     return checkCD;
711   }
712
713
714   fLabelArray = new Int_t[fProngs];
715
716   if (fGenLcOption==kCountLambdapi) { // Lc -> Lambda + pion OR cc
717
718     if (!(TMath::Abs(mcPartDaughter0->GetPdgCode())==3122 &&
719           TMath::Abs(mcPartDaughter1->GetPdgCode())==211) && 
720         !(TMath::Abs(mcPartDaughter0->GetPdgCode())==211  &&
721           TMath::Abs(mcPartDaughter1->GetPdgCode())==3122)) {
722       AliDebug(2, "The Lc MC doesn't decay in Lambda+pion (or cc), skipping!!");
723       delete [] fLabelArray;
724       fLabelArray = 0x0;
725       return checkCD;  
726     }
727
728     // it is Lc -> Lambda + pion OR cc
729     if (TMath::Abs(mcPartDaughter0->GetPdgCode())==3122) {
730       mcPartDaughter0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1)); // the bachelor
731       mcPartDaughter1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter0)); // the V0
732       Int_t daughterTemp = daughter0;
733       daughter0 = daughter1; // the bachelor label
734       daughter1 = daughterTemp; // the V0 label
735     }
736
737     if (mcPartDaughter1->GetNDaughters()!=2) {
738       AliDebug(2, "The Lambda MC particle doesn't decay in 2 particles, skipping!!");
739       delete [] fLabelArray;
740       fLabelArray = 0x0;
741       return checkCD;
742     }
743
744     Int_t daughter1D0 = mcPartDaughter1->GetDaughter(0);
745     Int_t daughter1D1 = mcPartDaughter1->GetDaughter(1);
746     if (daughter1D0<=0 || daughter1D1<=0) {
747       AliDebug(2, Form("The Lambda MC particle doesn't have correct daughters, skipping!!"));
748       delete [] fLabelArray;
749       fLabelArray = 0x0;
750       return checkCD;
751     }
752
753     AliAODMCParticle* mcPartDaughter1D0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1D0));
754     AliAODMCParticle* mcPartDaughter1D1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1D1));
755     if (!mcPartDaughter1D0 || !mcPartDaughter1D1) {
756       AliError("The Lambda daughter particles not found in MC array");
757       delete [] fLabelArray;
758       fLabelArray = 0x0;
759       return checkCD;
760     }
761
762     if (!(TMath::Abs(mcPartDaughter1D0->GetPdgCode())==211   &&
763           TMath::Abs(mcPartDaughter1D1->GetPdgCode())==2212) &&
764         !(TMath::Abs(mcPartDaughter1D0->GetPdgCode())==2212  &&
765           TMath::Abs(mcPartDaughter1D1->GetPdgCode())==211)) {
766       AliDebug(2, "The Lambda MC doesn't decay in pi+proton (or cc), skipping!!");
767       delete [] fLabelArray;
768       fLabelArray = 0x0;
769       return checkCD;
770     }
771
772     // Lambda -> p+pi OR cc
773
774     fLabelArray[0] = daughter0;//mcPartDaughter0->GetLabel(); // bachelor
775
776     if (fmcPartCandidate->Charge()>0) {
777
778       if (mcPartDaughter1D0->GetPdgCode()==2212) {
779         fLabelArray[1] = daughter1D0;//mcPartDaughter1D0->GetLabel(); // proton
780         fLabelArray[2] = daughter1D1;//mcPartDaughter1D1->GetLabel(); // pion
781       } else if (mcPartDaughter1D1->GetPdgCode()==2212) {
782         fLabelArray[1] = daughter1D1;//mcPartDaughter1D1->GetLabel(); // proton
783         fLabelArray[2] = daughter1D0;//mcPartDaughter1D0->GetLabel(); // pion
784       }
785
786     } else if (fmcPartCandidate->Charge()<0) {
787
788       if (mcPartDaughter1D0->GetPdgCode()==211) {
789         fLabelArray[1] = daughter1D0;//mcPartDaughter1D0->GetLabel(); // pion
790         fLabelArray[2] = daughter1D1;//mcPartDaughter1D1->GetLabel(); // proton
791       } else if (mcPartDaughter1D1->GetPdgCode()==211) {
792         fLabelArray[1] = daughter1D1;//mcPartDaughter1D1->GetLabel(); // pion
793         fLabelArray[2] = daughter1D0;//mcPartDaughter1D0->GetLabel(); // proton
794       }
795
796     }
797
798   } else if (fGenLcOption==kCountK0Sp) { // Lc -> K0bar + proton OR cc
799
800     if (!(TMath::Abs(mcPartDaughter0->GetPdgCode())==311   &&
801           TMath::Abs(mcPartDaughter1->GetPdgCode())==2212) &&
802         !(TMath::Abs(mcPartDaughter0->GetPdgCode())==2212  &&
803           TMath::Abs(mcPartDaughter1->GetPdgCode())==311)) {
804       AliDebug(2, "The Lc MC doesn't decay in K0bar+proton (or cc), skipping!!");
805       delete [] fLabelArray;
806       fLabelArray = 0x0;
807       return checkCD;  
808     }
809
810     if (TMath::Abs(mcPartDaughter0->GetPdgCode())==311) {
811       mcPartDaughter0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter1)); // the bachelor
812       mcPartDaughter1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter0)); // the V0
813       Int_t daughterTemp = daughter0;
814       daughter0 = daughter1; // the bachelor label
815       daughter1 = daughterTemp; // the V0 label
816     }
817
818     if (mcPartDaughter1->GetNDaughters()!=1) {
819       AliDebug(2, "The K0/K0bar MC particle doesn't decay in 1 particles, skipping!!");
820       delete [] fLabelArray;
821       fLabelArray = 0x0;
822       return checkCD;
823     }
824
825     Int_t daughter = mcPartDaughter1->GetDaughter(0);
826     if (daughter<=0) {
827       AliDebug(2, Form("The K0/K0bar MC particle doesn't have correct daughter, skipping!!"));
828       delete [] fLabelArray;
829       fLabelArray = 0x0;
830       return checkCD;
831     }
832
833     AliAODMCParticle* mcPartDaughter = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughter));
834     if (!mcPartDaughter) {
835       AliError("The K0/K0bar daughter particle not found in MC array");
836       delete [] fLabelArray;
837       fLabelArray = 0x0;
838       return checkCD;
839     }
840
841     if (!(TMath::Abs(mcPartDaughter->GetPdgCode())==310)) {
842       AliDebug(2, "The K0/K0bar MC doesn't go in K0S, skipping!!");
843       delete [] fLabelArray;
844       fLabelArray = 0x0;
845       return checkCD;
846     }
847
848     if (mcPartDaughter->GetNDaughters()!=2) {
849       AliDebug(2, "The K0S MC doesn't decay in 2 particles, skipping!!");
850       delete [] fLabelArray;
851       fLabelArray = 0x0;
852       return checkCD;
853     }
854
855     Int_t daughterD0 = mcPartDaughter->GetDaughter(0);
856     Int_t daughterD1 = mcPartDaughter->GetDaughter(1);
857     if (daughterD0<=0 || daughterD1<=0) {
858       AliDebug(2, Form("The K0S MC particle doesn't have correct daughters, skipping!!"));
859       delete [] fLabelArray;
860       fLabelArray = 0x0;
861       return checkCD;
862     }
863
864     AliAODMCParticle* mcPartDaughterD0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughterD0));
865     AliAODMCParticle* mcPartDaughterD1 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughterD1));
866     if (!mcPartDaughterD0 || !mcPartDaughterD1) {
867       AliError("The K0S daughter particles not found in MC array");
868       delete [] fLabelArray;
869       fLabelArray = 0x0;
870       return checkCD;
871     }
872
873     if (! ( TMath::Abs(mcPartDaughterD0->GetPdgCode())==211 &&
874             TMath::Abs(mcPartDaughterD1->GetPdgCode())==211 ) ) {
875       AliDebug(2, "The K0S MC doesn't decay in pi+ pi-, skipping!!");
876       delete [] fLabelArray;
877       fLabelArray = 0x0;
878       return checkCD;
879     }
880
881     // K0S -> pi+ pi-
882
883     fLabelArray[0] = daughter0;//mcPartDaughter0->GetLabel(); // bachelor
884
885     if (mcPartDaughterD0->GetPdgCode()==211) {
886       fLabelArray[1] = daughterD0;//mcPartDaughterD0->GetLabel(); // pi+
887       fLabelArray[2] = daughterD1;//mcPartDaughterD1->GetLabel(); // pi-
888       AliDebug(2,Form(" daughter0=%d ------ daughter1=%d ------ dg0->GetLabel()=%d ------ dg1->GetLabel()=%d ",daughterD0,daughterD1,mcPartDaughterD0->GetLabel(),mcPartDaughterD1->GetLabel()));
889     } else if (mcPartDaughterD1->GetPdgCode()==211) {
890       fLabelArray[1] = daughterD1;//mcPartDaughterD1->GetLabel(); // pi+
891       fLabelArray[2] = daughterD0;//mcPartDaughterD0->GetLabel(); // pi-
892       AliDebug(2,Form(" daughter0=%d ------ daughter1=%d ------ dg0->GetLabel()=%d ------ dg1->GetLabel()=%d ",daughterD1,daughterD0,mcPartDaughterD1->GetLabel(),mcPartDaughterD0->GetLabel()));
893     }
894   }
895
896   AliDebug(2,Form(" label0=%d, label1=%d, label2=%d",fLabelArray[0],fLabelArray[1],fLabelArray[2]));
897   
898   SetAccCut(); // setting the pt and eta acceptance cuts
899
900   checkCD = kTRUE;
901   return checkCD;
902
903 }
904 //____________________________________________
905 Bool_t AliCFVertexingHFLctoV0bachelor::FillVectorFromMCarray(AliAODMCParticle *mcPartDaughterBachelor,
906                                                              AliAODMCParticle *mcPartDaughterK0,
907                                                              Double_t *vectorMC)
908 {
909   // fill the vector
910
911   Bool_t bGenValues = kFALSE;
912
913   AliAODMCParticle *mcPartV0DaughterPos = dynamic_cast<AliAODMCParticle*>(fmcArray->At(fLabelArray[1]));
914   AliAODMCParticle *mcPartV0DaughterNeg = dynamic_cast<AliAODMCParticle*>(fmcArray->At(fLabelArray[2]));
915   AliAODMCParticle *mcPartDaughterV0 = 0x0;
916
917   if(!mcPartV0DaughterPos || !mcPartV0DaughterNeg) return bGenValues;
918
919   if (TMath::Abs(mcPartDaughterK0->GetPdgCode())==311) {
920     Int_t daughterK0 = mcPartDaughterK0->GetDaughter(0);
921     if (daughterK0<=0) {
922       AliDebug(2, Form("The K0/K0bar particle doesn't have correct daughter, skipping!!"));
923       return bGenValues;
924     }
925     mcPartDaughterV0 = dynamic_cast<AliAODMCParticle*>(fmcArray->At(daughterK0));
926     if (!mcPartDaughterV0) {
927       AliDebug(2,"The K0/K0bar daughter particle not found in MC array");
928       return bGenValues;
929     }
930     if (TMath::Abs(mcPartDaughterV0->GetPdgCode())!=310) {
931       AliDebug(2,"The K0/K0bar daughter particle is not a K0S");
932       return bGenValues;
933     }
934   } else if (TMath::Abs(mcPartDaughterK0->GetPdgCode())==3122) {
935     mcPartDaughterV0 = dynamic_cast<AliAODMCParticle*>(mcPartDaughterK0);
936     if (!mcPartDaughterV0) {
937       AliDebug(2,"The Lambda particle not found in MC array");
938       return bGenValues;
939     }
940   }
941
942   if (!mcPartDaughterV0) {
943     AliDebug(2,"V0 particle not found in MC array");
944     return bGenValues;
945   }
946
947   //Double_t cTLc = Ctau(fmcPartCandidate); // by default wrt Primary Vtx
948   //Double_t pTbach = mcPartDaughterBachelor->Pt(); // get the bachelor pT
949   Double_t pbach = mcPartDaughterBachelor->P(); // get the bachelor p
950
951   Double_t vtx1[3] = {0,0,0};   // primary vertex               
952   Bool_t hasPrimVtx = fmcPartCandidate->XvYvZv(vtx1);  // cm
953
954   // getting vertex from daughters
955   Double_t vtx1daughter0[3] = {0,0,0};   // secondary vertex from daughter 0
956   Bool_t hasSecVtx1 = mcPartDaughterBachelor->XvYvZv(vtx1daughter0);  //cm
957   Double_t vtx1daughter1[3] = {0,0,0};   // secondary vertex from daughter 1
958   Bool_t hasSecVtx2 = mcPartDaughterV0->XvYvZv(vtx1daughter1);  //cm
959   if (!hasPrimVtx || !hasSecVtx1 || !hasSecVtx2) {
960     AliDebug(2,"At least one of Prim.vtx, V0vtx, BachelorVtx doesn't exist!");
961     //return bGenValues;
962   }
963
964   if (TMath::Abs(vtx1daughter0[0]-vtx1daughter1[0])>1E-5 ||
965       TMath::Abs(vtx1daughter0[1]-vtx1daughter1[1])>1E-5 ||
966       TMath::Abs(vtx1daughter0[2]-vtx1daughter1[2])>1E-5) {
967     AliError("Daughters have different secondary vertex, skipping the track");
968     //return bGenValues;
969   }
970
971   // getting the momentum from the daughters
972   Double_t px1[2] = {mcPartDaughterBachelor->Px(), mcPartDaughterV0->Px()};
973   Double_t py1[2] = {mcPartDaughterBachelor->Py(), mcPartDaughterV0->Py()};
974   Double_t pz1[2] = {mcPartDaughterBachelor->Pz(), mcPartDaughterV0->Pz()};
975
976   Int_t nprongs = 2;
977   Short_t charge = mcPartDaughterBachelor->Charge();
978   Double_t d0[2] = {0.,0.};
979   AliAODRecoDecayHF* decayLc = new AliAODRecoDecayHF(vtx1,vtx1daughter0,nprongs,charge,px1,py1,pz1,d0);
980   Double_t cosPAwrtPrimVtxLc = decayLc->CosPointingAngle();
981   delete decayLc;
982
983   // getting vertex from daughters
984   Double_t vtx2daughter0[3] = {0,0,0};   // secondary vertex from daughter 0
985   Bool_t hasSecVtx3 = mcPartV0DaughterPos->XvYvZv(vtx2daughter0);  //cm
986   Double_t vtx2daughter1[3] = {0,0,0};   // secondary vertex from daughter 1
987   Bool_t hasSecVtx4 = mcPartV0DaughterNeg->XvYvZv(vtx2daughter1);  //cm
988   if (!hasSecVtx3 || !hasSecVtx4) {
989     AliDebug(2,"At least one of V0Posvtx, V0Negtx doesn't exist!");
990     //return bGenValues;
991   }
992
993   if (TMath::Abs(vtx2daughter0[0]-vtx2daughter1[0])>1E-5 ||
994       TMath::Abs(vtx2daughter0[1]-vtx2daughter1[1])>1E-5 ||
995       TMath::Abs(vtx2daughter0[2]-vtx2daughter1[2])>1E-5) {
996     AliError("Daughters have different secondary vertex, skipping the track");
997     //return bGenValues;
998   }
999
1000   // getting the momentum from the daughters
1001   Double_t px[2] = {mcPartV0DaughterPos->Px(), mcPartV0DaughterNeg->Px()};              
1002   Double_t py[2] = {mcPartV0DaughterPos->Py(), mcPartV0DaughterNeg->Py()};              
1003   Double_t pz[2] = {mcPartV0DaughterPos->Pz(), mcPartV0DaughterNeg->Pz()};
1004
1005   nprongs = 2;
1006   charge = 0;
1007   AliAODRecoDecayHF* decay = new AliAODRecoDecayHF(vtx1,vtx2daughter0,nprongs,charge,px,py,pz,d0);
1008   Double_t cosPAwrtPrimVtxV0 = decay->CosPointingAngle();
1009   /*
1010   Double_t cTV0 = 0.; //ct
1011   if (fGenLcOption==kCountK0Sp) {
1012     cTV0 = decay->Ct(310); // by default wrt Primary Vtx
1013   } else if (fGenLcOption==kCountLambdapi) {
1014     cTV0 = decay->Ct(3122); // by default wrt Primary Vtx
1015   }
1016   */
1017
1018   Double_t invMass = 0.; //invMass
1019   if (fGenLcOption==kCountK0Sp) {
1020     invMass = decay->InvMass2Prongs(0,1,211,211);
1021   } else if (fGenLcOption==kCountLambdapi) {
1022     if (fmcPartCandidate->GetPdgCode() == 4122)
1023       invMass = decay->InvMass2Prongs(0,1,2212,211);
1024     else if (fmcPartCandidate->GetPdgCode() ==-4122)
1025       invMass = decay->InvMass2Prongs(0,1,211,2212);
1026   }
1027   delete decay;
1028
1029   vectorMC[0]  = fmcPartCandidate->Pt();
1030   vectorMC[1]  = fmcPartCandidate->Y() ;
1031   vectorMC[2]  = fmcPartCandidate->Phi();
1032   vectorMC[3]  = 0; // dummy value x MC, onTheFlyStatus
1033   vectorMC[4]  = fzMCVertex;
1034   vectorMC[5]  = fCentValue; // reconstructed centrality
1035   vectorMC[6]  = 1; // dummy value x MC, fFake
1036   vectorMC[7]  = fMultiplicity; // reconstructed multiplicity
1037
1038   if (fConfiguration==AliCFTaskVertexingHF::kSnail) {
1039     //vectorMC[8]  = pTbach;
1040     vectorMC[8]  = pbach;
1041     vectorMC[9]  = mcPartDaughterV0->Pt();
1042     vectorMC[10] = mcPartDaughterV0->Y();
1043     vectorMC[11] = mcPartDaughterV0->Phi();
1044     vectorMC[12] = invMass;
1045     vectorMC[13] = 0; // dummy value x MC, V0 DCA
1046     vectorMC[14] = cosPAwrtPrimVtxV0;
1047     vectorMC[15] = cosPAwrtPrimVtxLc;
1048     //vectorMC[16] = cTV0*1.E4; // in micron
1049     //vectorMC[17] = cTLc*1.E4; // in micron
1050   }
1051
1052   bGenValues = kTRUE;
1053   return bGenValues;
1054
1055 }