[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/ExRootResult.h"
|
---|
| 35 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 36 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
| 37 |
|
---|
| 38 | #include "TMath.h"
|
---|
| 39 | #include "TString.h"
|
---|
| 40 | #include "TFormula.h"
|
---|
| 41 | #include "TRandom3.h"
|
---|
| 42 | #include "TObjArray.h"
|
---|
| 43 | #include "TDatabasePDG.h"
|
---|
| 44 | #include "TLorentzVector.h"
|
---|
| 45 |
|
---|
| 46 | #include <algorithm>
|
---|
| 47 | #include <stdexcept>
|
---|
| 48 | #include <iostream>
|
---|
| 49 | #include <sstream>
|
---|
| 50 |
|
---|
| 51 | using namespace std;
|
---|
| 52 |
|
---|
[36b4099] | 53 | namespace {
|
---|
| 54 | // integer power (faster than TMath::Pow() + cast to integer)
|
---|
| 55 | int ipow(int base, int exp)
|
---|
| 56 | {
|
---|
| 57 | int result = 1;
|
---|
| 58 | while (exp)
|
---|
| 59 | {
|
---|
| 60 | if (exp & 1)
|
---|
| 61 | result *= base;
|
---|
| 62 | exp >>= 1;
|
---|
| 63 | base *= base;
|
---|
| 64 | }
|
---|
| 65 |
|
---|
| 66 | return result;
|
---|
| 67 | }
|
---|
| 68 |
|
---|
| 69 | // standalone function to extract the i-th digit from a number (counting from 0 = rightmost, etc..)
|
---|
| 70 | int digit(int val, int i)
|
---|
| 71 | {
|
---|
| 72 | int y = ipow(10,i);
|
---|
| 73 | int z = val/y;
|
---|
| 74 | int val2 = val / (y * 10);
|
---|
| 75 | return (z - val2*10);
|
---|
| 76 | }
|
---|
| 77 |
|
---|
| 78 | // return the first two digits if this is a "fundamental" particle
|
---|
| 79 | // ID = 100 is a special case (internal generator ID's are 81-100)
|
---|
| 80 | // also, 101 and 102 are now used (by HepPID) for geantinos
|
---|
| 81 | int fundamentalID(int pdgCode)
|
---|
| 82 | {
|
---|
| 83 | pdgCode = abs(pdgCode);
|
---|
| 84 | if( ( digit(pdgCode, 9) == 1 ) && ( digit(pdgCode, 8) == 0 ) ) { return 0; }
|
---|
| 85 | if( digit(pdgCode, 2) == 0 && digit(pdgCode, 3) == 0) {
|
---|
| 86 | return pdgCode%10000;
|
---|
| 87 | } else if( pdgCode <= 102 ) {
|
---|
| 88 | return pdgCode;
|
---|
| 89 | } else {
|
---|
| 90 | return 0;
|
---|
| 91 | }
|
---|
| 92 | }
|
---|
| 93 |
|
---|
| 94 | bool hasBottom(int pdgCode)
|
---|
| 95 | {
|
---|
| 96 | if ( (pdgCode/10000000) > 0)
|
---|
| 97 | return false;
|
---|
| 98 | if (pdgCode <= 100)
|
---|
| 99 | return false;
|
---|
| 100 | if (fundamentalID(pdgCode) <= 100 && fundamentalID(pdgCode) > 0)
|
---|
| 101 | return false;
|
---|
| 102 | if( digit(pdgCode, 3) == 5 || digit(pdgCode, 2) == 5 || digit(pdgCode, 1) == 5 )
|
---|
| 103 | return true;
|
---|
| 104 | return false;
|
---|
| 105 | }
|
---|
| 106 | }
|
---|
| 107 |
|
---|
| 108 |
|
---|
[d7d2da3] | 109 | //------------------------------------------------------------------------------
|
---|
| 110 |
|
---|
| 111 | StatusPidFilter::StatusPidFilter() :
|
---|
| 112 | fItInputArray(0)
|
---|
| 113 | {
|
---|
| 114 | }
|
---|
| 115 |
|
---|
| 116 | //------------------------------------------------------------------------------
|
---|
| 117 |
|
---|
| 118 | StatusPidFilter::~StatusPidFilter()
|
---|
| 119 | {
|
---|
| 120 | }
|
---|
| 121 |
|
---|
| 122 | //------------------------------------------------------------------------------
|
---|
| 123 |
|
---|
| 124 | void StatusPidFilter::Init()
|
---|
| 125 | {
|
---|
| 126 | // PT threshold
|
---|
| 127 |
|
---|
| 128 | fPTMin = GetDouble("PTMin", 0.5);
|
---|
| 129 |
|
---|
| 130 | // import input array
|
---|
| 131 | fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
|
---|
| 132 | fItInputArray = fInputArray->MakeIterator();
|
---|
| 133 |
|
---|
| 134 | // create output array
|
---|
| 135 |
|
---|
| 136 | fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
|
---|
| 137 | }
|
---|
| 138 |
|
---|
| 139 | //------------------------------------------------------------------------------
|
---|
| 140 |
|
---|
| 141 | void StatusPidFilter::Finish()
|
---|
| 142 | {
|
---|
| 143 | if(fItInputArray) delete fItInputArray;
|
---|
| 144 | }
|
---|
| 145 |
|
---|
| 146 | //------------------------------------------------------------------------------
|
---|
| 147 |
|
---|
| 148 | void StatusPidFilter::Process()
|
---|
| 149 | {
|
---|
| 150 | Candidate *candidate;
|
---|
| 151 | Int_t status, pdgCode;
|
---|
[41da326] | 152 | Bool_t pass;
|
---|
[d7d2da3] | 153 |
|
---|
| 154 | fItInputArray->Reset();
|
---|
| 155 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
| 156 | {
|
---|
| 157 | status = candidate->Status;
|
---|
| 158 | pdgCode = TMath::Abs(candidate->PID);
|
---|
| 159 |
|
---|
[41da326] | 160 | pass = kFALSE;
|
---|
[d7d2da3] | 161 |
|
---|
[868daac] | 162 | // hard scattering particles (first condition for Py6, second for Py8)
|
---|
[a0f2226] | 163 | if(status == 3) pass = kTRUE;
|
---|
| 164 | if(status > 20 && status < 30 ) pass = kTRUE;
|
---|
[41da326] | 165 |
|
---|
[d5104a4] | 166 | // electrons, muons, taus and neutrinos
|
---|
| 167 | if(pdgCode > 10 && pdgCode < 17) pass = kTRUE;
|
---|
[41da326] | 168 |
|
---|
| 169 | // heavy quarks
|
---|
[868daac] | 170 | if(pdgCode == 4 ||pdgCode == 5 || pdgCode == 6) pass = kTRUE;
|
---|
[41da326] | 171 |
|
---|
| 172 | // Gauge bosons and other fundamental bosons
|
---|
[d5104a4] | 173 | if(pdgCode > 22 && pdgCode < 43) pass = kTRUE;
|
---|
[41da326] | 174 |
|
---|
[36b4099] | 175 | // logic ported from HepPDF: http://lcgapp.cern.ch/project/simu/HepPDT/HepPDT.2.05.02/html/ParticleID_8cc-source.html#l00081
|
---|
| 176 | bool is_b_hadron = hasBottom(pdgCode);
|
---|
[9c980d2] | 177 | bool is_b_quark = (pdgCode == 5);
|
---|
[36b4099] | 178 |
|
---|
| 179 | if (is_b_hadron)
|
---|
| 180 | pass = kTRUE;
|
---|
| 181 |
|
---|
[9c980d2] | 182 | // fPTMin not applied to b_hadrons / b_quarks to allow for b-enriched sample stitching
|
---|
| 183 | if(!pass || (candidate->Momentum.Pt() < fPTMin && !(is_b_hadron || is_b_quark)) ) continue;
|
---|
[d7d2da3] | 184 |
|
---|
| 185 | fOutputArray->Add(candidate);
|
---|
| 186 | }
|
---|
| 187 | }
|
---|
| 188 |
|
---|