Fork me on GitHub

source: svn/trunk/modules/LeptonDressing.cc@ 1247

Last change on this file since 1247 was 915, checked in by Pavel Demin, 12 years ago

add LeptonDressing module

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id Revision Date
File size: 2.9 KB
Line 
1
2/** \class LeptonDressing
3 *
4 *
5 *
6 *
7 *
8 * \author P. Demin && A. Mertens - UCL, Louvain-la-Neuve
9 *
10 */
11
12#include "modules/LeptonDressing.h"
13
14#include "classes/DelphesClasses.h"
15#include "classes/DelphesFactory.h"
16#include "classes/DelphesFormula.h"
17
18#include "ExRootAnalysis/ExRootResult.h"
19#include "ExRootAnalysis/ExRootFilter.h"
20#include "ExRootAnalysis/ExRootClassifier.h"
21
22#include "TMath.h"
23#include "TString.h"
24#include "TFormula.h"
25#include "TRandom3.h"
26#include "TObjArray.h"
27#include "TDatabasePDG.h"
28#include "TLorentzVector.h"
29
30#include <algorithm>
31#include <stdexcept>
32#include <iostream>
33#include <sstream>
34
35using namespace std;
36
37//------------------------------------------------------------------------------
38
39LeptonDressing::LeptonDressing() :
40 fItDressingInputArray(0), fItCandidateInputArray(0)
41{
42}
43
44//------------------------------------------------------------------------------
45
46LeptonDressing::~LeptonDressing()
47{
48}
49
50//------------------------------------------------------------------------------
51
52void LeptonDressing::Init()
53{
54 fDeltaR = GetDouble("DeltaRMax", 0.4);
55
56 // import input array(s)
57
58 fDressingInputArray = ImportArray(GetString("DressingInputArray", "Calorimeter/photons"));
59 fItDressingInputArray = fDressingInputArray->MakeIterator();
60
61 fCandidateInputArray = ImportArray(GetString("CandidateInputArray", "UniqueObjectFinder/electrons"));
62 fItCandidateInputArray = fCandidateInputArray->MakeIterator();
63
64 // create output array
65
66 fOutputArray = ExportArray(GetString("OutputArray", "electrons"));
67}
68
69//------------------------------------------------------------------------------
70
71void LeptonDressing::Finish()
72{
73 if(fItCandidateInputArray) delete fItCandidateInputArray;
74 if(fItDressingInputArray) delete fItDressingInputArray;
75}
76
77//------------------------------------------------------------------------------
78
79void LeptonDressing::Process()
80{
81 Candidate *candidate, *dressing, *mother;
82 TLorentzVector momentum;
83
84 // loop over all input candidate
85 fItCandidateInputArray->Reset();
86 while((candidate = static_cast<Candidate*>(fItCandidateInputArray->Next())))
87 {
88 const TLorentzVector &candidateMomentum = candidate->Momentum;
89
90 // loop over all input tracks
91 fItDressingInputArray->Reset();
92 momentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
93 while((dressing = static_cast<Candidate*>(fItDressingInputArray->Next())))
94 {
95 const TLorentzVector &dressingMomentum = dressing->Momentum;
96 if (dressingMomentum.Pt() > 0.1)
97 {
98 if(candidateMomentum.DeltaR(dressingMomentum) <= fDeltaR)
99 {
100 momentum += dressingMomentum;
101 }
102 }
103 }
104
105 mother = candidate;
106 candidate = static_cast<Candidate*>(candidate->Clone());
107
108 candidate->Momentum += momentum;
109 candidate->AddCandidate(mother);
110
111 fOutputArray->Add(candidate);
112 }
113}
114
115//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.