Fork me on GitHub

source: git/modules/ConstituentFilter.cc@ a0b6d15

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