Fork me on GitHub

source: git/modules/LeptonDressing.cc@ b438187

Last change on this file since b438187 was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

  • Property mode set to 100644
File size: 3.7 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[d7d2da3]19/** \class LeptonDressing
[1fa50c2]20 *
[d7d2da3]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"
[341014c]34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootResult.h"
[d7d2da3]36
37#include "TDatabasePDG.h"
[341014c]38#include "TFormula.h"
[d7d2da3]39#include "TLorentzVector.h"
[341014c]40#include "TMath.h"
41#include "TObjArray.h"
42#include "TRandom3.h"
43#include "TString.h"
[d7d2da3]44
[341014c]45#include <algorithm>
[d7d2da3]46#include <iostream>
47#include <sstream>
[341014c]48#include <stdexcept>
[d7d2da3]49
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54LeptonDressing::LeptonDressing() :
[341014c]55 fItDressingInputArray(0), fItCandidateInputArray(0)
[d7d2da3]56{
57}
58
59//------------------------------------------------------------------------------
60
61LeptonDressing::~LeptonDressing()
62{
63}
64
65//------------------------------------------------------------------------------
66
67void 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();
[341014c]75
[d7d2da3]76 fCandidateInputArray = ImportArray(GetString("CandidateInputArray", "UniqueObjectFinder/electrons"));
77 fItCandidateInputArray = fCandidateInputArray->MakeIterator();
[341014c]78
[d7d2da3]79 // create output array
80
81 fOutputArray = ExportArray(GetString("OutputArray", "electrons"));
82}
83
84//------------------------------------------------------------------------------
85
86void LeptonDressing::Finish()
87{
88 if(fItCandidateInputArray) delete fItCandidateInputArray;
89 if(fItDressingInputArray) delete fItDressingInputArray;
90}
91
92//------------------------------------------------------------------------------
93
94void LeptonDressing::Process()
95{
96 Candidate *candidate, *dressing, *mother;
97 TLorentzVector momentum;
[341014c]98
[d7d2da3]99 // loop over all input candidate
100 fItCandidateInputArray->Reset();
[341014c]101 while((candidate = static_cast<Candidate *>(fItCandidateInputArray->Next())))
[d7d2da3]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);
[341014c]108 while((dressing = static_cast<Candidate *>(fItDressingInputArray->Next())))
[d7d2da3]109 {
110 const TLorentzVector &dressingMomentum = dressing->Momentum;
[341014c]111 if(dressingMomentum.Pt() > 0.1)
[d7d2da3]112 {
113 if(candidateMomentum.DeltaR(dressingMomentum) <= fDeltaR)
114 {
115 momentum += dressingMomentum;
116 }
117 }
118 }
119
120 mother = candidate;
[341014c]121 candidate = static_cast<Candidate *>(candidate->Clone());
[d7d2da3]122
123 candidate->Momentum += momentum;
124 candidate->AddCandidate(mother);
[341014c]125
[d7d2da3]126 fOutputArray->Add(candidate);
127 }
128}
129
130//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.