Fork me on GitHub

source: git/modules/TruthVertexFinder.cc@ 3b02069

Last change on this file since 3b02069 was 5862d1f, checked in by michele <michele.selvaggi@…>, 4 years ago

added TruthVertexFinder module

  • Property mode set to 100644
File size: 4.4 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/** \class TruthVertexFinder
20 *
21 * Merges particles from pile-up sample into event
22 *
23 * \author M. Selvaggi - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/TruthVertexFinder.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesPileUpReader.h"
32#include "classes/DelphesTF2.h"
33
34#include "ExRootAnalysis/ExRootClassifier.h"
35#include "ExRootAnalysis/ExRootFilter.h"
36#include "ExRootAnalysis/ExRootResult.h"
37
38#include "TDatabasePDG.h"
39#include "TFormula.h"
40#include "TLorentzVector.h"
41#include "TMath.h"
42#include "TObjArray.h"
43#include "TRandom3.h"
44#include "TString.h"
45
46#include <algorithm>
47#include <iostream>
48#include <sstream>
49#include <stdexcept>
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55TruthVertexFinder::TruthVertexFinder() :
56fItInputArray(0), fItOutputArray(0)
57{
58}
59
60//------------------------------------------------------------------------------
61
62TruthVertexFinder::~TruthVertexFinder()
63{
64}
65
66//------------------------------------------------------------------------------
67
68void TruthVertexFinder::Init()
69{
70
71 fResolution = GetDouble("Resolution", 1E-06); // resolution in meters
72 // import input array
73 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
74 fItInputArray = fInputArray->MakeIterator();
75
76 // create output arrays
77 fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
78 //fItOutputArray = fVertexOutputArray->MakeIterator();
79}
80
81//------------------------------------------------------------------------------
82
83void TruthVertexFinder::Finish()
84{
85 if(fItInputArray) delete fItInputArray;
86 if(fItOutputArray) delete fItOutputArray;
87}
88
89//------------------------------------------------------------------------------
90
91void TruthVertexFinder::Process()
92{
93 Int_t nvtx = -1;
94 Float_t pt;
95 Candidate *candidate, *vertex;
96 DelphesFactory *factory;
97
98
99 fItInputArray->Reset();
100
101 factory = GetFactory();
102 vertex = factory->NewCandidate();
103
104 TLorentzVector vertexPosition(0., 0., 0., 0.);
105
106 nvtx=0;
107 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
108 {
109
110 const TLorentzVector &candidatePosition = candidate->Position;
111 const TLorentzVector &candidateMomentum = candidate->Momentum;
112
113 pt = candidateMomentum.Pt();
114 // check whether vertex already included, if so add particle
115 Bool_t old_vertex=false;
116 fItOutputArray = fVertexOutputArray->MakeIterator();
117 fItOutputArray->Reset();
118 while((vertex = static_cast<Candidate *>(fItOutputArray->Next())))
119 {
120 const TLorentzVector &vertexPosition = vertex->Position;
121 // check whether spatial difference is < 1 um, in that case assume it is the same vertex
122 if ( TMath::Abs((candidatePosition.P() - vertexPosition.P())) < fResolution*1.E3)
123 {
124 old_vertex=true;
125 vertex->AddCandidate(candidate);
126 if (TMath::Abs(candidate->Charge) > 0)
127 {
128 vertex->ClusterNDF += 1;
129 vertex->GenSumPT2 += pt*pt;
130 }
131 }
132 }
133
134 // else fill new vertex
135 if (!old_vertex)
136 {
137 vertex = factory->NewCandidate();
138 vertex->Position = candidatePosition;
139 vertex->ClusterIndex = nvtx;
140
141 if (TMath::Abs(candidate->Charge) > 0)
142 {
143 vertex->ClusterNDF = 1;
144 vertex->GenSumPT2 = pt*pt;
145 }
146 else
147 {
148 vertex->ClusterNDF = 0;
149 vertex->GenSumPT2 = 0.;
150 }
151 fVertexOutputArray->Add(vertex);
152 nvtx++;
153 }
154 }
155}
156
157//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.