]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Gamma finder added, all the cuts are optimized
authorsgorbuno <sgorbuno@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 27 Oct 2009 07:24:20 +0000 (07:24 +0000)
committersgorbuno <sgorbuno@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 27 Oct 2009 07:24:20 +0000 (07:24 +0000)
HLT/global/AliHLTVertexer.cxx
HLT/global/physics/AliHLTV0HistoComponent.cxx
HLT/global/physics/AliHLTV0HistoComponent.h

index ba0031b20622219f192e9e5b4cf4d8a8b7418ce6..96787a1e02810b5de2560595678a225c18e8f03f 100644 (file)
@@ -132,7 +132,6 @@ void AliHLTVertexer::FindV0s(  )
 
   for( Int_t iTr = 0; iTr<nTracks; iTr++ ){ //* first daughter
 
-    if( fESD->GetTrack(iTr)->GetTPCNcls()<60  ) continue;
     AliESDTrackInfo &info = fTrackInfos[iTr];
     if( !info.fOK ) continue;    
     if( info.fParticle.GetQ() >0 ) continue;    
@@ -140,7 +139,6 @@ void AliHLTVertexer::FindV0s(  )
 
     for( Int_t jTr = 0; jTr<nTracks; jTr++ ){  //* second daughter
 
-      if( fESD->GetTrack(jTr)->GetTPCNcls()<60  ) continue;
       AliESDTrackInfo &jnfo = fTrackInfos[jTr];
       if( !jnfo.fOK ) continue;
       if( jnfo.fParticle.GetQ() < 0 ) continue;
@@ -148,11 +146,11 @@ void AliHLTVertexer::FindV0s(  )
    
       //* construct V0 mother
 
-      AliKFParticle V0( info.fParticle, jnfo.fParticle );     
+      AliKFParticle v0( info.fParticle, jnfo.fParticle );     
 
       //* check V0 Chi^2
       
-      if( V0.GetChi2()<0 || V0.GetChi2() > 9.*V0.GetNDF() ) continue;
+      if( v0.GetChi2()<0 || v0.GetChi2() > 9.*v0.GetNDF() ) continue;
 
       //* subtruct daughters from primary vertex 
 
@@ -166,44 +164,31 @@ void AliHLTVertexer::FindV0s(  )
        if( primVtxCopy.GetNContributors()<=2 ) continue;
        primVtxCopy -= jnfo.fParticle;
       }
-      //* Check V0 Chi^2 deviation from primary vertex 
+      //* Check v0 Chi^2 deviation from primary vertex 
 
-      if( V0.GetDeviationFromVertex( primVtxCopy ) >3. ) continue;
+      if( v0.GetDeviationFromVertex( primVtxCopy ) >3. ) continue;
 
       //* Add V0 to primary vertex to improve the primary vertex resolution
 
-      primVtxCopy += V0;      
+      primVtxCopy += v0;      
 
       //* Set production vertex for V0
 
-      V0.SetProductionVertex( primVtxCopy );
+      v0.SetProductionVertex( primVtxCopy );
 
       //* Check chi^2 for a case
 
-      if( V0.GetChi2()<0 || V0.GetChi2()> 9.*V0.GetNDF() ) continue;
-
-      // Abschtand in [cm]
-
-      double dx = V0.GetX()-primVtxCopy.GetX();
-      double dy = V0.GetY()-primVtxCopy.GetY();
-      double r = sqrt(dx*dx + dy*dy);
-      //if( r>30 ) continue;
-      if( r<.2 ) continue;
-
-     //* Get V0 decay length with estimated error
+      if( v0.GetChi2()<0 || v0.GetChi2()> 9.*v0.GetNDF() ) continue;
+      
+      //* Get V0 decay length with estimated error
 
       Double_t length, sigmaLength;
-      if( V0.GetDecayLength( length, sigmaLength ) ) continue;
+      if( v0.GetDecayLength( length, sigmaLength ) ) continue;
 
       //* Reject V0 if it decays too close[sigma] to the primary vertex
 
       if( length  <3.5*sigmaLength ) continue;
 
-      //* Get V0 invariant mass 
-
-      // Double_t mass, sigmaMass;
-      //if( V0.GetMass( mass, sigmaMass ) ) continue;   
-      
       //* add ESD v0 
       
       AliESDv0 v0ESD( *fESD->GetTrack( iTr ), iTr, *fESD->GetTrack( jTr ), jTr );  
index 32e34a3a8487549cd09ac5083eaa5671adef4f57..e75bc6fa0754828eb4b01cec7976b3c071e37f30 100644 (file)
@@ -48,9 +48,12 @@ ClassImp(AliHLTV0HistoComponent)
 
 AliHLTV0HistoComponent::AliHLTV0HistoComponent()
 :
+  fGamma(0),
   fKShort(0),
   fLambda(0),
   fAP(0),
+  fGammaXY(0),
+  fNGammas(0),
   fNKShorts(0),
   fNLambdas(0)
 {
@@ -108,19 +111,32 @@ int AliHLTV0HistoComponent::DoInit( int argc, const char** argv )
 {
   // init
   
-  fKShort = new TH1F("hKShort","HLT KShort inv mass",80,0.4,.6); 
-  fKShort->SetFillColor(kBlue);
+  fGamma = new TH1F("hGamma","HLT:  #gamma inv mass",50,-.05,.2); 
+  fGamma->SetFillColor(kGreen);
+  fGamma->SetStats(0);
+  fKShort = new TH1F("hKShort","HLT:  K_{s}^{0} inv mass",80,0.4,.6); 
+  fKShort->SetFillColor(kGreen);
   fKShort->SetStats(0);
-  fLambda = new TH1F("hLambda","HLT Lambda inv mass",50,1.0,1.4); 
-  fLambda->SetFillColor(kBlue);
+  fLambda = new TH1F("hLambda","HLT:  #Lambda^{0} inv mass",50,1.0,1.4); 
+  fLambda->SetFillColor(kGreen);
   fLambda->SetStats(0);
-  fAP = new TH2F("hAP","HLT Armenteros-Podolanski plot",60,-1.,1.,60,0.,0.3);
+  fAP = new TH2F("hAP","HLT Armenteros-Podolanski plot",60,-1.,1.,60,0.,0.3);
   fAP->SetMarkerStyle(8);
   fAP->SetMarkerSize(0.4);
   fAP->SetYTitle("p_{t}[GeV]");
   fAP->SetXTitle("(p^{+}_{L}-p^{-}_{L})/(p^{+}_{L}+p^{-}_{L})");
   fAP->SetStats(0);
   fAP->SetOption("COLZ");
+
+  fGammaXY = new TH2F("hGammaXY","HLT:  #gamma conversions",100,-100.,100.,100,-100.,100.);
+  fGammaXY->SetMarkerStyle(8);
+  fGammaXY->SetMarkerSize(0.4);
+  fGammaXY->SetYTitle("X[cm]");
+  fGammaXY->SetXTitle("Y[cm]");
+  fGammaXY->SetStats(0);
+  fGammaXY->SetOption("COLZ");
+
+  fNGammas = 0;
   fNKShorts = 0;
   fNLambdas = 0;
 
@@ -143,9 +159,11 @@ int AliHLTV0HistoComponent::DoInit( int argc, const char** argv )
 int AliHLTV0HistoComponent::DoDeinit()
 {
   // see header file for class documentation
+  delete fGamma;
   delete fKShort;
   delete fLambda;
   delete fAP;
+  delete fGammaXY;
   return 0;
 }
 
@@ -161,83 +179,138 @@ int AliHLTV0HistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
     event->GetStdContent();
     Int_t nV0 = event->GetNumberOfV0s();
     AliKFParticle::SetField( event->GetMagneticField() );
-            
+
+    const double kKsMass = 0.49767;
+    const double kLambdaMass = 1.11568;
+
     for (Int_t iv=0; iv<nV0; iv++) {
-      AliESDv0 *v0 = event->GetV0(iv);
        
-      AliESDtrack *t1=event->GetTrack(v0->GetPindex());
-      AliESDtrack *t2=event->GetTrack(v0->GetNindex());      
+      AliESDtrack *t1=event->GetTrack( event->GetV0(iv)->GetNindex());
+      AliESDtrack *t2=event->GetTrack( event->GetV0(iv)->GetPindex());      
 
-      AliKFParticle kf1( *t1->GetInnerParam(), 211 );
-      AliKFParticle kf2( *t2->GetInnerParam(), 211 );
+      AliKFParticle kf1( *t1->GetInnerParam(), 11 );
+      AliKFParticle kf2( *t2->GetInnerParam(), 11 );
 
       AliKFVertex primVtx( *event->GetPrimaryVertexTracks() );
       double dev1 = kf1.GetDeviationFromVertex( primVtx );
       double dev2 = kf2.GetDeviationFromVertex( primVtx );
       
-      AliKFParticle V0( kf1, kf2 );
-      primVtx+=V0;
-      V0.SetProductionVertex( primVtx );
+      AliKFParticle v0( kf1, kf2 );
+      primVtx+=v0;
+      v0.SetProductionVertex( primVtx );
 
-      double dx = V0.GetX()-primVtx.GetX();
-      double dy = V0.GetY()-primVtx.GetY();
-      double r = sqrt(dx*dx + dy*dy);
-      
-      AliKFParticle kf01 = kf1, kf02 = kf2;
-      kf01.SetProductionVertex(V0);
-      kf02.SetProductionVertex(V0);
-      kf01.TransportToProductionVertex();
-      kf02.TransportToProductionVertex();
+      // Gamma finder
+      bool isGamma = 0;
+      {
+       double mass, error;       
+       v0.GetMass(mass,error);
+       if( fGamma ) fGamma->Fill( mass );
+
+       if( TMath::Abs(mass)<0.05 ){
+         AliKFParticle gamma = v0;
+         gamma.SetMassConstraint(0);
+         if( fGammaXY ) fGammaXY->Fill(gamma.GetX(), gamma.GetY());
+         isGamma = 1;
+         fNGammas++;
+       }            
+      }
       
-      double px1=kf01.GetPx(), py1=kf01.GetPy(), pz1=kf01.GetPz();
-      double px2=kf02.GetPx(), py2=kf02.GetPy(), pz2=kf02.GetPz();
-      double px = px1+px2, py = py1+py2, pz = pz1+pz2;
-      double p = sqrt(px*px+py*py+pz*pz);
-      double l1 = (px*px1 + py*py1 + pz*pz1)/p;
-      double l2 = (px*px2 + py*py2 + pz*pz2)/p;
-      double pt = sqrt(px1*px1+py1*py1+pz1*pz1 - l1*l1);
-      double ap = (l1-l2)/(l1+l2);
+      if( isGamma ) continue;
+
+      double dx = v0.GetX()-primVtx.GetX();
+      double dy = v0.GetY()-primVtx.GetY();
+      double r = sqrt(dx*dx + dy*dy);
+
+      if( r>50. ) continue;
 
       // AP plot
-      //if( dev1>3.5 && dev2>3.5 ){
+
+      double ap=0;
+
+      {
+       AliKFParticle kf01 = kf1, kf02 = kf2;
+       kf01.SetProductionVertex(v0);
+       kf02.SetProductionVertex(v0);
+       kf01.TransportToProductionVertex();
+       kf02.TransportToProductionVertex();      
+       double px1=kf01.GetPx(), py1=kf01.GetPy(), pz1=kf01.GetPz();
+       double px2=kf02.GetPx(), py2=kf02.GetPy(), pz2=kf02.GetPz();
+       double px = px1+px2, py = py1+py2, pz = pz1+pz2;
+       double p = sqrt(px*px+py*py+pz*pz);
+       double l1 = (px*px1 + py*py1 + pz*pz1)/p;
+       double l2 = (px*px2 + py*py2 + pz*pz2)/p;
+       double pt = sqrt(px1*px1+py1*py1+pz1*pz1 - l1*l1);
+       ap = (l1-l2)/(l1+l2);
        if( fAP ) fAP->Fill( ap, pt );
-       //}
+      } 
 
-      // kShort finder 
-            
-      if( fKShort ) fKShort->Fill( V0.GetMass() );
-      fNKShorts++;
+      double length = v0.GetDecayLength();
+
+      // KShort finder
+      
+      bool isKs = 0;
+      
+      if( t1->GetTPCNcls()>=60 && t2->GetTPCNcls()>=60 && length>=1.5 ){
+       
+       AliKFParticle piN( *t1->GetInnerParam(), 211 ); 
+       AliKFParticle piP( *t2->GetInnerParam(), 211 ); 
+       
+       AliKFParticle Ks( piN, piP );
+       Ks.SetProductionVertex( primVtx );
+       
+       double mass, error;
+       Ks.GetMass( mass, error);
+       if( fKShort ) fKShort->Fill( mass );
+       
+       if( TMath::Abs( mass - kKsMass )<3.5*error ){  
+         isKs = 1;
+         fNKShorts++;
+       }
+      }
+      
+      if( isKs ) continue;
 
       // Lambda finder 
-      if( r>=2. && fabs( ap )>= 0.4 ){                    
-       fNLambdas++;
+     
+      if( t1->GetTPCNcls()>=60 && t2->GetTPCNcls()>=60 
+         && dev1>=3. && dev2>=3. 
+         && length>=4. 
+         ){
+       
        AliKFParticle kP, kpi;
-       if( ap<0 ){ 
-         kP = AliKFParticle( *t1->GetInnerParam(), 2212 );
-         kpi = AliKFParticle( *t2->GetInnerParam(), 211 );
-       } else {
+       if( ap<0 ){ 
          kP = AliKFParticle( *t2->GetInnerParam(), 2212 );
          kpi = AliKFParticle( *t1->GetInnerParam(), 211 );
+       } else {
+         kP = AliKFParticle( *t1->GetInnerParam(), 2212 );
+         kpi = AliKFParticle( *t2->GetInnerParam(), 211 );
        }
-       AliKFParticle Lambda = AliKFParticle( kP, kpi );
-       Lambda.SetProductionVertex( primVtx );
-       
-       if( fLambda ) fLambda->Fill( Lambda.GetMass() );
+
+       AliKFParticle lambda = AliKFParticle( kP, kpi );
+       lambda.SetProductionVertex( primVtx );  
+       double mass, error;
+       lambda.GetMass( mass, error);
+       if( fLambda ) fLambda->Fill( mass );
+
+       if( TMath::Abs( mass - kLambdaMass )<3.5*error ){  
+         fNLambdas++;
+       }       
       }
+
     }// V0's
   
-    if( fKShort ){
-      PushBack( (TObject*) fKShort, kAliHLTDataTypeHistogram,0);
-    }
-    if( fAP ){
-      PushBack( (TObject*) fAP, kAliHLTDataTypeHistogram,0);
-    }
-    if( fLambda ){
-      PushBack( (TObject*) fLambda, kAliHLTDataTypeHistogram, 0);
-    }
+    if( fGamma ) PushBack( (TObject*) fGamma, kAliHLTDataTypeHistogram,0);
+    
+    if( fKShort ) PushBack( (TObject*) fKShort, kAliHLTDataTypeHistogram,0);
+    
+    if( fLambda ) PushBack( (TObject*) fLambda, kAliHLTDataTypeHistogram, 0);
+    
+    if( fAP ) PushBack( (TObject*) fAP, kAliHLTDataTypeHistogram,0);    
+
+    if( fGammaXY ) PushBack( (TObject*) fGammaXY, kAliHLTDataTypeHistogram,0);
   }  
   
-  HLTInfo("V0Histo found %d K0*'s and %d Lambdas total", fNKShorts, fNLambdas );
+  HLTInfo("Found %d Gammas, %d KShorts and %d Lambdas total", fNGammas, fNKShorts, fNLambdas );
   
   return 0;
 }
index 71d0c89dcca6a0781023b67003204d985569f3eb..6c428b752e3b3ad56e44ca851c31a355ce9a8600 100644 (file)
@@ -68,11 +68,14 @@ private:
    */ 
   int Configure(const char* arguments);
   
-  TH1F *fKShort;; // KShort inv. mass
-  TH1F *fLambda; // Lambda inv. mass
-  TH2F *fAP; // Armenteros-Podolanski 
+  TH1F *fGamma;   // Gamma inv. mass
+  TH1F *fKShort;   // Ks inv. mass
+  TH1F *fLambda;   // Lambda inv. mass
+  TH2F *fAP;       // Armenteros-Podolanski 
+  TH2F *fGammaXY;  // XY distribution of gamma convertions
+  Int_t fNGammas; // n found total
   Int_t fNKShorts; // n found total
-  Int_t fNLambdas;  // n found total
+  Int_t fNLambdas; // n found total
 
   ClassDef(AliHLTV0HistoComponent, 0);