1 | /*
|
---|
2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
4 | *
|
---|
5 | * This program is free software: you can redistribute it and/or modify
|
---|
6 | * it under the terms of the GNU General Public License as published by
|
---|
7 | * the Free Software Foundation, either version 3 of the License, or
|
---|
8 | * (at your option) any later version.
|
---|
9 | *
|
---|
10 | * This program is distributed in the hope that it will be useful,
|
---|
11 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
12 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
13 | * GNU General Public License for more details.
|
---|
14 | *
|
---|
15 | * You should have received a copy of the GNU General Public License
|
---|
16 | * along with this program. If not, see <http://www.gnu.org/licenses/>.
|
---|
17 | */
|
---|
18 |
|
---|
19 |
|
---|
20 | /** \class PhotonID
|
---|
21 | *
|
---|
22 | * Applies complex photon Id. Reconstructed photon candidtes are first separated into matched and non-matched to gen particles.
|
---|
23 | * Non-matched pass the "fake" efficiency. Matched photons get further splitted into isolated and non-isolated (user can choose criterion for isolation)
|
---|
24 | * Isolated photons pass the "prompt" efficiency while the non-isolated pass the "non-prompt" efficiency
|
---|
25 | *
|
---|
26 | * \author M. Selvaggi CERN
|
---|
27 | *
|
---|
28 | */
|
---|
29 |
|
---|
30 | #include "modules/PhotonID.h"
|
---|
31 |
|
---|
32 | #include "classes/DelphesClasses.h"
|
---|
33 | #include "classes/DelphesFactory.h"
|
---|
34 | #include "classes/DelphesFormula.h"
|
---|
35 |
|
---|
36 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
37 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
38 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
39 |
|
---|
40 | #include "TMath.h"
|
---|
41 | #include "TString.h"
|
---|
42 | #include "TFormula.h"
|
---|
43 | #include "TRandom3.h"
|
---|
44 | #include "TObjArray.h"
|
---|
45 | #include "TDatabasePDG.h"
|
---|
46 | #include "TLorentzVector.h"
|
---|
47 |
|
---|
48 | #include <algorithm>
|
---|
49 | #include <stdexcept>
|
---|
50 | #include <iostream>
|
---|
51 | #include <sstream>
|
---|
52 |
|
---|
53 | using namespace std;
|
---|
54 |
|
---|
55 | //------------------------------------------------------------------------------
|
---|
56 |
|
---|
57 | PhotonID::PhotonID() :
|
---|
58 | fPromptFormula(0), fNonPromptFormula(0), fFakeFormula(0), fItInputPhotonArray(0), fItInputGenArray(0)
|
---|
59 | {
|
---|
60 | fPromptFormula = new DelphesFormula;
|
---|
61 | fNonPromptFormula = new DelphesFormula;
|
---|
62 | fFakeFormula = new DelphesFormula;
|
---|
63 | }
|
---|
64 |
|
---|
65 | //------------------------------------------------------------------------------
|
---|
66 |
|
---|
67 | PhotonID::~PhotonID()
|
---|
68 | {
|
---|
69 | if(fPromptFormula) delete fPromptFormula;
|
---|
70 | if(fNonPromptFormula) delete fNonPromptFormula;
|
---|
71 | if(fFakeFormula) delete fFakeFormula;
|
---|
72 | }
|
---|
73 |
|
---|
74 | //------------------------------------------------------------------------------
|
---|
75 |
|
---|
76 | void PhotonID::Init()
|
---|
77 | {
|
---|
78 |
|
---|
79 | // read PhotonID formulae
|
---|
80 | fPromptFormula->Compile(GetString("PromptFormula", "1.0"));
|
---|
81 | fNonPromptFormula->Compile(GetString("NonPromptFormula", "1.0"));
|
---|
82 | fFakeFormula->Compile(GetString("FakeFormula", "1.0"));
|
---|
83 |
|
---|
84 | // import input arrays
|
---|
85 | fInputPhotonArray = ImportArray(GetString("InputPhotonArray", "PhotonIsolation/photons"));
|
---|
86 | fItInputPhotonArray = fInputPhotonArray->MakeIterator();
|
---|
87 |
|
---|
88 | // use filtered collection for speed
|
---|
89 | fInputGenArray = ImportArray(GetString("InputGenArray", "GenParticleFilter/filteredParticles"));
|
---|
90 | fItInputGenArray = fInputGenArray->MakeIterator();
|
---|
91 |
|
---|
92 | // min pt to be considered, make sure this threshold is higher than threshold in particle filter
|
---|
93 | fPTMin = GetDouble("PTMin", 10.0);
|
---|
94 |
|
---|
95 | // to be tuned, since FS and delphes have different isolation profiles
|
---|
96 | fRelIsoMax = GetDouble("fRelIsoMax", 0.3);
|
---|
97 |
|
---|
98 | // create output array
|
---|
99 | fOutputArray = ExportArray(GetString("OutputArray", "photons"));
|
---|
100 | }
|
---|
101 |
|
---|
102 | //------------------------------------------------------------------------------
|
---|
103 |
|
---|
104 | void PhotonID::Finish()
|
---|
105 | {
|
---|
106 | if(fItInputPhotonArray) delete fItInputPhotonArray;
|
---|
107 | if(fItInputGenArray) delete fItInputGenArray;
|
---|
108 | }
|
---|
109 |
|
---|
110 | //------------------------------------------------------------------------------
|
---|
111 |
|
---|
112 | void PhotonID::Process()
|
---|
113 | {
|
---|
114 | Candidate *candidate, *mother;
|
---|
115 | Double_t pt, eta, phi, e;
|
---|
116 | Double_t relIso;
|
---|
117 | Bool_t isolated;
|
---|
118 |
|
---|
119 | //cout<< "---- new event ---------"<<endl;
|
---|
120 |
|
---|
121 | fItInputPhotonArray->Reset();
|
---|
122 | while((candidate = static_cast<Candidate*>(fItInputPhotonArray->Next())))
|
---|
123 | {
|
---|
124 |
|
---|
125 | mother = candidate;
|
---|
126 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
127 | candidate->AddCandidate(mother);
|
---|
128 |
|
---|
129 | const TLorentzVector &candidatePosition = candidate->Position;
|
---|
130 | const TLorentzVector &candidateMomentum = candidate->Momentum;
|
---|
131 | eta = candidatePosition.Eta();
|
---|
132 | phi = candidatePosition.Phi();
|
---|
133 | pt = candidateMomentum.Pt();
|
---|
134 | e = candidateMomentum.E();
|
---|
135 |
|
---|
136 | if (pt < fPTMin) continue;
|
---|
137 |
|
---|
138 | //cout<< " ---- photon -----: "<<pt<<","<<eta<<","<<phi<<endl;
|
---|
139 |
|
---|
140 | // find out if photon matches does not match photon in gen collection and apply fae efficiency
|
---|
141 | if (isFake(candidate) )
|
---|
142 | {
|
---|
143 | //cout<<" Fake!"<<endl;
|
---|
144 |
|
---|
145 | if(gRandom->Uniform() > fFakeFormula->Eval(pt, eta, phi, e)) continue;
|
---|
146 | //cout<<" passed"<<endl;
|
---|
147 | candidate->Status = 3;
|
---|
148 | fOutputArray->Add(candidate);
|
---|
149 | }
|
---|
150 |
|
---|
151 | // if matches photon in gen collection
|
---|
152 | else
|
---|
153 | {
|
---|
154 | relIso = candidate->IsolationVar;
|
---|
155 | isolated = (relIso < 0.3);
|
---|
156 | //cout<<" Prompt!: "<<relIso<<endl;
|
---|
157 |
|
---|
158 | // if isolated apply prompt formula
|
---|
159 | if (isolated)
|
---|
160 | {
|
---|
161 | //cout<<" isolated!: "<<relIso<<endl;
|
---|
162 | if(gRandom->Uniform() > fPromptFormula->Eval(pt, eta, phi, e)) continue;
|
---|
163 | //cout<<" passed"<<endl;
|
---|
164 | candidate->Status = 1;
|
---|
165 | fOutputArray->Add(candidate);
|
---|
166 |
|
---|
167 | }
|
---|
168 |
|
---|
169 | // if non-isolated apply non-prompt formula
|
---|
170 | else
|
---|
171 | {
|
---|
172 | //cout<<" non-isolated!: "<<relIso<<endl;
|
---|
173 | if(gRandom->Uniform() > fNonPromptFormula->Eval(pt, eta, phi, e)) continue;
|
---|
174 | //cout<<" passed"<<endl;
|
---|
175 | candidate->Status = 2;
|
---|
176 | fOutputArray->Add(candidate);
|
---|
177 | }
|
---|
178 | }
|
---|
179 | }
|
---|
180 | }
|
---|
181 |
|
---|
182 | //------------------------------------------------------------------------------
|
---|
183 |
|
---|
184 | Bool_t PhotonID::isFake(const Candidate *obj){
|
---|
185 |
|
---|
186 | const TLorentzVector &mom_rec = obj->Momentum;
|
---|
187 |
|
---|
188 | Bool_t matches = false;
|
---|
189 | fItInputGenArray->Reset();
|
---|
190 | Candidate *gen;
|
---|
191 |
|
---|
192 | while((gen = static_cast<Candidate*>(fItInputGenArray->Next())))
|
---|
193 | {
|
---|
194 | const TLorentzVector &mom_gen = gen->Momentum;
|
---|
195 | Int_t status = gen->Status;
|
---|
196 | Int_t pdgCode = TMath::Abs(gen->PID);
|
---|
197 | Float_t dPtOverPt = TMath::Abs((mom_gen.Pt() - mom_rec.Pt())/mom_rec.Pt());
|
---|
198 | Float_t deltaR = mom_gen.DeltaR(mom_rec);
|
---|
199 |
|
---|
200 | if (status != 1) continue;
|
---|
201 | if (pdgCode != 22) continue;
|
---|
202 | if (dPtOverPt > 0.5) continue;
|
---|
203 | if (deltaR > 0.1) continue;
|
---|
204 |
|
---|
205 | matches = true;
|
---|
206 | break;
|
---|
207 |
|
---|
208 | }
|
---|
209 |
|
---|
210 | return !matches;
|
---|
211 | }
|
---|