Fork me on GitHub

source: git/modules/TrackPileUpSubtractor.cc@ 5a697dde

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5a697dde was 781edf86, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

fix charged Pu subtraction

  • Property mode set to 100644
File size: 4.6 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
20/** \class TrackPileUpSubtractor
21 *
22 * Subtract pile-up contribution from tracks.
23 *
24 * \author P. Demin - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/TrackPileUpSubtractor.h"
29
30#include "classes/DelphesClasses.h"
31#include "classes/DelphesFactory.h"
32#include "classes/DelphesFormula.h"
33
34#include "ExRootAnalysis/ExRootResult.h"
35#include "ExRootAnalysis/ExRootFilter.h"
36#include "ExRootAnalysis/ExRootClassifier.h"
37
38#include "TMath.h"
39#include "TString.h"
40#include "TFormula.h"
41#include "TRandom3.h"
42#include "TObjArray.h"
43#include "TDatabasePDG.h"
44#include "TLorentzVector.h"
45
46#include <algorithm>
47#include <stdexcept>
48#include <iostream>
49#include <sstream>
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54
[4de073b]55TrackPileUpSubtractor::TrackPileUpSubtractor() :
56fFormula(0)
[d7d2da3]57{
[4de073b]58 fFormula = new DelphesFormula;
[d7d2da3]59}
60
61//------------------------------------------------------------------------------
62
63TrackPileUpSubtractor::~TrackPileUpSubtractor()
64{
[4de073b]65 if(fFormula) delete fFormula;
[d7d2da3]66}
67
68//------------------------------------------------------------------------------
69
70void TrackPileUpSubtractor::Init()
71{
[715ab7c]72 // import input array
[22dc7fd]73
74 fVertexInputArray = ImportArray(GetString("VertexInputArray", "PileUpMerger/vertices"));
75 fItVertexInputArray = fVertexInputArray->MakeIterator();
[715ab7c]76
[4de073b]77 // read resolution formula in m
78 fFormula->Compile(GetString("ZVertexResolution", "0.001"));
[d7d2da3]79
[e2d3977]80 fPTMin = GetDouble("PTMin", 0.);
[715ab7c]81
[d7d2da3]82 // import arrays with output from other modules
[715ab7c]83
[d7d2da3]84 ExRootConfParam param = GetParam("InputArray");
85 Long_t i, size;
86 const TObjArray *array;
87 TIterator *iterator;
88
89 size = param.GetSize();
90 for(i = 0; i < size/2; ++i)
91 {
92 array = ImportArray(param[i*2].GetString());
93 iterator = array->MakeIterator();
94
[498cda4]95 fInputMap[iterator] = ExportArray(param[i*2 + 1].GetString());
[d7d2da3]96 }
97}
98
99//------------------------------------------------------------------------------
100
101void TrackPileUpSubtractor::Finish()
102{
103 map< TIterator *, TObjArray * >::iterator itInputMap;
104 TIterator *iterator;
105
106 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
107 {
108 iterator = itInputMap->first;
109
110 if(iterator) delete iterator;
111 }
[22dc7fd]112
113 if(fItVertexInputArray) delete fItVertexInputArray;
[d7d2da3]114}
115
116//------------------------------------------------------------------------------
117
118void TrackPileUpSubtractor::Process()
119{
120 Candidate *candidate, *particle;
121 map< TIterator *, TObjArray * >::iterator itInputMap;
122 TIterator *iterator;
123 TObjArray *array;
[24d005f]124 Double_t z, zvtx=0;
[4de073b]125 Double_t pt, eta, phi, e;
[22dc7fd]126
[715ab7c]127
[22dc7fd]128 // find z position of primary vertex
[715ab7c]129
[22dc7fd]130 fItVertexInputArray->Reset();
131 while((candidate = static_cast<Candidate*>(fItVertexInputArray->Next())))
132 {
[54b6dfc]133 if(!candidate->IsPU)
[22dc7fd]134 {
[715ab7c]135 zvtx = candidate->Position.Z();
136 // break;
[22dc7fd]137 }
138 }
139
[d7d2da3]140 // loop over all input arrays
141 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
142 {
143 iterator = itInputMap->first;
144 array = itInputMap->second;
145
146 // loop over all candidates
147 iterator->Reset();
148 while((candidate = static_cast<Candidate*>(iterator->Next())))
149 {
150 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
[4de073b]151 const TLorentzVector &candidateMomentum = particle->Momentum;
152
153 eta = candidateMomentum.Eta();
154 pt = candidateMomentum.Pt();
155 phi = candidateMomentum.Phi();
156 e = candidateMomentum.E();
157
[d7d2da3]158 z = particle->Position.Z();
[715ab7c]159
[d7d2da3]160 // apply pile-up subtraction
161 // assume perfect pile-up subtraction for tracks outside fZVertexResolution
[715ab7c]162
[781edf86]163 if(candidate->Charge !=0 && candidate->IsPU && TMath::Abs(z-zvtx) > fFormula->Eval(pt, eta, phi, e)* 1.0e3)
[715ab7c]164 {
165 candidate->IsRecoPU = 1;
166 }
167 else
[e2d3977]168 {
[715ab7c]169 candidate->IsRecoPU = 0;
170 if(candidate->Momentum.Pt() > fPTMin) array->Add(candidate);
[e2d3977]171 }
[d7d2da3]172 }
173 }
174}
Note: See TracBrowser for help on using the repository browser.