Fork me on GitHub

source: git/modules/UniqueObjectFinder.cc@ 769f65b

ImprovedOutputFile Timing llp
Last change on this file since 769f65b was 70bb4cb, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

add possibility to calculate MET with smeared jets

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