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 LeptonDressing
|
---|
20 | *
|
---|
21 | *
|
---|
22 | *
|
---|
23 | * \author P. Demin && A. Mertens - UCL, Louvain-la-Neuve
|
---|
24 | *
|
---|
25 | */
|
---|
26 |
|
---|
27 | #include "modules/LeptonDressing.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 |
|
---|
50 | using namespace std;
|
---|
51 |
|
---|
52 | //------------------------------------------------------------------------------
|
---|
53 |
|
---|
54 | LeptonDressing::LeptonDressing() :
|
---|
55 | fItDressingInputArray(0), fItCandidateInputArray(0)
|
---|
56 | {
|
---|
57 | }
|
---|
58 |
|
---|
59 | //------------------------------------------------------------------------------
|
---|
60 |
|
---|
61 | LeptonDressing::~LeptonDressing()
|
---|
62 | {
|
---|
63 | }
|
---|
64 |
|
---|
65 | //------------------------------------------------------------------------------
|
---|
66 |
|
---|
67 | void LeptonDressing::Init()
|
---|
68 | {
|
---|
69 | fDeltaR = GetDouble("DeltaRMax", 0.4);
|
---|
70 |
|
---|
71 | // import input array(s)
|
---|
72 |
|
---|
73 | fDressingInputArray = ImportArray(GetString("DressingInputArray", "Calorimeter/photons"));
|
---|
74 | fItDressingInputArray = fDressingInputArray->MakeIterator();
|
---|
75 |
|
---|
76 | fCandidateInputArray = ImportArray(GetString("CandidateInputArray", "UniqueObjectFinder/electrons"));
|
---|
77 | fItCandidateInputArray = fCandidateInputArray->MakeIterator();
|
---|
78 |
|
---|
79 | // create output array
|
---|
80 |
|
---|
81 | fOutputArray = ExportArray(GetString("OutputArray", "electrons"));
|
---|
82 | }
|
---|
83 |
|
---|
84 | //------------------------------------------------------------------------------
|
---|
85 |
|
---|
86 | void LeptonDressing::Finish()
|
---|
87 | {
|
---|
88 | if(fItCandidateInputArray) delete fItCandidateInputArray;
|
---|
89 | if(fItDressingInputArray) delete fItDressingInputArray;
|
---|
90 | }
|
---|
91 |
|
---|
92 | //------------------------------------------------------------------------------
|
---|
93 |
|
---|
94 | void LeptonDressing::Process()
|
---|
95 | {
|
---|
96 | Candidate *candidate, *dressing, *mother;
|
---|
97 | TLorentzVector momentum;
|
---|
98 |
|
---|
99 | // loop over all input candidate
|
---|
100 | fItCandidateInputArray->Reset();
|
---|
101 | while((candidate = static_cast<Candidate *>(fItCandidateInputArray->Next())))
|
---|
102 | {
|
---|
103 | const TLorentzVector &candidateMomentum = candidate->Momentum;
|
---|
104 |
|
---|
105 | // loop over all input tracks
|
---|
106 | fItDressingInputArray->Reset();
|
---|
107 | momentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
|
---|
108 | while((dressing = static_cast<Candidate *>(fItDressingInputArray->Next())))
|
---|
109 | {
|
---|
110 | const TLorentzVector &dressingMomentum = dressing->Momentum;
|
---|
111 | if(dressingMomentum.Pt() > 0.1)
|
---|
112 | {
|
---|
113 | if(candidateMomentum.DeltaR(dressingMomentum) <= fDeltaR)
|
---|
114 | {
|
---|
115 | momentum += dressingMomentum;
|
---|
116 | }
|
---|
117 | }
|
---|
118 | }
|
---|
119 |
|
---|
120 | mother = candidate;
|
---|
121 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
122 |
|
---|
123 | candidate->Momentum += momentum;
|
---|
124 | candidate->AddCandidate(mother);
|
---|
125 |
|
---|
126 | fOutputArray->Add(candidate);
|
---|
127 | }
|
---|
128 | }
|
---|
129 |
|
---|
130 | //------------------------------------------------------------------------------
|
---|