1 |
ceballos |
1.1 |
// $Id $
|
2 |
|
|
|
3 |
|
|
#include "MitPhysics/SelMods/interface/WBFExampleAnalysisMod.h"
|
4 |
|
|
#include <TH1D.h>
|
5 |
|
|
#include <TH2D.h>
|
6 |
|
|
#include <TParameter.h>
|
7 |
|
|
#include "MitAna/DataUtil/interface/Debug.h"
|
8 |
|
|
#include "MitAna/DataTree/interface/Names.h"
|
9 |
|
|
#include "MitPhysics/Init/interface/ModNames.h"
|
10 |
|
|
#include "MitAna/DataCont/interface/ObjArray.h"
|
11 |
|
|
#include "MitCommon/MathTools/interface/MathUtils.h"
|
12 |
|
|
#include "MitAna/DataTree/interface/ParticleCol.h"
|
13 |
|
|
#include "TFile.h"
|
14 |
|
|
#include "TTree.h"
|
15 |
|
|
|
16 |
|
|
using namespace mithep;
|
17 |
|
|
ClassImp(mithep::WBFExampleAnalysisMod)
|
18 |
|
|
|
19 |
|
|
//--------------------------------------------------------------------------------------------------
|
20 |
|
|
WBFExampleAnalysisMod::WBFExampleAnalysisMod(const char *name, const char *title) :
|
21 |
|
|
BaseMod(name,title),
|
22 |
|
|
fMetName("NoDefaultNameSet"),
|
23 |
|
|
fCleanJetsName("NoDefaultNameSet"),
|
24 |
|
|
fVertexName(ModNames::gkGoodVertexesName),
|
25 |
|
|
fMet(0),
|
26 |
|
|
fVertices(0),
|
27 |
|
|
fJetPtMax(30),
|
28 |
|
|
fJetPtMin(30),
|
29 |
|
|
fDeltaEtaMin(4),
|
30 |
ceballos |
1.5 |
fDiJetMassMin(500)
|
31 |
ceballos |
1.1 |
{
|
32 |
|
|
// Constructor.
|
33 |
|
|
}
|
34 |
|
|
|
35 |
|
|
//--------------------------------------------------------------------------------------------------
|
36 |
|
|
void WBFExampleAnalysisMod::Begin()
|
37 |
|
|
{
|
38 |
|
|
// Run startup code on the client machine. For this module, we dont do
|
39 |
|
|
// anything here.
|
40 |
|
|
}
|
41 |
|
|
|
42 |
|
|
//--------------------------------------------------------------------------------------------------
|
43 |
|
|
void WBFExampleAnalysisMod::SlaveBegin()
|
44 |
|
|
{
|
45 |
|
|
// Run startup code on the computer (slave) doing the actual analysis. Here,
|
46 |
|
|
// we typically initialize histograms and other analysis objects and request
|
47 |
|
|
// branches. For this module, we request a branch of the MitTree.
|
48 |
|
|
|
49 |
|
|
//Create your histograms here
|
50 |
|
|
|
51 |
|
|
//*************************************************************************************************
|
52 |
|
|
// Selection Histograms
|
53 |
|
|
//*************************************************************************************************
|
54 |
ceballos |
1.2 |
AddTH1(fWBFSelection,"fWBFSelection", ";Cut Number;Number of Events", 7, -1.5, 5.5);
|
55 |
ceballos |
1.4 |
AddTH1(fWBFDeltaPhi, "fWBFDeltaPhi", ";Delta Phi;Number of Events", 90,0.,180 );
|
56 |
ceballos |
1.1 |
|
57 |
|
|
//***********************************************************************************************
|
58 |
|
|
// N-1 Histograms
|
59 |
|
|
//***********************************************************************************************
|
60 |
|
|
//All events
|
61 |
ceballos |
1.2 |
AddTH1(fWBFPtJetMax_NMinusOne, "fWBFPtJetMax_NMinusOne", ";Pt Jet Max;Number of Events",300,0.,300. );
|
62 |
|
|
AddTH1(fWBFPtJetMin_NMinusOne, "fWBFPtJetMin_NMinusOne", ";Pt Jet Min;Number of Events",300,0.,300. );
|
63 |
|
|
AddTH1(fWBFEta12_NMinusOne, "fWBFEta12_NMinusOne", ";eta1*eta2;Number of Events", 2,-1.5,1.5 );
|
64 |
|
|
AddTH1(fWBFdeltaEta_NMinusOne, "fWBFdeltaEta_NMinusOne", ";Delta Eta;Number of Events", 100,0.,10. );
|
65 |
|
|
AddTH1(fWBFdijetMass_NMinusOne,"fWBFdijetMass_NMinusOne",";DiJet Mass;Number of Events",300,0.,3000.);
|
66 |
|
|
AddTH1(fWBFZVar_NMinusOne, "fWBFZVar_NMinusOne", ";Z variable;Number of Events",50,0.,5. );
|
67 |
ceballos |
1.1 |
|
68 |
|
|
//***********************************************************************************************
|
69 |
|
|
// After all cuts Histograms
|
70 |
|
|
//***********************************************************************************************
|
71 |
|
|
AddTH1(fWBFSSMass_afterCuts,"fWBFSSMass_afterCuts",";SS dilepton Mass;Number of Events",200,0.,400);
|
72 |
ceballos |
1.2 |
AddTH1(fWBFSSDeltaPhi_afterCuts,"fWBFSSDeltaPhi_afterCuts",";SS DeltaPhi 2l;Number of Events",90,0.,180);
|
73 |
|
|
AddTH1(fWBFOSMass_afterCuts,"fWBFOSMass_afterCuts",";OS dilepton Mass;Number of Events",200,0.,400);
|
74 |
|
|
AddTH1(fWBFOSDeltaPhi_afterCuts,"fWBFOSDeltaPhi_afterCuts",";OS DeltaPhi 2l;Number of Events",90,0.,180);
|
75 |
ceballos |
1.1 |
|
76 |
|
|
}
|
77 |
|
|
|
78 |
|
|
//--------------------------------------------------------------------------------------------------
|
79 |
|
|
void WBFExampleAnalysisMod::Process()
|
80 |
|
|
{
|
81 |
|
|
//Obtain all the good objects from the event cleaning module
|
82 |
|
|
fVertices = GetObjThisEvt<VertexOArr>(fVertexName);
|
83 |
|
|
ParticleOArr *CleanLeptons = dynamic_cast<mithep::ParticleOArr*>
|
84 |
|
|
(FindObjThisEvt(ModNames::gkMergedLeptonsName));
|
85 |
|
|
ObjArray<Jet> *CleanJets = dynamic_cast<ObjArray<Jet>* >
|
86 |
|
|
(FindObjThisEvt(fCleanJetsName.Data()));
|
87 |
|
|
TParameter<Double_t> *NNLOWeight = GetObjThisEvt<TParameter<Double_t> >("NNLOWeight");
|
88 |
|
|
|
89 |
|
|
MetCol *met = dynamic_cast<ObjArray<Met>* >(FindObjThisEvt(fMetName));
|
90 |
|
|
const Met *caloMet = 0;
|
91 |
|
|
if (met) {
|
92 |
|
|
caloMet = met->At(0);
|
93 |
|
|
} else {
|
94 |
|
|
cout << "Error: Met Collection " << fMetName << " could not be loaded.\n";
|
95 |
|
|
return;
|
96 |
|
|
}
|
97 |
|
|
|
98 |
|
|
if(CleanJets->GetEntries() < 2) return;
|
99 |
|
|
|
100 |
|
|
//*********************************************************************************************
|
101 |
|
|
//Define Cuts
|
102 |
|
|
//*********************************************************************************************
|
103 |
ceballos |
1.2 |
const int nCuts = 6;
|
104 |
|
|
bool passCut[nCuts] = {false, false, false, false, false, false};
|
105 |
ceballos |
1.1 |
|
106 |
|
|
// ptjet max cut
|
107 |
|
|
if(CleanJets->At(0)->Pt() > fJetPtMax) passCut[0] = true;
|
108 |
|
|
|
109 |
|
|
// ptjet min cut
|
110 |
|
|
if(CleanJets->At(1)->Pt() > fJetPtMin) passCut[1] = true;
|
111 |
|
|
|
112 |
|
|
// this cut is always applied
|
113 |
|
|
if(CleanJets->At(0)->Eta()*CleanJets->At(1)->Eta() < 0) passCut[2] = true;
|
114 |
|
|
|
115 |
|
|
// deltaEta cut
|
116 |
|
|
double deltaEta = TMath::Abs(CleanJets->At(0)->Eta()-CleanJets->At(1)->Eta());
|
117 |
|
|
if(deltaEta > fDeltaEtaMin) passCut[3] = true;
|
118 |
|
|
|
119 |
|
|
// dijet mass cut
|
120 |
|
|
CompositeParticle dijet;
|
121 |
|
|
dijet.AddDaughter(CleanJets->At(0));
|
122 |
|
|
dijet.AddDaughter(CleanJets->At(1));
|
123 |
|
|
if(dijet.Mass() > fDiJetMassMin) passCut[4] = true;
|
124 |
ceballos |
1.2 |
|
125 |
ceballos |
1.3 |
// jet veto cut, use of Zeppenfeld variable
|
126 |
ceballos |
1.2 |
passCut[5] = true;
|
127 |
|
|
double zVarMin = 30.;
|
128 |
|
|
if(CleanJets->GetEntries() >=2 ){
|
129 |
|
|
for(UInt_t i=2; i<CleanJets->GetEntries(); i++){
|
130 |
|
|
if(CleanJets->At(i)->Pt() <= 30.0) return;
|
131 |
|
|
double zVar = TMath::Abs(CleanJets->At(i)->Eta()-(CleanJets->At(0)->Eta()+CleanJets->At(1)->Eta())/2)/
|
132 |
|
|
TMath::Abs(CleanJets->At(0)->Eta()-CleanJets->At(1)->Eta());
|
133 |
|
|
if(zVar < zVarMin) zVarMin = zVar;
|
134 |
|
|
}
|
135 |
|
|
if(zVarMin < 1) passCut[5] = false;
|
136 |
|
|
}
|
137 |
ceballos |
1.4 |
|
138 |
ceballos |
1.1 |
//*********************************************************************************************
|
139 |
|
|
//Make Selection Histograms. Number of events passing each level of cut
|
140 |
|
|
//*********************************************************************************************
|
141 |
|
|
bool passAllCuts = true;
|
142 |
|
|
for(int c=0; c<nCuts; c++) passAllCuts = passAllCuts & passCut[c];
|
143 |
ceballos |
1.4 |
|
144 |
|
|
if(passAllCuts){
|
145 |
|
|
double deltaPhiJJ = MathUtils::DeltaPhi(CleanJets->At(0)->Phi(),
|
146 |
|
|
CleanJets->At(1)->Phi())* 180.0 / TMath::Pi();
|
147 |
|
|
fWBFDeltaPhi->Fill(deltaPhiJJ,NNLOWeight->GetVal());
|
148 |
|
|
}
|
149 |
|
|
|
150 |
ceballos |
1.1 |
//Cut Selection Histograms
|
151 |
|
|
fWBFSelection->Fill(-1,NNLOWeight->GetVal());
|
152 |
|
|
for (int k=0;k<nCuts;k++) {
|
153 |
|
|
bool pass = true;
|
154 |
|
|
bool passPreviousCut = true;
|
155 |
|
|
for (int p=0;p<=k;p++) {
|
156 |
|
|
pass = (pass && passCut[p]);
|
157 |
|
|
if (p<k)
|
158 |
|
|
passPreviousCut = (passPreviousCut&& passCut[p]);
|
159 |
|
|
}
|
160 |
|
|
|
161 |
|
|
if (pass) {
|
162 |
|
|
fWBFSelection->Fill(k,NNLOWeight->GetVal());
|
163 |
|
|
}
|
164 |
|
|
}
|
165 |
|
|
|
166 |
|
|
//*********************************************************************************************
|
167 |
|
|
// N-1 Histograms
|
168 |
|
|
//*********************************************************************************************
|
169 |
|
|
// N-1 ptjet max
|
170 |
|
|
bool pass = true;
|
171 |
|
|
for (int k=0;k<nCuts;k++) {
|
172 |
|
|
if (k != 0) {
|
173 |
|
|
pass = (pass && passCut[k]);
|
174 |
|
|
}
|
175 |
|
|
}
|
176 |
|
|
if (pass) {
|
177 |
ceballos |
1.2 |
fWBFPtJetMax_NMinusOne->Fill(TMath::Min(CleanJets->At(0)->Pt(),299.999),NNLOWeight->GetVal());
|
178 |
ceballos |
1.1 |
}
|
179 |
|
|
|
180 |
|
|
// N-1 ptjet min
|
181 |
|
|
pass = true;
|
182 |
|
|
for (int k=0;k<nCuts;k++) {
|
183 |
|
|
if (k != 1) {
|
184 |
|
|
pass = (pass && passCut[k]);
|
185 |
|
|
}
|
186 |
|
|
}
|
187 |
|
|
if (pass) {
|
188 |
ceballos |
1.2 |
fWBFPtJetMin_NMinusOne->Fill(TMath::Min(CleanJets->At(1)->Pt(),299.999),NNLOWeight->GetVal());
|
189 |
|
|
}
|
190 |
|
|
|
191 |
|
|
// N-1 eta1*eta2
|
192 |
|
|
pass = true;
|
193 |
|
|
for (int k=0;k<nCuts;k++) {
|
194 |
|
|
if (k != 2) {
|
195 |
|
|
pass = (pass && passCut[k]);
|
196 |
|
|
}
|
197 |
|
|
}
|
198 |
|
|
if (pass) {
|
199 |
|
|
fWBFEta12_NMinusOne->Fill(CleanJets->At(0)->Eta()*CleanJets->At(1)->Eta()/
|
200 |
|
|
TMath::Abs(CleanJets->At(0)->Eta()*CleanJets->At(1)->Eta()),NNLOWeight->GetVal());
|
201 |
ceballos |
1.1 |
}
|
202 |
|
|
|
203 |
|
|
// N-1 deltaEta
|
204 |
|
|
pass = true;
|
205 |
|
|
for (int k=0;k<nCuts;k++) {
|
206 |
|
|
if (k != 3) {
|
207 |
|
|
pass = (pass && passCut[k]);
|
208 |
|
|
}
|
209 |
|
|
}
|
210 |
|
|
if (pass) {
|
211 |
|
|
fWBFdeltaEta_NMinusOne->Fill(TMath::Min(deltaEta,9.999),NNLOWeight->GetVal());
|
212 |
|
|
}
|
213 |
|
|
|
214 |
|
|
// N-1 dijet mass
|
215 |
|
|
pass = true;
|
216 |
|
|
for (int k=0;k<nCuts;k++) {
|
217 |
|
|
if (k != 4) {
|
218 |
|
|
pass = (pass && passCut[k]);
|
219 |
|
|
}
|
220 |
|
|
}
|
221 |
|
|
if (pass) {
|
222 |
ceballos |
1.2 |
fWBFdijetMass_NMinusOne->Fill(TMath::Min(dijet.Mass(),2999.999),NNLOWeight->GetVal());
|
223 |
|
|
}
|
224 |
|
|
|
225 |
|
|
// N-1 dijet mass
|
226 |
|
|
pass = true;
|
227 |
|
|
for (int k=0;k<nCuts;k++) {
|
228 |
|
|
if (k != 5) {
|
229 |
|
|
pass = (pass && passCut[k]);
|
230 |
|
|
}
|
231 |
|
|
}
|
232 |
|
|
if (pass) {
|
233 |
|
|
fWBFZVar_NMinusOne->Fill(TMath::Min(zVarMin,4.999),NNLOWeight->GetVal());
|
234 |
ceballos |
1.1 |
}
|
235 |
|
|
|
236 |
|
|
//*********************************************************************************************
|
237 |
|
|
//Plots after all Cuts
|
238 |
|
|
//*********************************************************************************************
|
239 |
|
|
if (passAllCuts) {
|
240 |
|
|
// Distributions for dileptons events
|
241 |
|
|
if (CleanLeptons->GetEntries() >= 2 &&
|
242 |
|
|
CleanLeptons->At(0)->Pt() > 20 && CleanLeptons->At(1)->Pt() > 10){
|
243 |
|
|
|
244 |
|
|
CompositeParticle dilepton;
|
245 |
|
|
dilepton.AddDaughter(CleanLeptons->At(0));
|
246 |
|
|
dilepton.AddDaughter(CleanLeptons->At(1));
|
247 |
|
|
double deltaPhiLeptons = MathUtils::DeltaPhi(CleanLeptons->At(0)->Phi(),
|
248 |
|
|
CleanLeptons->At(1)->Phi())* 180.0 / TMath::Pi();
|
249 |
|
|
|
250 |
ceballos |
1.2 |
if(CleanLeptons->At(0)->Charge() * CleanLeptons->At(1)->Charge() > 0){ // same-sign
|
251 |
ceballos |
1.1 |
fWBFSSMass_afterCuts->Fill(TMath::Min(dilepton.Mass(),399.999),NNLOWeight->GetVal());
|
252 |
|
|
fWBFSSDeltaPhi_afterCuts->Fill(deltaPhiLeptons,NNLOWeight->GetVal());
|
253 |
|
|
}
|
254 |
ceballos |
1.2 |
else { // opposite-sign
|
255 |
ceballos |
1.1 |
fWBFOSMass_afterCuts->Fill(TMath::Min(dilepton.Mass(),399.999),NNLOWeight->GetVal());
|
256 |
|
|
fWBFOSDeltaPhi_afterCuts->Fill(deltaPhiLeptons,NNLOWeight->GetVal());
|
257 |
|
|
}
|
258 |
|
|
}
|
259 |
|
|
}
|
260 |
|
|
|
261 |
|
|
return;
|
262 |
|
|
}
|
263 |
|
|
|
264 |
|
|
//--------------------------------------------------------------------------------------------------
|
265 |
|
|
void WBFExampleAnalysisMod::SlaveTerminate()
|
266 |
|
|
{
|
267 |
|
|
|
268 |
|
|
// Run finishing code on the computer (slave) that did the analysis. For this
|
269 |
|
|
// module, we dont do anything here.
|
270 |
|
|
|
271 |
|
|
}
|
272 |
|
|
//--------------------------------------------------------------------------------------------------
|
273 |
|
|
void WBFExampleAnalysisMod::Terminate()
|
274 |
|
|
{
|
275 |
|
|
// Run finishing code on the client computer. For this module, we dont do
|
276 |
|
|
// anything here.
|
277 |
|
|
|
278 |
|
|
}
|