Fork me on GitHub

source: git/modules/Merger.cc@ a0c065d

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