[700] | 1 |
|
---|
[814] | 2 | /** \class Merger
|
---|
| 3 | *
|
---|
| 4 | * Merges multiple input arrays into one output array
|
---|
| 5 | * and sums transverse momenta of all input objects.
|
---|
| 6 | *
|
---|
| 7 | * $Date: 2013-04-26 10:39:14 +0000 (Fri, 26 Apr 2013) $
|
---|
| 8 | * $Revision: 1099 $
|
---|
| 9 | *
|
---|
| 10 | *
|
---|
| 11 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
| 12 | *
|
---|
| 13 | */
|
---|
| 14 |
|
---|
[760] | 15 | #include "modules/Merger.h"
|
---|
[700] | 16 |
|
---|
| 17 | #include "classes/DelphesClasses.h"
|
---|
| 18 | #include "classes/DelphesFactory.h"
|
---|
[766] | 19 | #include "classes/DelphesFormula.h"
|
---|
[700] | 20 |
|
---|
[703] | 21 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
| 22 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 23 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
| 24 |
|
---|
| 25 | #include "TMath.h"
|
---|
| 26 | #include "TString.h"
|
---|
[700] | 27 | #include "TFormula.h"
|
---|
[703] | 28 | #include "TRandom3.h"
|
---|
| 29 | #include "TObjArray.h"
|
---|
| 30 | #include "TDatabasePDG.h"
|
---|
[700] | 31 | #include "TLorentzVector.h"
|
---|
| 32 |
|
---|
[703] | 33 | #include <algorithm>
|
---|
| 34 | #include <stdexcept>
|
---|
[700] | 35 | #include <iostream>
|
---|
[703] | 36 | #include <sstream>
|
---|
[700] | 37 |
|
---|
| 38 | using namespace std;
|
---|
| 39 |
|
---|
| 40 | //------------------------------------------------------------------------------
|
---|
| 41 |
|
---|
[760] | 42 | Merger::Merger()
|
---|
[700] | 43 | {
|
---|
| 44 | }
|
---|
| 45 |
|
---|
| 46 | //------------------------------------------------------------------------------
|
---|
| 47 |
|
---|
[760] | 48 | Merger::~Merger()
|
---|
[700] | 49 | {
|
---|
| 50 | }
|
---|
| 51 |
|
---|
| 52 | //------------------------------------------------------------------------------
|
---|
| 53 |
|
---|
[760] | 54 | void Merger::Init()
|
---|
[700] | 55 | {
|
---|
| 56 | // import arrays with output from other modules
|
---|
| 57 |
|
---|
| 58 | ExRootConfParam param = GetParam("InputArray");
|
---|
| 59 | Long_t i, size;
|
---|
| 60 | const TObjArray *array;
|
---|
| 61 | TIterator *iterator;
|
---|
| 62 |
|
---|
| 63 | size = param.GetSize();
|
---|
[894] | 64 | for(i = 0; i < size; ++i)
|
---|
[700] | 65 | {
|
---|
[894] | 66 | array = ImportArray(param[i].GetString());
|
---|
[700] | 67 | iterator = array->MakeIterator();
|
---|
| 68 |
|
---|
[894] | 69 | fInputList.push_back(iterator);
|
---|
[700] | 70 | }
|
---|
| 71 |
|
---|
[760] | 72 | // create output arrays
|
---|
[700] | 73 |
|
---|
[760] | 74 | fOutputArray = ExportArray(GetString("OutputArray", "candidates"));
|
---|
[700] | 75 |
|
---|
[760] | 76 | fMomentumOutputArray = ExportArray(GetString("MomentumOutputArray", "momentum"));
|
---|
[894] | 77 |
|
---|
| 78 | fEnergyOutputArray = ExportArray(GetString("EnergyOutputArray", "energy"));
|
---|
[700] | 79 | }
|
---|
| 80 |
|
---|
| 81 | //------------------------------------------------------------------------------
|
---|
| 82 |
|
---|
[760] | 83 | void Merger::Finish()
|
---|
[700] | 84 | {
|
---|
[894] | 85 | vector< TIterator * >::iterator itInputList;
|
---|
[700] | 86 | TIterator *iterator;
|
---|
| 87 |
|
---|
[894] | 88 | for(itInputList = fInputList.begin(); itInputList != fInputList.end(); ++itInputList)
|
---|
[700] | 89 | {
|
---|
[894] | 90 | iterator = *itInputList;
|
---|
[700] | 91 | if(iterator) delete iterator;
|
---|
| 92 | }
|
---|
| 93 | }
|
---|
| 94 |
|
---|
| 95 | //------------------------------------------------------------------------------
|
---|
| 96 |
|
---|
[760] | 97 | void Merger::Process()
|
---|
[700] | 98 | {
|
---|
| 99 | Candidate *candidate;
|
---|
[760] | 100 | TLorentzVector momentum;
|
---|
[894] | 101 | Double_t sumPT, sumE;
|
---|
| 102 | vector< TIterator * >::iterator itInputList;
|
---|
[700] | 103 | TIterator *iterator;
|
---|
[894] | 104 |
|
---|
[760] | 105 | DelphesFactory *factory = GetFactory();
|
---|
| 106 |
|
---|
| 107 | momentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
|
---|
[894] | 108 | sumPT = 0;
|
---|
| 109 | sumE = 0;
|
---|
[700] | 110 |
|
---|
| 111 | // loop over all input arrays
|
---|
[894] | 112 | for(itInputList = fInputList.begin(); itInputList != fInputList.end(); ++itInputList)
|
---|
[700] | 113 | {
|
---|
[894] | 114 | iterator = *itInputList;
|
---|
[700] | 115 |
|
---|
| 116 | // loop over all candidates
|
---|
| 117 | iterator->Reset();
|
---|
| 118 | while((candidate = static_cast<Candidate*>(iterator->Next())))
|
---|
| 119 | {
|
---|
| 120 | const TLorentzVector &candidateMomentum = candidate->Momentum;
|
---|
| 121 |
|
---|
[894] | 122 | momentum += candidateMomentum;
|
---|
| 123 | sumPT += candidateMomentum.Pt();
|
---|
| 124 | sumE += candidateMomentum.E();
|
---|
[700] | 125 |
|
---|
[894] | 126 | fOutputArray->Add(candidate);
|
---|
[700] | 127 | }
|
---|
| 128 | }
|
---|
[760] | 129 |
|
---|
| 130 | candidate = factory->NewCandidate();
|
---|
| 131 |
|
---|
| 132 | candidate->Position.SetXYZT(0.0, 0.0, 0.0, 0.0);
|
---|
| 133 | candidate->Momentum = momentum;
|
---|
| 134 |
|
---|
| 135 | fMomentumOutputArray->Add(candidate);
|
---|
[700] | 136 |
|
---|
[894] | 137 | candidate = factory->NewCandidate();
|
---|
| 138 |
|
---|
| 139 | candidate->Position.SetXYZT(0.0, 0.0, 0.0, 0.0);
|
---|
| 140 | candidate->Momentum.SetPtEtaPhiE(sumPT, 0.0, 0.0, sumE);
|
---|
| 141 |
|
---|
[1099] | 142 | fEnergyOutputArray->Add(candidate);
|
---|
| 143 | }
|
---|
[894] | 144 |
|
---|
[700] | 145 | //------------------------------------------------------------------------------
|
---|