Fork me on GitHub

source: git/modules/UniqueObjectFinder.cc@ 4e09c3a

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 4e09c3a was cab38f6, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

remove svn tags and fix formatting

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