Fork me on GitHub

source: git/modules/PdgCodeFilter.cc@ fc4e460

ImprovedOutputFile Timing dual_readout llp
Last change on this file since fc4e460 was be2222c, checked in by Michele <michele.selvaggi@…>, 10 years ago

"PdgCodeFilter module added. Updated card with NeutrinoFilter before GenJet clustering"

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