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