[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 |
|
---|
| 53 | //------------------------------------------------------------------------------
|
---|
| 54 |
|
---|
| 55 | StatusPidFilter::StatusPidFilter() :
|
---|
| 56 | fItInputArray(0)
|
---|
| 57 | {
|
---|
| 58 | }
|
---|
| 59 |
|
---|
| 60 | //------------------------------------------------------------------------------
|
---|
| 61 |
|
---|
| 62 | StatusPidFilter::~StatusPidFilter()
|
---|
| 63 | {
|
---|
| 64 | }
|
---|
| 65 |
|
---|
| 66 | //------------------------------------------------------------------------------
|
---|
| 67 |
|
---|
| 68 | void StatusPidFilter::Init()
|
---|
| 69 | {
|
---|
| 70 | // PT threshold
|
---|
| 71 |
|
---|
| 72 | fPTMin = GetDouble("PTMin", 0.5);
|
---|
| 73 |
|
---|
| 74 | // import input array
|
---|
| 75 | fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
|
---|
| 76 | fItInputArray = fInputArray->MakeIterator();
|
---|
| 77 |
|
---|
| 78 | // create output array
|
---|
| 79 |
|
---|
| 80 | fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
|
---|
| 81 | }
|
---|
| 82 |
|
---|
| 83 | //------------------------------------------------------------------------------
|
---|
| 84 |
|
---|
| 85 | void StatusPidFilter::Finish()
|
---|
| 86 | {
|
---|
| 87 | if(fItInputArray) delete fItInputArray;
|
---|
| 88 | }
|
---|
| 89 |
|
---|
| 90 | //------------------------------------------------------------------------------
|
---|
| 91 |
|
---|
| 92 | void StatusPidFilter::Process()
|
---|
| 93 | {
|
---|
| 94 | Candidate *candidate;
|
---|
| 95 | Int_t status, pdgCode;
|
---|
[41da326] | 96 | Bool_t pass;
|
---|
[d7d2da3] | 97 |
|
---|
| 98 | fItInputArray->Reset();
|
---|
| 99 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
| 100 | {
|
---|
| 101 | status = candidate->Status;
|
---|
| 102 | pdgCode = TMath::Abs(candidate->PID);
|
---|
| 103 |
|
---|
[41da326] | 104 | pass = kFALSE;
|
---|
[d7d2da3] | 105 |
|
---|
[41da326] | 106 | // status == 3
|
---|
| 107 | if(status == 3) pass = kTRUE;
|
---|
| 108 |
|
---|
[d5104a4] | 109 | // electrons, muons, taus and neutrinos
|
---|
| 110 | if(pdgCode > 10 && pdgCode < 17) pass = kTRUE;
|
---|
[41da326] | 111 |
|
---|
| 112 | // heavy quarks
|
---|
| 113 | if(pdgCode == 5 || pdgCode == 6) pass = kTRUE;
|
---|
| 114 |
|
---|
| 115 | // Gauge bosons and other fundamental bosons
|
---|
[d5104a4] | 116 | if(pdgCode > 22 && pdgCode < 43) pass = kTRUE;
|
---|
[41da326] | 117 |
|
---|
| 118 | if(!pass || candidate->Momentum.Pt() <= fPTMin) continue;
|
---|
[d7d2da3] | 119 |
|
---|
| 120 | fOutputArray->Add(candidate);
|
---|
| 121 | }
|
---|
| 122 | }
|
---|
| 123 |
|
---|