Fork me on GitHub

source: git/modules/FastJetGridMedianEstimator.cc@ df35033

ImprovedOutputFile Timing dual_readout llp
Last change on this file since df35033 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: 5.7 KB
RevLine 
[01f457a]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[01f457a]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.
[1fa50c2]9 *
[01f457a]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.
[1fa50c2]14 *
[01f457a]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
[6fb1a5d]19
20/** \class FastJetGridMedianEstimator
21 *
22 * Computes median energy density per event using a fixed grid.
23 *
24 * \author M. Selvaggi - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/FastJetGridMedianEstimator.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#include <vector>
51#include <utility>
52
53#include "fastjet/PseudoJet.hh"
54#include "fastjet/JetDefinition.hh"
55#include "fastjet/ClusterSequence.hh"
56#include "fastjet/Selector.hh"
57#include "fastjet/RectangularGrid.hh"
58#include "fastjet/ClusterSequenceArea.hh"
59#include "fastjet/tools/JetMedianBackgroundEstimator.hh"
60
61#include "fastjet/tools/GridMedianBackgroundEstimator.hh"
62
63#include "fastjet/plugins/SISCone/fastjet/SISConePlugin.hh"
64#include "fastjet/plugins/CDFCones/fastjet/CDFMidPointPlugin.hh"
65#include "fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh"
66
67#include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh"
68#include "fastjet/contribs/Nsubjettiness/Njettiness.hh"
69#include "fastjet/contribs/Nsubjettiness/NjettinessPlugin.hh"
70#include "fastjet/contribs/Nsubjettiness/WinnerTakeAllRecombiner.hh"
71
72using namespace std;
73using namespace fastjet;
74using namespace fastjet::contrib;
75
76
77//------------------------------------------------------------------------------
78
79FastJetGridMedianEstimator::FastJetGridMedianEstimator() :
80 fItInputArray(0)
81{
82
83}
84
85//------------------------------------------------------------------------------
86
87FastJetGridMedianEstimator::~FastJetGridMedianEstimator()
88{
89
90}
91
92//------------------------------------------------------------------------------
93
94void FastJetGridMedianEstimator::Init()
95{
96 // read eta ranges
97
98 ExRootConfParam param = GetParam("GridRange");
99 Long_t i, size;
100
101 fGrid.clear();
102 size = param.GetSize();
103 for(i = 0; i < size/4; ++i)
104 {
105 fGrid[make_pair(param[i*4].GetDouble(), param[i*4 + 1].GetDouble())] = make_pair(param[i*4 + 2].GetDouble(), param[i*4 + 3].GetDouble());
106 //cout<<param[i*4].GetDouble()<<","<< param[i*4 + 1].GetDouble()<<","<< param[i*4 + 2].GetDouble()<<","<< param[i*4 + 3].GetDouble()<<endl;
107
108 }
109
110
111 //cout<<fGrid[make_pair(0.0,2.5)].first<<","<<fGrid[make_pair(0.0,2.5)].second<<endl;
112
113 // import input array
114
115 fInputArray = ImportArray(GetString("InputArray", "Calorimeter/towers"));
116 fItInputArray = fInputArray->MakeIterator();
117
118 fRhoOutputArray = ExportArray(GetString("RhoOutputArray", "rho"));
119
120}
121
122//------------------------------------------------------------------------------
123
124void FastJetGridMedianEstimator::Finish()
125{
126 if(fItInputArray) delete fItInputArray;
127}
128
129//------------------------------------------------------------------------------
130
131void FastJetGridMedianEstimator::Process()
132{
133 Candidate *candidate;
134 TLorentzVector momentum;
135 Double_t deta, dphi, detaMin, detaMax;
136 Int_t number;
137 Double_t rho = 0;
138 PseudoJet jet;
139 vector<PseudoJet> inputList, outputList;
140
141 std::map< std::pair< Double_t , Double_t > , std::pair< Double_t , Double_t > >::iterator itGrid;
142
143 DelphesFactory *factory = GetFactory();
144
145 inputList.clear();
146
147 // loop over input objects
148 fItInputArray->Reset();
149 number = 0;
150 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
151 {
152 momentum = candidate->Momentum;
153 jet = PseudoJet(momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
154 jet.set_user_index(number);
155 inputList.push_back(jet);
156 ++number;
157 }
158
159
160 // compute rho and store it
161
162 // cout<<"caio"<<endl;
163 for(itGrid = fGrid.begin(); itGrid != fGrid.end(); ++itGrid)
164 {
165 //Selector select_rapidity = SelectorAbsRapRange(itEtaRangeMap->first, itEtaRangeMap->second);
166 // JetMedianBackgroundEstimator estimator(select_rapidity, *fDefinition, *fAreaDefinition);
167
168 //cout<<itGrid->first.first<<endl;
169
170 detaMin = (itGrid->first).first;
171 detaMax = (itGrid->first).second;
172 deta = (itGrid->second).first;
173 dphi = (itGrid->second).second;
174
175 //cout<<detaMin<<","<<detaMax<<","<<deta<<","<<dphi<<endl;
176
177
178 RectangularGrid grid(detaMin, detaMax, deta, dphi);
179 //cout<<grid.is_initialised()<<endl;
180 //cout<<grid.rapmin()<<","<<grid.rapmax()<<","<<grid.drap()<<","<<grid.dphi()<<endl;
181
182 GridMedianBackgroundEstimator estimator(grid);
183
184 estimator.set_particles(inputList);
185
186 //cout<<estimator.description()<<endl;
187
188 rho = estimator.rho();
189 //cout<<rho<<endl;
190
191
192 candidate = factory->NewCandidate();
193 candidate->Momentum.SetPtEtaPhiE(rho, 0.0, 0.0, rho);
194 candidate->Edges[0] = detaMin;
195 candidate->Edges[1] = detaMax;
196 fRhoOutputArray->Add(candidate);
197 }
198
199}
Note: See TracBrowser for help on using the repository browser.