Fork me on GitHub

source: git/modules/Merger.cc@ a643a7f

Last change on this file since a643a7f 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: 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/** \class Merger
20 *
21 * Merges multiple input arrays into one output array
22 * and sums transverse momenta of all input objects.
23 *
24 * \author P. Demin - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/Merger.h"
29
30#include "classes/DelphesClasses.h"
31#include "classes/DelphesFactory.h"
32#include "classes/DelphesFormula.h"
33
34#include "ExRootAnalysis/ExRootClassifier.h"
35#include "ExRootAnalysis/ExRootFilter.h"
36#include "ExRootAnalysis/ExRootResult.h"
37
38#include "TDatabasePDG.h"
39#include "TFormula.h"
40#include "TLorentzVector.h"
41#include "TMath.h"
42#include "TObjArray.h"
43#include "TRandom3.h"
44#include "TString.h"
45
46#include <algorithm>
47#include <iostream>
48#include <sstream>
49#include <stdexcept>
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55Merger::Merger()
56{
57}
58
59//------------------------------------------------------------------------------
60
61Merger::~Merger()
62{
63}
64
65//------------------------------------------------------------------------------
66
67void Merger::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; ++i)
78 {
79 array = ImportArray(param[i].GetString());
80 iterator = array->MakeIterator();
81
82 fInputList.push_back(iterator);
83 }
84
85 // create output arrays
86
87 fOutputArray = ExportArray(GetString("OutputArray", "candidates"));
88
89 fMomentumOutputArray = ExportArray(GetString("MomentumOutputArray", "momentum"));
90
91 fEnergyOutputArray = ExportArray(GetString("EnergyOutputArray", "energy"));
92}
93
94//------------------------------------------------------------------------------
95
96void Merger::Finish()
97{
98 vector<TIterator *>::iterator itInputList;
99 TIterator *iterator;
100
101 for(itInputList = fInputList.begin(); itInputList != fInputList.end(); ++itInputList)
102 {
103 iterator = *itInputList;
104 if(iterator) delete iterator;
105 }
106}
107
108//------------------------------------------------------------------------------
109
110void Merger::Process()
111{
112 Candidate *candidate;
113 TLorentzVector momentum;
114 Double_t sumPT, sumE;
115 vector<TIterator *>::iterator itInputList;
116 TIterator *iterator;
117
118 DelphesFactory *factory = GetFactory();
119
120 momentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
121 sumPT = 0;
122 sumE = 0;
123
124 // loop over all input arrays
125 for(itInputList = fInputList.begin(); itInputList != fInputList.end(); ++itInputList)
126 {
127 iterator = *itInputList;
128
129 // loop over all candidates
130 iterator->Reset();
131 while((candidate = static_cast<Candidate *>(iterator->Next())))
132 {
133 const TLorentzVector &candidateMomentum = candidate->Momentum;
134
135 momentum += candidateMomentum;
136 sumPT += candidateMomentum.Pt();
137 sumE += candidateMomentum.E();
138
139 fOutputArray->Add(candidate);
140 }
141 }
142
143 candidate = factory->NewCandidate();
144
145 candidate->Position.SetXYZT(0.0, 0.0, 0.0, 0.0);
146 candidate->Momentum = momentum;
147
148 fMomentumOutputArray->Add(candidate);
149
150 candidate = factory->NewCandidate();
151
152 candidate->Position.SetXYZT(0.0, 0.0, 0.0, 0.0);
153 candidate->Momentum.SetPtEtaPhiE(sumPT, 0.0, 0.0, sumE);
154
155 fEnergyOutputArray->Add(candidate);
156}
157
158//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.