Fork me on GitHub

source: git/modules/TrackPileUpSubtractor.cc@ 8b04b31

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 8b04b31 was 7c0fcd5, checked in by Pavel Demin <demin@…>, 10 years ago

delete duplicate license file and prepend GPLv3 header to all source code files

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