]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add invariant mass histograms to check from MC the mass from conversions, check neede...
authorgconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 23 Sep 2011 18:02:02 +0000 (18:02 +0000)
committergconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 23 Sep 2011 18:02:02 +0000 (18:02 +0000)
PWG4/PartCorrDep/AliAnaPi0.cxx
PWG4/PartCorrDep/AliAnaPi0.h

index dd26370b94b0942ef438510ed6f4efcff5a89e00..78d12b7b9bee65b46217067b63d4f0d0b7bf64f4 100755 (executable)
@@ -87,7 +87,8 @@ fhPrimEtaPhi(0x0),           fhPrimEtaAccPhi(0x0),         fhPrimPi0PtOrigin(0x0
 fhMCOrgMass(),               fhMCOrgAsym(),                fhMCOrgDeltaEta(),            fhMCOrgDeltaPhi(),
 fhMCPi0MassPtRec(),          fhMCPi0MassPtTrue(),          fhMCPi0PtTruePtRec(),         
 fhMCEtaMassPtRec(),          fhMCEtaMassPtTrue(),          fhMCEtaPtTruePtRec(),
-fhMCPi0PtOrigin(0x0),        fhMCEtaPtOrigin(0x0)
+fhMCPi0PtOrigin(0x0),        fhMCEtaPtOrigin(0x0),
+fhReMCFromConversion(0),     fhReMCFromNotConversion(0),   fhReMCFromMixConversion(0)
 {
 //Default Ctor
  InitParameters();
@@ -625,6 +626,25 @@ TList * AliAnaPi0::GetCreateOutputObjects()
   
   //Histograms filled only if MC data is requested     
   if(IsDataMC()){
+    
+    fhReMCFromConversion = new TH2F("hReMCFromConversion","Invariant mass of 2 clusters originated in conversions",
+                         nptbins,ptmin,ptmax,nmassbins,massmin,massmax);
+    fhReMCFromConversion->SetXTitle("p_{T} (GeV/c)");
+    fhReMCFromConversion->SetYTitle("m_{#gamma,#gamma} (GeV/c^{2})");
+    outputContainer->Add(fhReMCFromConversion) ;
+    
+    fhReMCFromNotConversion = new TH2F("hReMCNotFromConversion","Invariant mass of 2 clusters not originated in conversions",
+                                    nptbins,ptmin,ptmax,nmassbins,massmin,massmax);
+    fhReMCFromNotConversion->SetXTitle("p_{T} (GeV/c)");
+    fhReMCFromNotConversion->SetYTitle("m_{#gamma,#gamma} (GeV/c^{2})");
+    outputContainer->Add(fhReMCFromNotConversion) ;
+    
+    fhReMCFromMixConversion = new TH2F("hReMCFromMixConversion","Invariant mass of 2 clusters one from conversion and the other not",
+                                    nptbins,ptmin,ptmax,nmassbins,massmin,massmax);
+    fhReMCFromMixConversion->SetXTitle("p_{T} (GeV/c)");
+    fhReMCFromMixConversion->SetYTitle("m_{#gamma,#gamma} (GeV/c^{2})");
+    outputContainer->Add(fhReMCFromMixConversion) ;
+    
     //Pi0
     fhPrimPi0Pt     = new TH1F("hPrimPi0Pt","Primary pi0 pt, Y<1",nptbins,ptmin,ptmax) ;
     fhPrimPi0AccPt  = new TH1F("hPrimPi0AccPt","Primary pi0 pt with both photons in acceptance",nptbins,ptmin,ptmax) ;
@@ -1920,7 +1940,25 @@ void AliAnaPi0::MakeAnalysisFillHistograms()
         // MC data
         //---------
         //Do some MC checks on the origin of the pair, is there any common ancestor and if there is one, who?
-        if(IsDataMC()) FillMCVersusRecDataHistograms(p1->GetLabel(), p2->GetLabel(),p1->Pt(), p2->Pt(),ncell1, ncell2, m, pt, a,deta, dphi);    
+        if(IsDataMC()) {
+          
+          if(GetMCAnalysisUtils()->CheckTagBit(p1->GetTag(),AliMCAnalysisUtils::kMCConversion) && 
+             GetMCAnalysisUtils()->CheckTagBit(p2->GetTag(),AliMCAnalysisUtils::kMCConversion))
+          {
+            fhReMCFromConversion->Fill(pt,m);
+          }
+          else if(!GetMCAnalysisUtils()->CheckTagBit(p1->GetTag(),AliMCAnalysisUtils::kMCConversion) && 
+                  !GetMCAnalysisUtils()->CheckTagBit(p2->GetTag(),AliMCAnalysisUtils::kMCConversion))
+          {
+            fhReMCFromNotConversion->Fill(pt,m);
+          }
+          else
+          {
+            fhReMCFromMixConversion->Fill(pt,m);
+          }
+                  
+          FillMCVersusRecDataHistograms(p1->GetLabel(), p2->GetLabel(),p1->Pt(), p2->Pt(),ncell1, ncell2, m, pt, a,deta, dphi); 
+        }
         
         //-----------------------
         //Multi cuts analysis 
index 1814d4a16c1aa5a101c22be429f4c758e20023c8..33ac1f476d47904ed45d9d637781afcff7f90d3b 100755 (executable)
@@ -228,7 +228,7 @@ class AliAnaPi0 : public AliAnaPartCorrBaseClass {
   //Multiple cuts: Assymmetry, pt, n cells, PID
   TH2F **  fhRePtNCellAsymCuts ;       //![fNPtCuts*fNAsymCuts*fNCellNCuts*] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
   TH2F **  fhMiPtNCellAsymCuts ;       //![fNPtCuts*fNAsymCuts*fNCellNCuts] Mixed two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
-  TH2F ** fhRePtNCellAsymCutsSM[12] ;  //![fNPtCuts*fNAsymCuts*fNCellNCutsfNModules] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry for each module
+  TH2F **  fhRePtNCellAsymCutsSM[12] ; //![fNPtCuts*fNAsymCuts*fNCellNCutsfNModules] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry for each module
  
   TH2F **  fhRePIDBits ;               //![fNPIDBits]  REAL two-photon invariant mass distribution for different PID bits
   TH3F **  fhRePtMult ;                //![fNAsymCuts] REAL two-photon invariant mass distribution for different track multiplicity and assymetry cuts
@@ -294,7 +294,13 @@ class AliAnaPi0 : public AliAnaPartCorrBaseClass {
   TH2F *   fhMCPi0PtOrigin ;           //! Mass of reoconstructed pi0 pairs  in calorimeter vs mother
   TH2F *   fhMCEtaPtOrigin ;           //! Mass of reoconstructed pi0 pairs  in calorimeter vs mother
 
-  ClassDef(AliAnaPi0,20)
+  TH2F *   fhReMCFromConversion ;      //! Invariant mass of 2 clusters originated in conversions
+  TH2F *   fhReMCFromNotConversion ;   //! Invariant mass of 2 clusters not originated in conversions
+  TH2F *   fhReMCFromMixConversion ;   //! Invariant mass of 2 clusters one from conversion and the other not
+
+
+  
+  ClassDef(AliAnaPi0,21)
 } ;