]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliAnalysisTaskDiJets.cxx
- AliPhysicsSelection: protected writing of fHistBunchCrossing and fHistTriggerPattern
[u/mrichter/AliRoot.git] / JETAN / AliAnalysisTaskDiJets.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 #include "AliAnalysisTaskDiJets.h"
19 #include "AliAODEvent.h"
20 #include "AliAODJet.h"
21 #include "AliAODDiJet.h"
22 #include <TClonesArray.h>
23 #include <TLorentzVector.h>
24 #include <TStyle.h>
25
26
27 ClassImp(AliAnalysisTaskDiJets)
28
29 ////////////////////////////////////////////////////////////////////////
30
31 AliAnalysisTaskDiJets::AliAnalysisTaskDiJets():
32     AliAnalysisTaskSE(),
33     fDiJets(0),
34     fDiJetsIn(0),
35     fFillAOD(kFALSE),
36     fHistList(0),
37     fH1DeltaPt(0),
38     fH1DeltaPhi(0),
39     fH1PhiImbal(0),
40     fH1Asym(0),
41     fH2Pt2vsPt1(0),
42     fH2DifvsSum(0)
43 {
44   // Default constructor
45 }
46
47 //----------------------------------------------------------------------
48 AliAnalysisTaskDiJets::AliAnalysisTaskDiJets(const char* name):
49     AliAnalysisTaskSE(name),
50     fDiJets(0),
51     fDiJetsIn(0),
52     fFillAOD(kFALSE),
53     fHistList(0),
54     fH1DeltaPt(0),
55     fH1DeltaPhi(0),
56     fH1PhiImbal(0),
57     fH1Asym(0),
58     fH2Pt2vsPt1(0),
59     fH2DifvsSum(0)
60 {
61   // Default constructor
62     DefineOutput(1, TList::Class());  
63 }
64
65 //----------------------------------------------------------------------
66 void AliAnalysisTaskDiJets::UserCreateOutputObjects()
67 {
68 // Create the output container
69 //
70     if (fDebug > 1) printf("AnalysisTaskDiJets::CreateOutPutData() \n");
71     fDiJets = new TClonesArray("AliAODDiJet", 0);
72     if (fFillAOD){
73       fDiJets->SetName("dijets");
74       AddAODBranch("TClonesArray", &fDiJets);
75         }
76
77     if (!fHistList) fHistList = new TList();
78     Float_t pi=TMath::Pi();
79     gStyle->SetPalette(1);
80
81     fH1DeltaPt  = new TH1F("DeltaPt","Difference between the jets' Pt;#Deltap_{T} (GeV/c);Entries",150,0.,150.);
82     fH1DeltaPt->SetMarkerSize(0.6);
83     fH1DeltaPt->SetMarkerColor(4);
84     fH1DeltaPt->SetMarkerStyle(21);
85     fH1DeltaPt->SetOption("E");
86
87     fH1DeltaPhi = new TH1F("DeltaPhi","Difference in the azimuthal angle;#Delta#phi;Entries",100,0.,pi);
88     fH1DeltaPhi->SetMarkerSize(0.6);
89     fH1DeltaPhi->SetMarkerColor(4);
90     fH1DeltaPhi->SetMarkerStyle(21);
91     fH1DeltaPhi->SetOption("E");
92
93     fH1PhiImbal = new TH1F("PhiImb","Phi imbalance;#phi;Entries",100,-pi,pi);
94     fH1PhiImbal->SetMarkerSize(0.6);
95     fH1PhiImbal->SetMarkerColor(4);
96     fH1PhiImbal->SetMarkerStyle(21);
97     fH1PhiImbal->SetOption("E");
98
99     fH1Asym     = new TH1F("Asym","Pt asymmetry;#Deltap_{T}/(p_{T,1}+p_{T,2});Entries",50,0.,1.);
100     fH1Asym->SetMarkerSize(0.6);
101     fH1Asym->SetMarkerColor(4);
102     fH1Asym->SetMarkerStyle(21);
103     fH1Asym->SetOption("E");
104
105     fH2Pt2vsPt1 = new TH2F("Pt2vsPt1","Pt2 versus Pt1;p_{T,1} (GeV/c);p_{T,2} (GeV/c)",250,0.,250.,250,0.,250.);
106     fH2Pt2vsPt1->SetOption("cont0");
107
108     fH2DifvsSum = new TH2F("DifvsSum","Pt difference vs Pt sum;p_{T,1}+p_{T,2} (GeV/c);#Deltap_{T} (GeV/c)",400,0.,400.,150,0.,150.);
109     fH2DifvsSum->SetOption("cont0");
110
111     fHistList->Add(fH1DeltaPt);
112     fHistList->Add(fH1DeltaPhi);
113     fHistList->Add(fH1PhiImbal);
114     fHistList->Add(fH1Asym);
115     fHistList->Add(fH2Pt2vsPt1);
116     fHistList->Add(fH2DifvsSum);
117 }
118
119 //----------------------------------------------------------------------
120 void AliAnalysisTaskDiJets::Init()
121 {
122     // Initialization
123     if (fDebug > 1) printf("AnalysisTaskDiJets::Init() \n");
124 }
125
126 //----------------------------------------------------------------------
127 void AliAnalysisTaskDiJets::UserExec(Option_t */*option*/)
128 {
129 // Execute analysis for current event
130 //
131     fDiJets->Delete();
132     AliAODEvent* aod   = dynamic_cast<AliAODEvent*> (InputEvent());
133
134     if(!aod){
135       // We do not have an input AOD, look in the output
136       aod = AODEvent();
137       if(!aod){
138         if (fDebug >1) printf("%s:%d AODEvent not found in the Output",(char*)__FILE__,__LINE__);
139         return;
140       }
141     }
142
143     TClonesArray* jets = aod->GetJets();
144
145     // N.B. if we take the aod from the output this is always
146     // empty and since it is the same as fDiJets 
147     fDiJetsIn = (TClonesArray*) (aod->GetList()->FindObject("dijets"));
148
149     if (fDiJetsIn) {
150       if (fDebug >1) printf("Found %d dijets in old list \n", fDiJetsIn->GetEntries());
151       AliAODJet* jj1, *jj2;
152       AliAODDiJet* testJ;
153
154       if (fDiJetsIn->GetEntries() > 0) {
155         testJ = (AliAODDiJet*) (fDiJetsIn->At(0));
156         jj1 = testJ->Jet(0);
157         jj1->Print("");
158         jj2 = testJ->Jet(1);
159         jj2->Print("");
160       }
161     }
162
163     Int_t nj = jets->GetEntriesFast();
164     if (fDebug >1) printf("There are %5d jets in the event \n", nj);
165
166     if (nj < 2){
167       PostData(1, fHistList);
168       return;
169     }
170     AliAODJet* jet1 = (AliAODJet*) (jets->At(0));
171     TLorentzVector v1 = *(jet1->MomentumVector());
172     AliAODJet* jet2 = (AliAODJet*) (jets->At(1));
173     TLorentzVector v2 = *(jet2->MomentumVector());
174     TLorentzVector v = v1 + v2;
175     Int_t ndi = fDiJets->GetEntriesFast();
176     TClonesArray &lref = *fDiJets;
177     new(lref[ndi]) AliAODDiJet(v);
178     AliAODDiJet* dijet = (AliAODDiJet*) (fDiJets->At(ndi));
179     dijet->SetJetRefs(jet1, jet2);
180
181     fH1DeltaPt->Fill(jet1->Pt()-jet2->Pt());
182     fH1DeltaPhi->Fill(dijet->DeltaPhi());
183     fH1PhiImbal->Fill(dijet->PhiImbalance());
184     fH1Asym->Fill((jet1->Pt()-jet2->Pt())/(jet1->Pt()+jet2->Pt()));
185     fH2Pt2vsPt1->Fill(jet1->Pt(),jet2->Pt());
186     fH2DifvsSum->Fill(jet1->Pt()+jet2->Pt(),jet1->Pt()-jet2->Pt());
187
188     PostData(1, fHistList);
189     return;
190 }
191
192 //----------------------------------------------------------------------
193 void AliAnalysisTaskDiJets::Terminate(Option_t */*option*/)
194 {
195 // Terminate analysis
196 //
197     if (fDebug > 1) printf("AnalysisDiJets: Terminate() \n");
198 }
199