Fork me on GitHub

source: git/modules/TrackPileUpSubtractor.cc@ 255a666

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 255a666 was 715ab7c, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

fix formatting in TrackPileUpSubtractor

  • Property mode set to 100644
File size: 4.3 KB
Line 
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
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
55TrackPileUpSubtractor::TrackPileUpSubtractor()
56{
57}
58
59//------------------------------------------------------------------------------
60
61TrackPileUpSubtractor::~TrackPileUpSubtractor()
62{
63}
64
65//------------------------------------------------------------------------------
66
67void TrackPileUpSubtractor::Init()
68{
69 // import input array
70
71 fVertexInputArray = ImportArray(GetString("VertexInputArray", "PileUpMerger/vertices"));
72 fItVertexInputArray = fVertexInputArray->MakeIterator();
73
74 fZVertexResolution = GetDouble("ZVertexResolution", 0.005)*1.0E3;
75
76 fPTMin = GetDouble("PTMin", 0.);
77
78 // import arrays with output from other modules
79
80 ExRootConfParam param = GetParam("InputArray");
81 Long_t i, size;
82 const TObjArray *array;
83 TIterator *iterator;
84
85 size = param.GetSize();
86 for(i = 0; i < size/2; ++i)
87 {
88 array = ImportArray(param[i*2].GetString());
89 iterator = array->MakeIterator();
90
91 fInputMap[iterator] = ExportArray(param[i*2 + 1].GetString());
92 }
93}
94
95//------------------------------------------------------------------------------
96
97void TrackPileUpSubtractor::Finish()
98{
99 map< TIterator *, TObjArray * >::iterator itInputMap;
100 TIterator *iterator;
101
102 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
103 {
104 iterator = itInputMap->first;
105
106 if(iterator) delete iterator;
107 }
108
109 if(fItVertexInputArray) delete fItVertexInputArray;
110}
111
112//------------------------------------------------------------------------------
113
114void TrackPileUpSubtractor::Process()
115{
116 Candidate *candidate, *particle;
117 map< TIterator *, TObjArray * >::iterator itInputMap;
118 TIterator *iterator;
119 TObjArray *array;
120 Double_t z, zvtx=0;
121
122
123 // find z position of primary vertex
124
125 fItVertexInputArray->Reset();
126 while((candidate = static_cast<Candidate*>(fItVertexInputArray->Next())))
127 {
128 if(!candidate->IsPU)
129 {
130 zvtx = candidate->Position.Z();
131 // break;
132 }
133 }
134
135 // loop over all input arrays
136 for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
137 {
138 iterator = itInputMap->first;
139 array = itInputMap->second;
140
141 // loop over all candidates
142 iterator->Reset();
143 while((candidate = static_cast<Candidate*>(iterator->Next())))
144 {
145 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
146 z = particle->Position.Z();
147
148 // apply pile-up subtraction
149 // assume perfect pile-up subtraction for tracks outside fZVertexResolution
150
151 if(candidate->IsPU && TMath::Abs(z-zvtx) > fZVertexResolution)
152 {
153 candidate->IsRecoPU = 1;
154 }
155 else
156 {
157 candidate->IsRecoPU = 0;
158 if(candidate->Momentum.Pt() > fPTMin) array->Add(candidate);
159 }
160 }
161 }
162}
163
164//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.