Fork me on GitHub

source: git/modules/PdgCodeFilter.cc@ ef8a06d

ImprovedOutputFile Timing dual_readout
Last change on this file since ef8a06d was 56562d0, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

added require charge possibility in PdgCodeCodeFilter

  • Property mode set to 100644
File size: 3.5 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 fInvert = GetBool("Invert", false);
77
78 fRequireStatus = GetBool("RequireStatus", false);
79 fStatus = GetInt("Status", 1);
80
81 fRequireCharge = GetBool("RequireCharge", false);
82 fCharge = GetInt("Charge", 1);
83
84 // import input array
85 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
86 fItInputArray = fInputArray->MakeIterator();
87
88 param = GetParam("PdgCode");
89 size = param.GetSize();
90
91 // read PdgCodes to be filtered out from the data card
92
93 fPdgCodes.clear();
94 for(i = 0; i < size; ++i)
95 {
96 fPdgCodes.push_back(param[i].GetInt());
97 }
98
99 // create output array
100 fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
101}
102
103//------------------------------------------------------------------------------
104
105void PdgCodeFilter::Finish()
106{
107 if(fItInputArray) delete fItInputArray;
108}
109
110//------------------------------------------------------------------------------
111
112void PdgCodeFilter::Process()
113{
114 Candidate *candidate;
115 Int_t pdgCode;
116 Bool_t pass;
117 Double_t pt;
118
119 fItInputArray->Reset();
120 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
121 {
122 pdgCode = candidate->PID;
123 const TLorentzVector &candidateMomentum = candidate->Momentum;
124 pt = candidateMomentum.Pt();
125
126 if(pt < fPTMin) continue;
127 if(fRequireStatus && (candidate->Status != fStatus)) continue;
128 if(fRequireCharge && (candidate->Charge != fCharge)) continue;
129
130 pass = kTRUE;
131 if(find(fPdgCodes.begin(), fPdgCodes.end(), pdgCode) != fPdgCodes.end()) pass = kFALSE;
132
133 if(fInvert) pass = !pass;
134 if(pass) fOutputArray->Add(candidate);
135 }
136}
137
Note: See TracBrowser for help on using the repository browser.