Fork me on GitHub

source: git/modules/PdgCodeFilter.cc@ dacd9c5

ImprovedOutputFile Timing dual_readout llp
Last change on this file since dacd9c5 was 6cdc544, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix formatting

  • Property mode set to 100644
File size: 3.1 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 PdgCodeFilter
20 *
21 * Removes particles with specific PDG codes
22 *
23 * \author M. Selvaggi
24 *
25 */
26
27#include "modules/PdgCodeFilter.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "ExRootAnalysis/ExRootResult.h"
34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootClassifier.h"
36
37#include "TMath.h"
38#include "TString.h"
39#include "TFormula.h"
40#include "TRandom3.h"
41#include "TObjArray.h"
42#include "TDatabasePDG.h"
43#include "TLorentzVector.h"
44
45#include <algorithm>
46#include <stdexcept>
47#include <iostream>
48#include <sstream>
49
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54PdgCodeFilter::PdgCodeFilter() :
55 fItInputArray(0)
56{
57}
58
59//------------------------------------------------------------------------------
60
61PdgCodeFilter::~PdgCodeFilter()
62{
63}
64
65//------------------------------------------------------------------------------
66
67void PdgCodeFilter::Init()
68{
69
70 ExRootConfParam param;
71 Size_t i, size;
72
73 // PT threshold
74 fPTMin = GetDouble("PTMin", 0.0);
75
76 // import input array
77 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
78 fItInputArray = fInputArray->MakeIterator();
79
80 param = GetParam("PdgCode");
81 size = param.GetSize();
82
83 // read PdgCodes to be filtered out from the data card
84
85 fPdgCodes.clear();
86 for(i = 0; i < size; ++i)
87 {
88 fPdgCodes.push_back(param[i].GetInt());
89 }
90
91 // create output array
92 fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
93}
94
95//------------------------------------------------------------------------------
96
97void PdgCodeFilter::Finish()
98{
99 if(fItInputArray) delete fItInputArray;
100}
101
102//------------------------------------------------------------------------------
103
104void PdgCodeFilter::Process()
105{
106 Candidate *candidate;
107 Int_t pdgCode;
108 Bool_t pass;
109 Double_t pt;
110
111 fItInputArray->Reset();
112 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
113 {
114 pdgCode = candidate->PID;
115 const TLorentzVector &candidateMomentum = candidate->Momentum;
116 pt = candidateMomentum.Pt();
117
118 pass = kTRUE;
119
120 if(pt < fPTMin) pass = kFALSE;
121 if(find(fPdgCodes.begin(), fPdgCodes.end(), pdgCode) != fPdgCodes.end()) pass = kFALSE;
122
123 if(pass) fOutputArray->Add(candidate);
124 }
125}
126
Note: See TracBrowser for help on using the repository browser.