Fork me on GitHub

source: git/modules/StatusPidFilter.cc@ 400597a

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 400597a was 7c0fcd5, checked in by Pavel Demin <demin@…>, 10 years ago

delete duplicate license file and prepend GPLv3 header to all source code files

  • Property mode set to 100644
File size: 3.2 KB
Line 
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/** \class StatusPidFilter
20 *
21 * Removes all generated particles except electrons, muons, taus,
22 * and particles with status == 3.
23 *
24 * $Date$
25 * $Revision$
26 *
27 *
28 * \author J. Hirschauer - FNAL
29 *
30 */
31
32#include "modules/StatusPidFilter.h"
33
34#include "classes/DelphesClasses.h"
35#include "classes/DelphesFactory.h"
36#include "classes/DelphesFormula.h"
37
38#include "ExRootAnalysis/ExRootResult.h"
39#include "ExRootAnalysis/ExRootFilter.h"
40#include "ExRootAnalysis/ExRootClassifier.h"
41
42#include "TMath.h"
43#include "TString.h"
44#include "TFormula.h"
45#include "TRandom3.h"
46#include "TObjArray.h"
47#include "TDatabasePDG.h"
48#include "TLorentzVector.h"
49
50#include <algorithm>
51#include <stdexcept>
52#include <iostream>
53#include <sstream>
54
55using namespace std;
56
57//------------------------------------------------------------------------------
58
59StatusPidFilter::StatusPidFilter() :
60 fItInputArray(0)
61{
62}
63
64//------------------------------------------------------------------------------
65
66StatusPidFilter::~StatusPidFilter()
67{
68}
69
70//------------------------------------------------------------------------------
71
72void StatusPidFilter::Init()
73{
74 // PT threshold
75
76 fPTMin = GetDouble("PTMin", 0.5);
77
78 // import input array
79 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
80 fItInputArray = fInputArray->MakeIterator();
81
82 // create output array
83
84 fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
85}
86
87//------------------------------------------------------------------------------
88
89void StatusPidFilter::Finish()
90{
91 if(fItInputArray) delete fItInputArray;
92}
93
94//------------------------------------------------------------------------------
95
96void StatusPidFilter::Process()
97{
98 Candidate *candidate;
99 Int_t status, pdgCode;
100 Bool_t pass;
101
102 fItInputArray->Reset();
103 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
104 {
105 status = candidate->Status;
106 pdgCode = TMath::Abs(candidate->PID);
107
108 pass = kFALSE;
109
110 // status == 3
111 if(status == 3) pass = kTRUE;
112
113 // electrons, muons, taus and neutrinos
114 if(pdgCode > 10 && pdgCode < 17) pass = kTRUE;
115
116 // heavy quarks
117 if(pdgCode == 5 || pdgCode == 6) pass = kTRUE;
118
119 // Gauge bosons and other fundamental bosons
120 if(pdgCode > 22 && pdgCode < 43) pass = kTRUE;
121
122 if(!pass || candidate->Momentum.Pt() <= fPTMin) continue;
123
124 fOutputArray->Add(candidate);
125 }
126}
127
Note: See TracBrowser for help on using the repository browser.