Fork me on GitHub

source: git/modules/ConstituentFilter.cc@ 6182516

ImprovedOutputFile Timing
Last change on this file since 6182516 was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

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