Fork me on GitHub

source: git/modules/UniqueObjectFinder.cc@ 3e2bb2b

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 3e2bb2b was 93da593, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

replace map with vector in UniqueObjectFinder

  • Property mode set to 100644
File size: 4.1 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 // import arrays with output from other modules
70
71 ExRootConfParam param = GetParam("InputArray");
72 Long_t i, size;
73 const TObjArray *array;
74 TIterator *iterator;
75
76 fInputMap.clear();
77
78 size = param.GetSize();
79 for(i = 0; i < size/2; ++i)
80 {
81 array = ImportArray(param[i*2].GetString());
82 iterator = array->MakeIterator();
83
84 fInputMap.push_back(make_pair(iterator, ExportArray(param[i*2 + 1].GetString())));
85 }
86}
87
88//------------------------------------------------------------------------------
89
90void UniqueObjectFinder::Finish()
91{
92 vector< pair< TIterator *, TObjArray * > >::iterator itInputMap;
93 TIterator *iterator;
94
95 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
96 {
97 iterator = itInputMap->first;
98
99 if(iterator) delete iterator;
100 }
101}
102
103//------------------------------------------------------------------------------
104
105void UniqueObjectFinder::Process()
106{
107 Candidate *candidate;
108 vector< pair< TIterator *, TObjArray * > >::iterator itInputMap;
109 TIterator *iterator;
110 TObjArray *array;
111
112 // loop over all input arrays
113 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
114 {
115 iterator = itInputMap->first;
116 array = itInputMap->second;
117
118 // loop over all candidates
119 iterator->Reset();
120 while((candidate = static_cast<Candidate*>(iterator->Next())))
121 {
122 if(Unique(candidate, itInputMap))
123 {
124 array->Add(candidate);
125 }
126 }
127 }
128}
129
130//------------------------------------------------------------------------------
131
132Bool_t UniqueObjectFinder::Unique(Candidate *candidate, vector< pair< TIterator *, TObjArray * > >::iterator itInputMap)
133{
134 Candidate *previousCandidate;
135 vector< pair< TIterator *, TObjArray * > >::iterator previousItInputMap;
136 TObjArray *array;
137
138 // loop over previous arrays
139 for(previousItInputMap = fInputMap.begin(); previousItInputMap != itInputMap; ++previousItInputMap)
140 {
141 array = previousItInputMap->second;
142 TIter iterator(array);
143
144 // loop over all candidates
145 iterator.Reset();
146 while((previousCandidate = static_cast<Candidate*>(iterator.Next())))
147 {
148 if(candidate->Overlaps(previousCandidate))
149 {
150 return kFALSE;
151 }
152 }
153 }
154
155 return kTRUE;
156}
157
158//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.