Fork me on GitHub

source: git/modules/FastJetGridMedianEstimator.cc@ 4298fdb

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 4298fdb was 973b92a, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

update FastJet library to 3.1.3 and Nsubjettiness library to 2.2.1

  • Property mode set to 100644
File size: 5.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 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/ExtraRecombiners.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 ExRootConfParam param;
97 Long_t i, size;
98 Double_t drap, dphi, rapMin, rapMax;
99
100 // read rapidity ranges
101
102 param = GetParam("GridRange");
103 size = param.GetSize();
104
105 fEstimators.clear();
106 for(i = 0; i < size/4; ++i)
107 {
108 rapMin = param[i*4].GetDouble();
109 rapMax = param[i*4 + 1].GetDouble();
110 drap = param[i*4 + 2].GetDouble();
111 dphi = param[i*4 + 3].GetDouble();
112 fEstimators.push_back(new GridMedianBackgroundEstimator(rapMin, rapMax, drap, dphi));
113 }
114
115 // import input array
116
117 fInputArray = ImportArray(GetString("InputArray", "Calorimeter/towers"));
118 fItInputArray = fInputArray->MakeIterator();
119
120 fRhoOutputArray = ExportArray(GetString("RhoOutputArray", "rho"));
121}
122
123//------------------------------------------------------------------------------
124
125void FastJetGridMedianEstimator::Finish()
126{
127 vector< GridMedianBackgroundEstimator * >::iterator itEstimators;
128
129 for(itEstimators = fEstimators.begin(); itEstimators != fEstimators.end(); ++itEstimators)
130 {
131 if(*itEstimators) delete *itEstimators;
132 }
133
134 if(fItInputArray) delete fItInputArray;
135}
136
137//------------------------------------------------------------------------------
138
139void FastJetGridMedianEstimator::Process()
140{
141 Candidate *candidate;
142 TLorentzVector momentum;
143 Int_t number;
144 Double_t rho = 0;
145 PseudoJet jet;
146 vector< PseudoJet > inputList, outputList;
147
148 vector< GridMedianBackgroundEstimator * >::iterator itEstimators;;
149
150 DelphesFactory *factory = GetFactory();
151
152 inputList.clear();
153
154 // loop over input objects
155 fItInputArray->Reset();
156 number = 0;
157 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
158 {
159 momentum = candidate->Momentum;
160 jet = PseudoJet(momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
161 jet.set_user_index(number);
162 inputList.push_back(jet);
163 ++number;
164 }
165
166 // compute rho and store it
167
168 for(itEstimators = fEstimators.begin(); itEstimators != fEstimators.end(); ++itEstimators)
169 {
170 (*itEstimators)->set_particles(inputList);
171
172 rho = (*itEstimators)->rho();
173
174 candidate = factory->NewCandidate();
175 candidate->Momentum.SetPtEtaPhiE(rho, 0.0, 0.0, rho);
176 candidate->Edges[0] = (*itEstimators)->rapmin();
177 candidate->Edges[1] = (*itEstimators)->rapmax();
178 fRhoOutputArray->Add(candidate);
179 }
180}
Note: See TracBrowser for help on using the repository browser.