]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/AliAnalysisTaskAODvsESD.cxx
Generation of muon filtered AOD files from standard AOD files (Roberta)
[u/mrichter/AliRoot.git] / PWG3 / muon / AliAnalysisTaskAODvsESD.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TNtuple.h"
4 #include "TBranch.h"
5 #include "TH1F.h"
6 #include "TCanvas.h"
7 #include <TLorentzVector.h>
8 #include <TVector3.h>
9
10 #include "AliAnalysisTask.h"
11 #include "AliAnalysisManager.h"
12 #include "AliESDInputHandler.h"
13 #include "AliESDEvent.h"
14 #include "AliESDMuonTrack.h"
15 #include "AliAODInputHandler.h"
16 #include "AliAODHandler.h"
17 #include "AliAODEvent.h"
18 #include "AliAnalysisTaskAODvsESD.h"
19
20 ClassImp(AliAnalysisTaskAODvsESD)
21
22 //________________________________________________________________________
23 AliAnalysisTaskAODvsESD::AliAnalysisTaskAODvsESD(const char *name) 
24 : AliAnalysisTask(name, ""), fESD(0), fAOD(0), fList(0), fMuonNtuple(0), fMuonNtupleAOD(0), fInvMass(0), fInvMassAOD(0)
25 {
26   // Constructor
27   // Define input and output slots here
28   // Input slot #0 works with a TChain
29   DefineInput(0, TChain::Class());
30   // Output slot #0 writes NTuple/histos into a TList
31   DefineOutput(0, TList::Class());  
32 }
33
34 //________________________________________________________________________
35 void AliAnalysisTaskAODvsESD::ConnectInputData(Option_t *) 
36 {
37   // Connect ESD here
38   TTree* esdTree = dynamic_cast<TTree*> (GetInputData(0));
39   if (!esdTree) {
40     Printf("ERROR: Could not read chain from input slot 0");
41   } else {
42     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());   
43     if (!esdH) {
44       Printf("ERROR: Could not get ESDInputHandler");
45     } else
46       fESD = esdH->GetEvent();
47   }  
48
49   // Connect AOD here
50   AliAODHandler *aodH = dynamic_cast<AliAODHandler*> (AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler());  
51   if (!aodH) {
52     Printf("ERROR: Could not get AODHandler");
53   } else
54     fAOD = aodH->GetAOD();
55 }
56
57 //________________________________________________________________________
58 void AliAnalysisTaskAODvsESD::CreateOutputObjects()
59 {
60   // This method has to be called INSIDE the user redefined CreateOutputObjects
61   // method, before creating each object corresponding to the output containers
62   // that are to be written to a file. This need to be done in general for the big output
63   // objects that may not fit memory during processing. 
64   OpenFile(0); 
65
66   // Create Ntuples 
67   // For single muons: thetaX, thetaY, ptInv, eta, phi, theta, px, py, pz,ch
68   fMuonNtuple = new TNtuple("fMuonNtuple","Muon information","thX:thY:ptI:eta:phi:theta:px:py:pz");
69   fMuonNtupleAOD = new TNtuple("fMuonNtupleAOD","Muon information","etaAOD:phiAOD:thetaAOD:pxAOD:pyAOD:pzAOD");
70   // Create histos for inv mass
71   fInvMass = new TH1F("fInvMass","Inv. mass from ESDs",140,0,7);
72   fInvMassAOD = new TH1F("fInvMassAOD","Inv. mass from AOD",140,0,7);
73
74   // Add Ntuples to the list
75   fList = new TList();
76   fList->Add(fMuonNtuple);
77   fList->Add(fMuonNtupleAOD);
78   fList->Add(fInvMass);
79   fList->Add(fInvMassAOD);
80 }
81
82 //________________________________________________________________________
83 void AliAnalysisTaskAODvsESD::Exec(Option_t *) 
84 {
85   // Main loop, called for each event
86   if (!fESD) {
87     Printf("ERROR: fESD not available");
88     return;
89   }
90   if (!fAOD) {
91     Printf("ERROR: fAOD not available");
92     return;
93   }
94   
95   // ESDs
96   // for muons
97   Float_t muonMass = 0.105658369;
98   Float_t thetaX1=0;    Float_t thetaY1=0;
99   Float_t ptInv1=0;     Float_t pYZ1=0;
100   Float_t pxM1=0;       Float_t pyM1=0;        Float_t pzM1=0;
101   Float_t etaM1=0;      Float_t phiM1=0;       Float_t thetaM1=0;
102   Int_t chargeM1=0;      Float_t energyM1=0;    TLorentzVector lvM1;
103   Float_t thetaX2=0;    Float_t thetaY2=0;
104   Float_t ptInv2=0;     Float_t pYZ2=0;
105   Float_t pxM2=0;       Float_t pyM2=0;        Float_t pzM2=0;
106   Int_t chargeM2=0;      Float_t energyM2=0;    TLorentzVector lvM2;
107   //for J/psi
108   Float_t invM=0;
109   TLorentzVector lvJpsi;
110  
111   int nmt = fESD->GetNumberOfMuonTracks();
112   // first loop over muon tracks in event
113   for (Int_t mTrack = 0; mTrack < nmt; mTrack++) {
114     thetaX1=0; thetaY1=0; ptInv1=0; pYZ1=0; chargeM1=0;
115     pxM1=0;    pyM1=0;    pzM1=0;   energyM1=0;
116     AliESDMuonTrack* muonTrack1 = fESD->GetMuonTrack(mTrack);
117     thetaX1 = muonTrack1->GetThetaX();
118     thetaY1 = muonTrack1->GetThetaY();
119     ptInv1 = TMath::Abs(muonTrack1->GetInverseBendingMomentum());
120     pYZ1 = 1/ptInv1;
121     pzM1  = - pYZ1 / TMath::Sqrt(1.0 + TMath::Tan(thetaY1)*TMath::Tan(thetaY1));
122     pxM1  = pzM1 * TMath::Tan(thetaX1);
123     pyM1  = pzM1 * TMath::Tan(thetaY1);
124     energyM1 = TMath::Sqrt(muonMass*muonMass + pxM1*pxM1 + pyM1*pyM1 + pzM1*pzM1);
125     lvM1.SetPxPyPzE(pxM1,pyM1,pzM1,energyM1);
126     chargeM1 = Int_t(TMath::Sign(1.,muonTrack1->GetInverseBendingMomentum()));      
127     etaM1 = muonTrack1->Eta();
128     phiM1 = muonTrack1->Phi();
129     thetaM1 = muonTrack1->Theta();
130     fMuonNtuple->Fill(thetaX1,thetaY1,ptInv1,etaM1,phiM1,thetaM1,pxM1,pyM1,pzM1);
131     // second loop over muon tracks in event
132     for (Int_t mTrack2 = mTrack+1; mTrack2 < nmt; mTrack2++) {
133       thetaX2=0; thetaY2=0; ptInv2=0; pYZ2=0; chargeM2=0;
134       pxM2=0; pyM2=0; pzM2=0; energyM2=0;
135       AliESDMuonTrack* muonTrack2 = fESD->GetMuonTrack(mTrack2);
136       thetaX2 = muonTrack2->GetThetaX();
137       thetaY2 = muonTrack2->GetThetaY();
138       ptInv2 = TMath::Abs(muonTrack2->GetInverseBendingMomentum());
139       pYZ2 = 1/ptInv2;
140       pzM2  = - pYZ2 / TMath::Sqrt(1.0 + TMath::Tan(thetaY2)*TMath::Tan(thetaY2));
141       pxM2  = pzM2 * TMath::Tan(thetaX2);
142       pyM2  = pzM2 * TMath::Tan(thetaY2);
143       chargeM2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum()));      
144       energyM2 = TMath::Sqrt(muonMass*muonMass + pxM2*pxM2 + pyM2*pyM2 + pzM2*pzM2);
145       // if muons have opposite charge
146       if(chargeM1*chargeM2 == -1){
147         lvM2.SetPxPyPzE(pxM2, pyM2, pzM2,energyM2);
148         lvJpsi = lvM1 + lvM2;
149         invM = lvJpsi.M();
150         fInvMass->Fill(invM);
151       } // end if muons with opposite charge
152     } // end second loop over muon tracks in event
153   } // end first loop over muon tracks in event
154
155
156   // Created AOD
157   Float_t pxAodM1=0;   Float_t pyAodM1=0;     Float_t pzAodM1=0;
158   Float_t etaAodM1=0;  Float_t phiAodM1=0;    Float_t thetaAodM1=0;
159   Int_t chargeAodM1=0;  Float_t energyAodM1=0; TLorentzVector lvAodM1;
160   Float_t pxAodM2=0;   Float_t pyAodM2=0;     Float_t pzAodM2=0;
161   Float_t etaAodM2=0;  Float_t phiAodM2=0;    Float_t thetaAodM2=0;
162   Int_t chargeAodM2=0;  Float_t energyAodM2=0; TLorentzVector lvAodM2;
163   //for J/psi
164   Float_t invMAOD=0;
165   TLorentzVector lvJpsiAOD;
166
167   int nmtAOD = fAOD->GetNumberOfTracks();
168   // first loop over tracks
169   for (Int_t mTrack = 0; mTrack < nmtAOD; mTrack++) {
170     pxM1=0; pyM1=0; pzM1=0; 
171     AliAODTrack* muonTrack1 = fAOD->GetTrack(mTrack);
172     if(muonTrack1->IsMuonTrack()){
173       etaAodM1 = muonTrack1->Eta();
174       phiAodM1 = muonTrack1->Phi();
175       thetaAodM1 = muonTrack1->Theta();
176       pxAodM1 = muonTrack1->Px();
177       pyAodM1 = muonTrack1->Py();
178       pzAodM1 = muonTrack1->Pz();
179       chargeAodM1 = muonTrack1->Charge();
180       energyAodM1 = TMath::Sqrt(muonMass*muonMass + pxAodM1*pxAodM1 + pyAodM1*pyAodM1 + pzAodM1*pzAodM1);
181       lvAodM1.SetPxPyPzE(pxAodM1, pyAodM1, pzAodM1,energyAodM1);
182       fMuonNtupleAOD->Fill(etaAodM1,phiAodM1,thetaAodM1,pxAodM1,pyAodM1,pzAodM1);
183     }
184     // second loop over tracks in event
185     for (Int_t mTrack2 = mTrack+1; mTrack2 < nmtAOD; mTrack2++) {
186       chargeAodM2=0;
187       pxAodM2=0; pyAodM2=0; pzAodM2=0; energyAodM2=0;
188       AliAODTrack* muonTrack2 = fAOD->GetTrack(mTrack2);
189       if(muonTrack2->IsMuonTrack()){
190         etaAodM2 = muonTrack2->Eta();
191         phiAodM2 = muonTrack2->Phi();
192         thetaAodM2 = muonTrack2->Theta();
193         pxAodM2 = muonTrack2->Px();
194         pyAodM2 = muonTrack2->Py();
195         pzAodM2 = muonTrack2->Pz();
196         energyAodM2 = TMath::Sqrt(muonMass*muonMass + pxAodM2*pxAodM2 + pyAodM2*pyAodM2 + pzAodM2*pzAodM2);
197         chargeAodM2 = muonTrack2->Charge();
198         // if muons of opposite charge
199         if(chargeAodM1*chargeAodM2 == -1){
200           lvAodM2.SetPxPyPzE(pxAodM2, pyAodM2, pzAodM2,energyAodM2);
201           lvJpsiAOD = lvAodM1 + lvAodM2;
202           invMAOD = lvJpsiAOD.M();
203           fInvMassAOD->Fill(invMAOD);
204         }// end if muons with opposite charge
205       }// end if muon track
206     }// end second loop over tracks in event  
207   }// end first loop over tracks in event
208   
209   // Post final data. Write histo list to a file with option "RECREATE"
210   PostData(0,fList);
211 }      
212
213 //________________________________________________________________________
214 void AliAnalysisTaskAODvsESD::Terminate(const Option_t*)
215 {
216   // check if major differences between the two Ntuples (more comparisons can be added)
217   int n1 = fMuonNtuple->GetEntries();
218   int n2 = fMuonNtupleAOD->GetEntries();
219  
220   if(n1!=n2){
221     printf("ERROR: Different number of entries in single muon Ntuples\n");
222     return;
223   }
224   else
225     printf("Same number of entries in single muon Ntuples\n");
226   
227 //   TCanvas* cv1 = new TCanvas("cvn1","cvn1",500,350);
228 //   cv1->cd(1);
229 //   fInvMass->SetMarkerStyle(29);
230 //   fInvMass->SetMarkerColor(4);
231 //   fInvMass->Draw();
232 //   fInvMassAOD->SetMarkerStyle(28);
233 //   fInvMassAOD->SetMarkerColor(2);
234 //   fInvMassAOD->Draw("same"); 
235 }
236