Fork me on GitHub

source: git/modules/UniqueObjectFinder.cc@ d612dec

Last change on this file since d612dec was 77e9ae1, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

set Standard to Cpp03 in .clang-format

  • 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/** \class UniqueObjectFinder
20 *
21 * Finds uniquely identified photons, electrons and jets.
22 *
23 * \author P. Demin - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/UniqueObjectFinder.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
54UniqueObjectFinder::UniqueObjectFinder()
55{
56}
57
58//------------------------------------------------------------------------------
59
60UniqueObjectFinder::~UniqueObjectFinder()
61{
62}
63
64//------------------------------------------------------------------------------
65
66void UniqueObjectFinder::Init()
67{
68 // use GetUniqueID algorithm to find unique objects (faster than the default Overlaps method)
69 fUseUniqueID = GetBool("UseUniqueID", false);
70
71 // import arrays with output from other modules
72
73 ExRootConfParam param = GetParam("InputArray");
74 Long_t i, size;
75 const TObjArray *array;
76 TIterator *iterator;
77
78 fInputMap.clear();
79
80 size = param.GetSize();
81 for(i = 0; i < size / 2; ++i)
82 {
83 array = ImportArray(param[i * 2].GetString());
84 iterator = array->MakeIterator();
85
86 fInputMap.push_back(make_pair(iterator, ExportArray(param[i * 2 + 1].GetString())));
87 }
88}
89
90//------------------------------------------------------------------------------
91
92void UniqueObjectFinder::Finish()
93{
94 vector<pair<TIterator *, TObjArray *> >::iterator itInputMap;
95 TIterator *iterator;
96
97 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
98 {
99 iterator = itInputMap->first;
100
101 if(iterator) delete iterator;
102 }
103}
104
105//------------------------------------------------------------------------------
106
107void UniqueObjectFinder::Process()
108{
109 Candidate *candidate;
110 vector<pair<TIterator *, TObjArray *> >::iterator itInputMap;
111 TIterator *iterator;
112 TObjArray *array;
113
114 // loop over all input arrays
115 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
116 {
117 iterator = itInputMap->first;
118 array = itInputMap->second;
119
120 // loop over all candidates
121 iterator->Reset();
122 while((candidate = static_cast<Candidate *>(iterator->Next())))
123 {
124 if(Unique(candidate, itInputMap))
125 {
126 array->Add(candidate);
127 }
128 }
129 }
130}
131
132//------------------------------------------------------------------------------
133
134Bool_t UniqueObjectFinder::Unique(Candidate *candidate, vector<pair<TIterator *, TObjArray *> >::iterator itInputMap)
135{
136 Candidate *previousCandidate;
137 vector<pair<TIterator *, TObjArray *> >::iterator previousItInputMap;
138 TObjArray *array;
139
140 // loop over previous arrays
141 for(previousItInputMap = fInputMap.begin(); previousItInputMap != itInputMap; ++previousItInputMap)
142 {
143 array = previousItInputMap->second;
144 TIter iterator(array);
145
146 // loop over all candidates
147 iterator.Reset();
148 while((previousCandidate = static_cast<Candidate *>(iterator.Next())))
149 {
150 if(fUseUniqueID)
151 {
152 if(candidate->GetUniqueID() == previousCandidate->GetUniqueID())
153 {
154 return kFALSE;
155 }
156 }
157 else
158 {
159 if(candidate->Overlaps(previousCandidate))
160 {
161 return kFALSE;
162 }
163 }
164 }
165 }
166
167 return kTRUE;
168}
169
170//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.