Fork me on GitHub

source: git/modules/ConstituentFilter.cc@ 0bd8fc8

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

remove svn tags and fix formatting

  • Property mode set to 100644
File size: 4.4 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
20/** \class ConstituentFilter
21 *
22 * Drops all input objects that are not constituents of any jet.
23 *
24 * \author P. Demin - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/ConstituentFilter.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
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55ConstituentFilter::ConstituentFilter()
56{
57}
58
59//------------------------------------------------------------------------------
60
61ConstituentFilter::~ConstituentFilter()
62{
63}
64
65//------------------------------------------------------------------------------
66
67void ConstituentFilter::Init()
68{
69 ExRootConfParam param;
70 Long_t i, size;
71 const TObjArray *array;
72 TIterator *iterator;
73
74 fJetPTMin = GetDouble("JetPTMin", 0.0);
75
76 // import input array(s)
77
78 param = GetParam("JetInputArray");
79 size = param.GetSize();
80 for(i = 0; i < size; ++i)
81 {
82 array = ImportArray(param[i].GetString());
83 iterator = array->MakeIterator();
84
85 fInputList.push_back(iterator);
86 }
87
88 param = GetParam("ConstituentInputArray");
89 size = param.GetSize();
90 for(i = 0; i < size/2; ++i)
91 {
92 array = ImportArray(param[i*2].GetString());
93 iterator = array->MakeIterator();
94
95 fInputMap[iterator] = ExportArray(param[i*2 + 1].GetString());
96 }
97}
98
99//------------------------------------------------------------------------------
100
101void ConstituentFilter::Finish()
102{
103 map< TIterator *, TObjArray * >::iterator itInputMap;
104 vector< TIterator * >::iterator itInputList;
105 TIterator *iterator;
106
107 for(itInputList = fInputList.begin(); itInputList != fInputList.end(); ++itInputList)
108 {
109 iterator = *itInputList;
110 if(iterator) delete iterator;
111 }
112
113 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
114 {
115 iterator = itInputMap->first;
116 if(iterator) delete iterator;
117 }
118}
119
120//------------------------------------------------------------------------------
121
122void ConstituentFilter::Process()
123{
124 Candidate *jet, *constituent;
125 map< TIterator *, TObjArray * >::iterator itInputMap;
126 vector< TIterator * >::iterator itInputList;
127 TIterator *iterator;
128 TObjArray *array;
129
130 // loop over all jet input arrays
131 for(itInputList = fInputList.begin(); itInputList != fInputList.end(); ++itInputList)
132 {
133 iterator = *itInputList;
134
135 // loop over all jets
136 iterator->Reset();
137 while((jet = static_cast<Candidate*>(iterator->Next())))
138 {
139 TIter itConstituents(jet->GetCandidates());
140
141 if(jet->Momentum.Pt() <= fJetPTMin) continue;
142
143 // loop over all constituents
144 while((constituent = static_cast<Candidate*>(itConstituents.Next())))
145 {
146 // set the IsConstituent flag
147 constituent->IsConstituent = 1;
148 }
149 }
150 }
151
152 // loop over all constituent input arrays
153 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
154 {
155 iterator = itInputMap->first;
156 array = itInputMap->second;
157
158 // loop over all constituents
159 iterator->Reset();
160 while((constituent = static_cast<Candidate*>(iterator->Next())))
161 {
162 // check the IsConstituent flag
163 if(constituent->IsConstituent)
164 {
165 array->Add(constituent);
166 }
167 }
168 }
169}
170
171//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.