Fork me on GitHub

source: git/modules/FastJetGridMedianEstimator.cc@ be1a3be

Last change on this file since be1a3be was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

  • 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/** \class FastJetGridMedianEstimator
20 *
21 * Computes median energy density per event using a fixed grid.
22 *
23 * \author M. Selvaggi - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/FastJetGridMedianEstimator.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "ExRootAnalysis/ExRootClassifier.h"
34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootResult.h"
36
37#include "TDatabasePDG.h"
38#include "TFormula.h"
39#include "TLorentzVector.h"
40#include "TMath.h"
41#include "TObjArray.h"
42#include "TRandom3.h"
43#include "TString.h"
44
45#include <algorithm>
46#include <iostream>
47#include <sstream>
48#include <stdexcept>
49#include <utility>
50#include <vector>
51
52#include "fastjet/ClusterSequence.hh"
53#include "fastjet/ClusterSequenceArea.hh"
54#include "fastjet/JetDefinition.hh"
55#include "fastjet/PseudoJet.hh"
56#include "fastjet/RectangularGrid.hh"
57#include "fastjet/Selector.hh"
58#include "fastjet/tools/JetMedianBackgroundEstimator.hh"
59
60#include "fastjet/tools/GridMedianBackgroundEstimator.hh"
61
62#include "fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh"
63#include "fastjet/plugins/CDFCones/fastjet/CDFMidPointPlugin.hh"
64#include "fastjet/plugins/SISCone/fastjet/SISConePlugin.hh"
65
66#include "fastjet/contribs/Nsubjettiness/ExtraRecombiners.hh"
67#include "fastjet/contribs/Nsubjettiness/Njettiness.hh"
68#include "fastjet/contribs/Nsubjettiness/NjettinessPlugin.hh"
69#include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh"
70
71using namespace std;
72using namespace fastjet;
73using namespace fastjet::contrib;
74
75//------------------------------------------------------------------------------
76
77FastJetGridMedianEstimator::FastJetGridMedianEstimator() :
78 fItInputArray(0)
79{
80}
81
82//------------------------------------------------------------------------------
83
84FastJetGridMedianEstimator::~FastJetGridMedianEstimator()
85{
86}
87
88//------------------------------------------------------------------------------
89
90void FastJetGridMedianEstimator::Init()
91{
92 ExRootConfParam param;
93 Long_t i, size;
94 Double_t drap, dphi, rapMin, rapMax;
95
96 // read rapidity ranges
97
98 param = GetParam("GridRange");
99 size = param.GetSize();
100
101 fEstimators.clear();
102 for(i = 0; i < size / 4; ++i)
103 {
104 rapMin = param[i * 4].GetDouble();
105 rapMax = param[i * 4 + 1].GetDouble();
106 drap = param[i * 4 + 2].GetDouble();
107 dphi = param[i * 4 + 3].GetDouble();
108 fEstimators.push_back(new GridMedianBackgroundEstimator(rapMin, rapMax, drap, dphi));
109 }
110
111 // import input array
112
113 fInputArray = ImportArray(GetString("InputArray", "Calorimeter/towers"));
114 fItInputArray = fInputArray->MakeIterator();
115
116 fRhoOutputArray = ExportArray(GetString("RhoOutputArray", "rho"));
117}
118
119//------------------------------------------------------------------------------
120
121void FastJetGridMedianEstimator::Finish()
122{
123 vector<GridMedianBackgroundEstimator *>::iterator itEstimators;
124
125 for(itEstimators = fEstimators.begin(); itEstimators != fEstimators.end(); ++itEstimators)
126 {
127 if(*itEstimators) delete *itEstimators;
128 }
129
130 if(fItInputArray) delete fItInputArray;
131}
132
133//------------------------------------------------------------------------------
134
135void FastJetGridMedianEstimator::Process()
136{
137 Candidate *candidate;
138 TLorentzVector momentum;
139 Int_t number;
140 Double_t rho = 0;
141 PseudoJet jet;
142 vector<PseudoJet> inputList, outputList;
143
144 vector<GridMedianBackgroundEstimator *>::iterator itEstimators;
145 ;
146
147 DelphesFactory *factory = GetFactory();
148
149 inputList.clear();
150
151 // loop over input objects
152 fItInputArray->Reset();
153 number = 0;
154 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
155 {
156 momentum = candidate->Momentum;
157 jet = PseudoJet(momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
158 jet.set_user_index(number);
159 inputList.push_back(jet);
160 ++number;
161 }
162
163 // compute rho and store it
164
165 for(itEstimators = fEstimators.begin(); itEstimators != fEstimators.end(); ++itEstimators)
166 {
167 (*itEstimators)->set_particles(inputList);
168
169 rho = (*itEstimators)->rho();
170
171 candidate = factory->NewCandidate();
172 candidate->Momentum.SetPtEtaPhiE(rho, 0.0, 0.0, rho);
173 candidate->Edges[0] = (*itEstimators)->rapmin();
174 candidate->Edges[1] = (*itEstimators)->rapmax();
175 fRhoOutputArray->Add(candidate);
176 }
177}
Note: See TracBrowser for help on using the repository browser.