[b443089] | 1 | /*
|
---|
| 2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
| 3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
[1fa50c2] | 4 | *
|
---|
[b443089] | 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.
|
---|
[1fa50c2] | 9 | *
|
---|
[b443089] | 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.
|
---|
[1fa50c2] | 14 | *
|
---|
[b443089] | 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 |
|
---|
[d7d2da3] | 19 | /** \class StatusPidFilter
|
---|
| 20 | *
|
---|
| 21 | * Removes all generated particles except electrons, muons, taus,
|
---|
| 22 | * and particles with status == 3.
|
---|
| 23 | *
|
---|
| 24 | * \author J. Hirschauer - FNAL
|
---|
| 25 | *
|
---|
| 26 | */
|
---|
| 27 |
|
---|
| 28 | #include "modules/StatusPidFilter.h"
|
---|
| 29 |
|
---|
| 30 | #include "classes/DelphesClasses.h"
|
---|
| 31 | #include "classes/DelphesFactory.h"
|
---|
| 32 | #include "classes/DelphesFormula.h"
|
---|
| 33 |
|
---|
| 34 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
[341014c] | 35 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 36 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
[d7d2da3] | 37 |
|
---|
| 38 | #include "TDatabasePDG.h"
|
---|
[341014c] | 39 | #include "TFormula.h"
|
---|
[d7d2da3] | 40 | #include "TLorentzVector.h"
|
---|
[341014c] | 41 | #include "TMath.h"
|
---|
| 42 | #include "TObjArray.h"
|
---|
| 43 | #include "TRandom3.h"
|
---|
| 44 | #include "TString.h"
|
---|
[d7d2da3] | 45 |
|
---|
| 46 | #include <algorithm>
|
---|
| 47 | #include <iostream>
|
---|
| 48 | #include <sstream>
|
---|
[341014c] | 49 | #include <stdexcept>
|
---|
[d7d2da3] | 50 |
|
---|
| 51 | using namespace std;
|
---|
| 52 |
|
---|
[341014c] | 53 | namespace
|
---|
| 54 | {
|
---|
| 55 | // integer power (faster than TMath::Pow() + cast to integer)
|
---|
| 56 | int ipow(int base, int exp)
|
---|
| 57 | {
|
---|
| 58 | int result = 1;
|
---|
| 59 | while(exp)
|
---|
[36b4099] | 60 | {
|
---|
[341014c] | 61 | if(exp & 1)
|
---|
| 62 | result *= base;
|
---|
| 63 | exp >>= 1;
|
---|
| 64 | base *= base;
|
---|
[36b4099] | 65 | }
|
---|
| 66 |
|
---|
[341014c] | 67 | return result;
|
---|
| 68 | }
|
---|
| 69 |
|
---|
| 70 | // standalone function to extract the i-th digit from a number (counting from 0 = rightmost, etc..)
|
---|
| 71 | int digit(int val, int i)
|
---|
| 72 | {
|
---|
| 73 | int y = ipow(10, i);
|
---|
| 74 | int z = val / y;
|
---|
| 75 | int val2 = val / (y * 10);
|
---|
| 76 | return (z - val2 * 10);
|
---|
| 77 | }
|
---|
| 78 |
|
---|
| 79 | // return the first two digits if this is a "fundamental" particle
|
---|
| 80 | // ID = 100 is a special case (internal generator ID's are 81-100)
|
---|
| 81 | // also, 101 and 102 are now used (by HepPID) for geantinos
|
---|
| 82 | int fundamentalID(int pdgCode)
|
---|
| 83 | {
|
---|
| 84 | pdgCode = abs(pdgCode);
|
---|
| 85 | if((digit(pdgCode, 9) == 1) && (digit(pdgCode, 8) == 0))
|
---|
[36b4099] | 86 | {
|
---|
[341014c] | 87 | return 0;
|
---|
[36b4099] | 88 | }
|
---|
[341014c] | 89 | if(digit(pdgCode, 2) == 0 && digit(pdgCode, 3) == 0)
|
---|
[36b4099] | 90 | {
|
---|
[341014c] | 91 | return pdgCode % 10000;
|
---|
[36b4099] | 92 | }
|
---|
[341014c] | 93 | else if(pdgCode <= 102)
|
---|
[36b4099] | 94 | {
|
---|
[341014c] | 95 | return pdgCode;
|
---|
[36b4099] | 96 | }
|
---|
[341014c] | 97 | else
|
---|
| 98 | {
|
---|
| 99 | return 0;
|
---|
| 100 | }
|
---|
| 101 | }
|
---|
[010117a] | 102 |
|
---|
[341014c] | 103 | bool hasBottom(int pdgCode)
|
---|
| 104 | {
|
---|
| 105 | if((pdgCode / 10000000) > 0)
|
---|
| 106 | return false;
|
---|
| 107 | if(pdgCode <= 100)
|
---|
| 108 | return false;
|
---|
| 109 | if(fundamentalID(pdgCode) <= 100 && fundamentalID(pdgCode) > 0)
|
---|
| 110 | return false;
|
---|
| 111 | if(digit(pdgCode, 3) == 5 || digit(pdgCode, 2) == 5 || digit(pdgCode, 1) == 5)
|
---|
| 112 | return true;
|
---|
| 113 | return false;
|
---|
| 114 | }
|
---|
[010117a] | 115 |
|
---|
[341014c] | 116 | bool isTauDaughter(int pdgCode, int M1, const TObjArray *fInputArray)
|
---|
| 117 | {
|
---|
| 118 | //not needed, just to speed up the code - can be further refined but gives only negligible improvement:
|
---|
| 119 | if(pdgCode == 15 || pdgCode < 11 || (pdgCode > 22 && pdgCode < 100) || pdgCode > 1000)
|
---|
| 120 | return false;
|
---|
[010117a] | 121 |
|
---|
[341014c] | 122 | if(M1 < 0)
|
---|
[010117a] | 123 | return false;
|
---|
| 124 |
|
---|
[341014c] | 125 | Candidate *mother;
|
---|
| 126 | mother = static_cast<Candidate *>(fInputArray->At(M1));
|
---|
| 127 | if(TMath::Abs(mother->PID) == 15)
|
---|
| 128 | return true;
|
---|
[27738e8] | 129 |
|
---|
[341014c] | 130 | return false;
|
---|
| 131 | }
|
---|
[27738e8] | 132 |
|
---|
[341014c] | 133 | bool isWDaughter(int M1, const TObjArray *fInputArray)
|
---|
| 134 | {
|
---|
| 135 | if(M1 < 0) return false;
|
---|
[27738e8] | 136 |
|
---|
[341014c] | 137 | Candidate *mother;
|
---|
| 138 | mother = static_cast<Candidate *>(fInputArray->At(M1));
|
---|
| 139 | if(TMath::Abs(mother->PID) == 24) return true;
|
---|
| 140 |
|
---|
| 141 | return false;
|
---|
[36b4099] | 142 | }
|
---|
| 143 |
|
---|
[341014c] | 144 | } // namespace
|
---|
[36b4099] | 145 |
|
---|
[d7d2da3] | 146 | //------------------------------------------------------------------------------
|
---|
| 147 |
|
---|
| 148 | StatusPidFilter::StatusPidFilter() :
|
---|
| 149 | fItInputArray(0)
|
---|
| 150 | {
|
---|
| 151 | }
|
---|
| 152 |
|
---|
| 153 | //------------------------------------------------------------------------------
|
---|
| 154 |
|
---|
| 155 | StatusPidFilter::~StatusPidFilter()
|
---|
| 156 | {
|
---|
| 157 | }
|
---|
| 158 |
|
---|
| 159 | //------------------------------------------------------------------------------
|
---|
| 160 |
|
---|
| 161 | void StatusPidFilter::Init()
|
---|
| 162 | {
|
---|
| 163 | // PT threshold
|
---|
| 164 | fPTMin = GetDouble("PTMin", 0.5);
|
---|
| 165 |
|
---|
[7a86a5c] | 166 | // keep or remove pileup particles
|
---|
| 167 | fRequireNotPileup = GetBool("RequireNotPileup", false);
|
---|
| 168 |
|
---|
[d7d2da3] | 169 | // import input array
|
---|
| 170 | fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
|
---|
| 171 | fItInputArray = fInputArray->MakeIterator();
|
---|
| 172 |
|
---|
| 173 | // create output array
|
---|
| 174 |
|
---|
| 175 | fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
|
---|
| 176 | }
|
---|
| 177 |
|
---|
| 178 | //------------------------------------------------------------------------------
|
---|
| 179 |
|
---|
| 180 | void StatusPidFilter::Finish()
|
---|
| 181 | {
|
---|
| 182 | if(fItInputArray) delete fItInputArray;
|
---|
| 183 | }
|
---|
| 184 |
|
---|
| 185 | //------------------------------------------------------------------------------
|
---|
| 186 |
|
---|
| 187 | void StatusPidFilter::Process()
|
---|
| 188 | {
|
---|
| 189 | Candidate *candidate;
|
---|
| 190 | Int_t status, pdgCode;
|
---|
[41da326] | 191 | Bool_t pass;
|
---|
[d7d2da3] | 192 |
|
---|
| 193 | fItInputArray->Reset();
|
---|
[341014c] | 194 | while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
|
---|
[d7d2da3] | 195 | {
|
---|
| 196 | status = candidate->Status;
|
---|
| 197 | pdgCode = TMath::Abs(candidate->PID);
|
---|
| 198 |
|
---|
[41da326] | 199 | pass = kFALSE;
|
---|
[d7d2da3] | 200 |
|
---|
[cbc56c5] | 201 | // Store all SUSY particles
|
---|
| 202 | if(pdgCode >= 1000001 && pdgCode <= 1000039) pass = kTRUE;
|
---|
| 203 |
|
---|
[868daac] | 204 | // hard scattering particles (first condition for Py6, second for Py8)
|
---|
[a0f2226] | 205 | if(status == 3) pass = kTRUE;
|
---|
[341014c] | 206 | if(status > 20 && status < 30) pass = kTRUE;
|
---|
[41da326] | 207 |
|
---|
[d5104a4] | 208 | // electrons, muons, taus and neutrinos
|
---|
| 209 | if(pdgCode > 10 && pdgCode < 17) pass = kTRUE;
|
---|
[41da326] | 210 |
|
---|
| 211 | // heavy quarks
|
---|
[341014c] | 212 | if(pdgCode == 4 || pdgCode == 5 || pdgCode == 6) pass = kTRUE;
|
---|
[41da326] | 213 |
|
---|
| 214 | // Gauge bosons and other fundamental bosons
|
---|
[d5104a4] | 215 | if(pdgCode > 22 && pdgCode < 43) pass = kTRUE;
|
---|
[41da326] | 216 |
|
---|
[341014c] | 217 | //Stable photons
|
---|
| 218 | if(pdgCode == 22 && status == 1) pass = kTRUE;
|
---|
[f3598d0] | 219 |
|
---|
[36b4099] | 220 | // logic ported from HepPDF: http://lcgapp.cern.ch/project/simu/HepPDT/HepPDT.2.05.02/html/ParticleID_8cc-source.html#l00081
|
---|
| 221 | bool is_b_hadron = hasBottom(pdgCode);
|
---|
[341014c] | 222 | bool is_b_quark = (pdgCode == 5);
|
---|
[36b4099] | 223 |
|
---|
[010117a] | 224 | bool is_tau_daughter = isTauDaughter(pdgCode, candidate->M1, fInputArray);
|
---|
| 225 |
|
---|
[341014c] | 226 | if(is_b_hadron)
|
---|
[36b4099] | 227 | pass = kTRUE;
|
---|
| 228 |
|
---|
[341014c] | 229 | if(is_tau_daughter)
|
---|
[010117a] | 230 | pass = kTRUE;
|
---|
| 231 |
|
---|
[27738e8] | 232 | bool is_W_daughter = isWDaughter(candidate->M1, fInputArray);
|
---|
[341014c] | 233 | if(is_W_daughter)
|
---|
[27738e8] | 234 | pass = kTRUE;
|
---|
| 235 |
|
---|
[9c980d2] | 236 | // fPTMin not applied to b_hadrons / b_quarks to allow for b-enriched sample stitching
|
---|
[010117a] | 237 | // fPTMin not applied to tau decay products to allow visible-tau four momentum determination
|
---|
[341014c] | 238 | if(!pass || (candidate->Momentum.Pt() < fPTMin && !(is_b_hadron || is_b_quark || is_tau_daughter || is_W_daughter))) continue;
|
---|
[d7d2da3] | 239 |
|
---|
[7a86a5c] | 240 | // not pileup particles
|
---|
[341014c] | 241 | if(fRequireNotPileup && (candidate->IsPU > 0)) continue;
|
---|
[7a86a5c] | 242 |
|
---|
[d7d2da3] | 243 | fOutputArray->Add(candidate);
|
---|
| 244 | }
|
---|
| 245 | }
|
---|