1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 //_________________________________________________________________________
18 // Class for reading data (Kinematics) in order to do prompt gamma
19 // or other particle identification and correlations
20 // Separates generated particles into charged (CTS)
21 // and neutral (PHOS or EMCAL acceptance)
22 // Now, it only works with data stored in Kinematics.root and
23 // not in filtered Kinematics branch in AODs
25 //*-- Author: Gustavo Conesa (LNF-INFN)
26 //////////////////////////////////////////////////////////////////////////////
29 // --- ROOT system ---
31 #include <TParticle.h>
32 #include <TDatabasePDG.h>
35 #include "TParticle.h"
36 //#include "Riostream.h"
38 //---- ANALYSIS system ----
39 #include "AliCaloTrackMCReader.h"
40 #include "AliGenEventHeader.h"
42 #include "AliVCluster.h"
43 #include "AliAODTrack.h"
44 #include "AliAODEvent.h"
45 #include "AliFiducialCut.h"
46 #include "AliMCAnalysisUtils.h"
48 ClassImp(AliCaloTrackMCReader)
50 //____________________________________________________________________________
51 AliCaloTrackMCReader::AliCaloTrackMCReader() :
52 AliCaloTrackReader(), fDecayPi0(0),
53 fNeutralParticlesArray(0x0), fChargedParticlesArray(0x0),
54 fStatusArray(0x0), fKeepAllStatus(0), fCheckOverlap(0),
55 fEMCALOverlapAngle(0),fPHOSOverlapAngle(0), fIndex2ndPhoton(0),
56 fOnlyGeneratorParticles(kTRUE)
60 //Initialize parameters
63 //_________________________________
64 AliCaloTrackMCReader::~AliCaloTrackMCReader() {
67 if(fChargedParticlesArray) delete fChargedParticlesArray ;
68 if(fNeutralParticlesArray) delete fNeutralParticlesArray ;
69 if(fStatusArray) delete fStatusArray ;
73 //____________________________________________________________________________
74 void AliCaloTrackMCReader::GetVertex(Double_t v[3]) const {
75 //Return vertex position
78 GetGenEventHeader()->PrimaryVertex(pv);
85 //_______________________________________________________________
86 void AliCaloTrackMCReader::InitParameters()
89 //Initialize the parameters of the analysis.
93 fChargedParticlesArray = new TArrayI(1);
94 fChargedParticlesArray->SetAt(11,0);
95 //Int_t pdgarray[]={12,14,16};// skip neutrinos
96 //fNeutralParticlesArray = new TArrayI(3, pdgarray);
97 fNeutralParticlesArray = new TArrayI(3);
98 fNeutralParticlesArray->SetAt(12,0); fNeutralParticlesArray->SetAt(14,1); fNeutralParticlesArray->SetAt(16,2);
99 fStatusArray = new TArrayI(1);
100 fStatusArray->SetAt(1,0);
102 fOnlyGeneratorParticles = kTRUE;
103 fKeepAllStatus = kTRUE;
105 fCheckOverlap = kFALSE;
106 fEMCALOverlapAngle = 2.5 * TMath::DegToRad();
107 fPHOSOverlapAngle = 0.5 * TMath::DegToRad();
108 fIndex2ndPhoton = -1;
112 fReadAODMCParticles = kFALSE; //This class only works with Kinematics.root input.
114 //For this reader we own the objects of the arrays
115 fAODCTS ->SetOwner(kTRUE);
116 fAODEMCAL->SetOwner(kTRUE);
117 fAODPHOS ->SetOwner(kTRUE);
121 //____________________________________________________________________________
122 void AliCaloTrackMCReader::CheckOverlap(const Float_t anglethres, const Int_t imom, Int_t & iPrimary, Int_t & index, TLorentzVector & mom, Int_t & pdg) {
123 //Check overlap of decay photons
124 if( fIndex2ndPhoton==iPrimary ){
128 else fIndex2ndPhoton=-1;
133 TLorentzVector ph1, ph2;
134 TParticle *meson = GetStack()->Particle(imom);
135 Int_t mepdg = meson->GetPdgCode();
136 Int_t idaug1 = meson->GetFirstDaughter();
137 if((mepdg == 111 || mepdg == 221 ) && meson->GetNDaughters() == 2){ //Check only decay in 2 photons
138 TParticle * d1 = GetStack()->Particle(idaug1);
139 TParticle *d2 = GetStack()->Particle(idaug1+1);
140 if(d1->GetPdgCode() == 22 && d2->GetPdgCode() == 22 ){
143 //printf("angle %2.2f\n",ph1.Angle(ph2.Vect()));
145 if(anglethres > ph1.Angle(ph2.Vect())){
149 meson->Momentum(mom);
150 //printf("Overlap:: pt %2.2f, phi %2.2f, eta %2.2f\n",mom.Pt(),mom.Phi(),mom.Eta());
151 if(iPrimary == idaug1) iPrimary++; //skip next photon in list
154 //Do not check overlapping for next decay photon from same meson
155 if(iPrimary == idaug1) {fIndex2ndPhoton = idaug1+1;
160 }//Meson Decay with 2 photon daughters
163 //____________________________________________________________________________
164 void AliCaloTrackMCReader::FillCalorimeters(Int_t & iParticle, TParticle* particle, TLorentzVector &momentum) {
165 //Fill AODCaloClusters or TParticles lists of PHOS or EMCAL
167 if(fFillPHOS && momentum.Pt() > fPHOSPtMin){
169 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS")) return;
171 Int_t index = iParticle ;
172 Int_t pdg = TMath::Abs(particle->GetPdgCode());
174 CheckOverlap(fPHOSOverlapAngle,particle->GetFirstMother(),index, iParticle, momentum, pdg);
176 Char_t ttype= AliVCluster::kPHOSNeutral;
177 Int_t labels[] = {index};
178 Float_t x[] = {momentum.X(), momentum.Y(), momentum.Z()};
179 //Create object and write it to file
180 AliAODCaloCluster *calo = new AliAODCaloCluster(index,1,labels,momentum.E(), x, NULL, ttype, 0);
182 SetCaloClusterPID(pdg,calo) ;
183 if(fDebug > 3 && momentum.Pt() > 0.2)
184 printf("AliCaloTrackMCReader::FillCalorimeters() - PHOS : Selected cluster %s E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
185 particle->GetName(),momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
186 fAODPHOS->Add(calo);//reference the selected object to the list
190 if(fFillEMCAL && momentum.Pt() > fEMCALPtMin){
192 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) return;
194 Int_t index = iParticle ;
195 Int_t pdg = TMath::Abs(particle->GetPdgCode());
198 CheckOverlap(fEMCALOverlapAngle,particle->GetFirstMother(),iParticle, index, momentum, pdg);
200 Char_t ttype= AliVCluster::kEMCALClusterv1;
201 Int_t labels[] = {index};
202 Float_t x[] = {momentum.X(), momentum.Y(), momentum.Z()};
203 //Create object and write it to file
204 AliAODCaloCluster *calo = new AliAODCaloCluster(iParticle,1,labels,momentum.E(), x, NULL, ttype, 0);
206 SetCaloClusterPID(pdg,calo) ;
207 if(fDebug > 3 && momentum.Pt() > 0.2)
208 printf("AliCaloTrackMCReader::FillCalorimeters() - EMCAL : Selected cluster %s E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
209 particle->GetName(),momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
210 fAODEMCAL->Add(calo);//reference the selected object to the list
214 //____________________________________________________________________________
215 Bool_t AliCaloTrackMCReader::FillInputEvent(const Int_t iEntry, const char * /*currentFileName*/){
216 //Fill the event counter and input lists that are needed, called by the analysis maker.
218 fEventNumber = iEntry;
219 //fCurrentFileName = TString(currentFileName);
221 //In case of analysis of events with jets, skip those with jet pt > 5 pt hard
222 if(fComparePtHardAndJetPt && GetStack()) {
223 if(!ComparePtHardAndJetPt()) return kFALSE ;
229 Int_t iParticle = 0 ;
230 Double_t charge = 0.;
231 Int_t nparticles = GetStack()->GetNtrack() ;
232 if(fOnlyGeneratorParticles) nparticles=GetStack()->GetNprimary();
233 for (iParticle = 0 ; iParticle < nparticles ; iParticle++) {
234 TParticle * particle = GetStack()->Particle(iParticle);
235 TLorentzVector momentum;
238 Int_t pdg = particle->GetPdgCode();
240 //Keep particles with a given status
241 if(KeepParticleWithStatus(particle->GetStatusCode()) && (particle->Pt() > 0) ){
243 //Skip bizarre particles, they crash when charge is calculated
244 // if(TMath::Abs(pdg) == 3124 || TMath::Abs(pdg) > 10000000) continue ;
246 charge = TDatabasePDG::Instance()->GetParticle(pdg)->Charge();
247 particle->Momentum(momentum);
248 //---------- Charged particles ----------------------
251 if(TMath::Abs(momentum.Eta())< fTrackMultEtaCut) fTrackMult++;
253 if(fFillCTS && (momentum.Pt() > fCTSPtMin)){
254 //Particles in CTS acceptance
256 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) continue;
258 if(TMath::Abs(pdg) == 11 && GetStack()->Particle(particle->GetFirstMother())->GetPdgCode()==22) continue ;
260 if(fDebug > 3 && momentum.Pt() > 0.2)
261 printf("AliCaloTrackMCReader::FillInputEvent() - CTS : Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
262 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
264 x[0] = particle->Vx(); x[1] = particle->Vy(); x[2] = particle->Vz();
265 p[0] = particle->Px(); p[1] = particle->Py(); p[2] = particle->Pz();
266 //Create object and write it to file
267 AliAODTrack *aodTrack = new AliAODTrack(0, iParticle, p, kTRUE, x, kFALSE,NULL, 0, 0,
270 kFALSE, // No fit performed
271 kFALSE, // No fit performed
272 AliAODTrack::kPrimary,
274 SetTrackChargeAndPID(pdg, aodTrack);
275 fAODCTS->Add(aodTrack);//reference the selected object to the list
277 //Keep some charged particles in calorimeters lists
278 if((fFillPHOS || fFillEMCAL) && KeepChargedParticles(pdg)) FillCalorimeters(iParticle, particle, momentum);
282 //-------------Neutral particles ----------------------
283 else if(charge == 0 && (fFillPHOS || fFillEMCAL)){
284 //Skip neutrinos or other neutral particles
285 //if(SkipNeutralParticles(pdg) || particle->IsPrimary()) continue ; // protection added (MG)
286 if(SkipNeutralParticles(pdg)) continue ;
287 //Fill particle/calocluster arrays
289 FillCalorimeters(iParticle, particle, momentum);
292 //Sometimes pi0 are stable for the generator, if needed decay it by hand
294 if(momentum.Pt() > fPHOSPtMin || momentum.Pt() > fEMCALPtMin){
295 TLorentzVector lvGamma1, lvGamma2 ;
296 //Double_t angle = 0;
299 MakePi0Decay(momentum,lvGamma1,lvGamma2);//,angle);
301 //Check if Pi0 is in the acceptance of the calorimeters, if aperture angle is small, keep it
302 TParticle * pPhoton1 = new TParticle(22,1,iParticle,0,0,0,lvGamma1.Px(),lvGamma1.Py(),
303 lvGamma1.Pz(),lvGamma1.E(),0,0,0,0);
304 TParticle * pPhoton2 = new TParticle(22,1,iParticle,0,0,0,lvGamma2.Px(),lvGamma2.Py(),
305 lvGamma2.Pz(),lvGamma2.E(),0,0,0,0);
306 //Fill particle/calocluster arrays
307 FillCalorimeters(iParticle,pPhoton1,lvGamma1);
308 FillCalorimeters(iParticle,pPhoton2,lvGamma2);
311 else FillCalorimeters(iParticle,particle, momentum); //Add the rest
314 } //particle with correct status
317 fIndex2ndPhoton = -1; //In case of overlapping studies, reset for each event
323 //________________________________________________________________
324 Bool_t AliCaloTrackMCReader::KeepParticleWithStatus(Int_t status) const {
325 //Check if status is equal to one of the list
326 //These particles will be used in analysis.
328 for(Int_t i= 0; i < fStatusArray->GetSize(); i++)
329 if(status == fStatusArray->At(i)) return kTRUE ;
338 //________________________________________________________________
339 Bool_t AliCaloTrackMCReader::KeepChargedParticles(Int_t pdg) const {
340 //Check if pdg is equal to one of the charged particles list
341 //These particles will be added to the calorimeters lists.
343 for(Int_t i= 0; i < fChargedParticlesArray->GetSize(); i++)
344 if(TMath::Abs(pdg) == fChargedParticlesArray->At(i)) return kTRUE ;
350 //________________________________________________________________
351 void AliCaloTrackMCReader::Print(const Option_t * opt) const
353 //Print some relevant parameters set for the analysis
357 AliCaloTrackReader::Print(opt);
359 printf("**** Print **** %s %s ****\n", GetName(), GetTitle() ) ;
361 printf("Decay Pi0? : %d\n", fDecayPi0) ;
362 printf("Check Overlap in Calo? : %d\n", fCheckOverlap) ;
363 printf("Keep all status? : %d\n", fKeepAllStatus) ;
365 if(!fKeepAllStatus) printf("Keep particles with status : ");
366 for(Int_t i= 0; i < fStatusArray->GetSize(); i++)
367 printf(" %d ; ", fStatusArray->At(i));
370 printf("Skip neutral particles in calo : ");
371 for(Int_t i= 0; i < fNeutralParticlesArray->GetSize(); i++)
372 printf(" %d ; ", fNeutralParticlesArray->At(i));
375 printf("Keep charged particles in calo : ");
376 for(Int_t i= 0; i < fChargedParticlesArray->GetSize(); i++)
377 printf(" %d ; ", fChargedParticlesArray->At(i));
382 //____________________________________________________________________________
383 void AliCaloTrackMCReader::MakePi0Decay(TLorentzVector &p0, TLorentzVector &p1,
384 TLorentzVector &p2) const {//, Double_t &angle) {
385 // Perform isotropic decay pi0 -> 2 photons
386 // p0 is pi0 4-momentum (inut)
387 // p1 and p2 are photon 4-momenta (output)
388 // cout<<"Boost vector"<<endl;
389 Double_t mPi0 = TDatabasePDG::Instance()->GetParticle(111)->Mass();
390 TVector3 b = p0.BoostVector();
391 //cout<<"Parameters"<<endl;
392 //Double_t mPi0 = p0.M();
393 Double_t phi = TMath::TwoPi() * gRandom->Rndm();
394 Double_t cosThe = 2 * gRandom->Rndm() - 1;
395 Double_t cosPhi = TMath::Cos(phi);
396 Double_t sinPhi = TMath::Sin(phi);
397 Double_t sinThe = TMath::Sqrt(1-cosThe*cosThe);
398 Double_t ePi0 = mPi0/2.;
399 //cout<<"ePi0 "<<ePi0<<endl;
400 //cout<<"Components"<<endl;
401 p1.SetPx(+ePi0*cosPhi*sinThe);
402 p1.SetPy(+ePi0*sinPhi*sinThe);
403 p1.SetPz(+ePi0*cosThe);
405 //cout<<"p1: "<<p1.Px()<<" "<<p1.Py()<<" "<<p1.Pz()<<" "<<p1.E()<<endl;
406 //cout<<"p1 Mass: "<<p1.Px()*p1.Px()+p1.Py()*p1.Py()+p1.Pz()*p1.Pz()-p1.E()*p1.E()<<endl;
407 p2.SetPx(-ePi0*cosPhi*sinThe);
408 p2.SetPy(-ePi0*sinPhi*sinThe);
409 p2.SetPz(-ePi0*cosThe);
411 //cout<<"p2: "<<p2.Px()<<" "<<p2.Py()<<" "<<p2.Pz()<<" "<<p2.E()<<endl;
412 //cout<<"p2 Mass: "<<p2.Px()*p2.Px()+p2.Py()*p2.Py()+p2.Pz()*p2.Pz()-p2.E()*p2.E()<<endl;
413 //cout<<"Boost "<<b.X()<<" "<<b.Y()<<" "<<b.Z()<<endl;
415 //cout<<"p1: "<<p1.Px()<<" "<<p1.Py()<<" "<<p1.Pz()<<" "<<p1.E()<<endl;
417 //cout<<"p2: "<<p2.Px()<<" "<<p2.Py()<<" "<<p2.Pz()<<" "<<p2.E()<<endl;
418 //cout<<"angle"<<endl;
419 //angle = p1.Angle(p2.Vect());
423 //____________________________________________________________________________
424 void AliCaloTrackMCReader::SetInputOutputMCEvent(AliVEvent* /*esd*/, AliAODEvent* aod, AliMCEvent* mc) {
425 // Connect the data pointer
431 //________________________________________________________________
432 Bool_t AliCaloTrackMCReader::SkipNeutralParticles(Int_t pdg) const {
433 //Check if pdg is equal to one of the neutral particles list
434 //These particles will be skipped from analysis.
436 for(Int_t i= 0; i < fNeutralParticlesArray->GetSize(); i++)
437 if(TMath::Abs(pdg) == fNeutralParticlesArray->At(i)) return kTRUE ;
444 //____________________________________________________________________
445 void AliCaloTrackMCReader::SetTrackChargeAndPID(const Int_t pdgCode, AliAODTrack *track) const {
446 //Give a PID weight for tracks equal to 1 depending on the particle type
448 Float_t pid[10] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0.};
454 pid[AliAODTrack::kUnknown] = 1.;
459 track->SetCharge(-1);
460 pid[AliAODTrack::kElectron] = 1.;
465 track->SetCharge(+1);
466 pid[AliAODTrack::kElectron] = 1.;
471 track->SetCharge(-1);
472 pid[AliAODTrack::kMuon] = 1.;
477 track->SetCharge(+1);
478 pid[AliAODTrack::kMuon] = 1.;
484 pid[AliAODTrack::kUnknown] = 1.;
489 track->SetCharge(+1);
490 pid[AliAODTrack::kPion] = 1.;
495 track->SetCharge(-1);
496 pid[AliAODTrack::kPion] = 1.;
502 pid[AliAODTrack::kUnknown] = 1.;
507 track->SetCharge(+1);
508 pid[AliAODTrack::kKaon] = 1.;
513 track->SetCharge(-1);
514 pid[AliAODTrack::kKaon] = 1.;
520 pid[AliAODTrack::kUnknown] = 1.;
525 track->SetCharge(+1);
526 pid[AliAODTrack::kProton] = 1.;
530 case -2212: // anti-p
531 track->SetCharge(-1);
532 pid[AliAODTrack::kProton] = 1.;
538 pid[AliAODTrack::kUnknown] = 1.;
544 pid[AliAODTrack::kUnknown] = 1.;
548 case -311: // anti-K0
550 pid[AliAODTrack::kUnknown] = 1.;
556 pid[AliAODTrack::kUnknown] = 1.;
562 pid[AliAODTrack::kUnknown] = 1.;
567 track->SetCharge(+1);
568 pid[AliAODTrack::kUnknown] = 1.;
573 track->SetCharge(-1);
574 pid[AliAODTrack::kUnknown] = 1.;
579 track->SetCharge(-1);
580 pid[AliAODTrack::kUnknown] = 1.;
586 pid[AliAODTrack::kUnknown] = 1.;
591 track->SetCharge(-1);
592 pid[AliAODTrack::kUnknown] = 1.;
597 track->SetCharge(-1);
598 pid[AliAODTrack::kUnknown] = 1.;
604 pid[AliAODTrack::kUnknown] = 1.;
608 case -3122: // anti-Lambda
610 pid[AliAODTrack::kUnknown] = 1.;
614 case -3222: // anti-Sigma-
615 track->SetCharge(-1);
616 pid[AliAODTrack::kUnknown] = 1.;
620 case -3212: // anti-Sigma0
622 pid[AliAODTrack::kUnknown] = 1.;
626 case -3112: // anti-Sigma+
627 track->SetCharge(+1);
628 pid[AliAODTrack::kUnknown] = 1.;
632 case -3322: // anti-Xi0
634 pid[AliAODTrack::kUnknown] = 1.;
638 case -3312: // anti-Xi+
639 track->SetCharge(+1);
642 case -3334: // anti-Omega+
643 track->SetCharge(+1);
644 pid[AliAODTrack::kUnknown] = 1.;
649 track->SetCharge(+1);
650 pid[AliAODTrack::kUnknown] = 1.;
655 track->SetCharge(-1);
656 pid[AliAODTrack::kUnknown] = 1.;
662 pid[AliAODTrack::kUnknown] = 1.;
666 case -421: // anti-D0
668 pid[AliAODTrack::kUnknown] = 1.;
673 track->SetCharge(-99);
674 pid[AliAODTrack::kUnknown] = 1.;
683 //____________________________________________________________________
684 void AliCaloTrackMCReader::SetCaloClusterPID(const Int_t pdgCode, AliVCluster *calo) const {
685 //Give a PID weight for CaloClusters equal to 1 depending on the particle type
687 Float_t pid[13] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.};
692 pid[AliVCluster::kPhoton] = 1.;
697 pid[AliVCluster::kElectron] = 1.;
702 pid[AliVCluster::kElectron] = 1.;
707 pid[AliVCluster::kCharged] = 1.;
712 pid[AliVCluster::kCharged] = 1.;
717 pid[AliVCluster::kPi0] = 1.;
722 pid[AliVCluster::kCharged] = 1.;
727 pid[AliVCluster::kCharged] = 1.;
732 pid[AliVCluster::kKaon0] = 1.;
733 pid[AliVCluster::kNeutral] = 1;
738 pid[AliVCluster::kCharged] = 1.;
743 pid[AliVCluster::kCharged] = 1.;
748 pid[AliVCluster::kNeutron] = 1.;
749 pid[AliVCluster::kNeutral] = 1.;
754 pid[AliVCluster::kCharged] = 1.;
758 case -2212: // anti-p
759 pid[AliVCluster::kCharged] = 1.;
764 pid[AliVCluster::kKaon0] = 1.;
765 pid[AliVCluster::kNeutral] = 1.;
770 pid[AliVCluster::kKaon0] = 1.;
771 pid[AliVCluster::kNeutral] = 1.;
775 case -311: // anti-K0
776 pid[AliVCluster::kKaon0] = 1.;
777 pid[AliVCluster::kNeutral] = 1.;
782 pid[AliVCluster::kNeutral] = 1.;
787 pid[AliVCluster::kUnknown] = 1.;
792 pid[AliVCluster::kUnknown] = 1.;
797 pid[AliVCluster::kUnknown] = 1.;
802 pid[AliVCluster::kUnknown] = 1.;
807 pid[AliVCluster::kUnknown] = 1.;
812 pid[AliVCluster::kUnknown] = 1.;
817 pid[AliVCluster::kUnknown] = 1.;
822 pid[AliVCluster::kNeutron] = 1.;
823 pid[AliVCluster::kNeutral] = 1.;
827 case -3122: // anti-Lambda
828 pid[AliVCluster::kUnknown] = 1.;
832 case -3222: // anti-Sigma-
833 pid[AliVCluster::kUnknown] = 1.;
837 case -3212: // anti-Sigma0
838 pid[AliVCluster::kUnknown] = 1.;
842 case -3112: // anti-Sigma+
843 pid[AliVCluster::kUnknown] = 1.;
847 case -3322: // anti-Xi0
848 pid[AliVCluster::kUnknown] = 1.;
852 case -3312: // anti-Xi+
853 pid[AliVCluster::kUnknown] = 1.;
857 case -3334: // anti-Omega+
858 pid[AliVCluster::kUnknown] = 1.;
863 pid[AliVCluster::kUnknown] = 1.;
868 pid[AliVCluster::kUnknown] = 1.;
873 pid[AliVCluster::kUnknown] = 1.;
877 case -421: // anti-D0
878 pid[AliVCluster::kUnknown] = 1.;
883 pid[AliVCluster::kUnknown] = 1.;