Fork me on GitHub

source: git/modules/VertexSorter.cc@ 4569b59

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 4569b59 was 4569b59, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

added try catch for BS collection

  • Property mode set to 100644
File size: 7.1 KB
Line 
1/** \class VertexSorter
2 *
3 * Merges particles from pile-up sample into event
4 *
5 *
6 * $Date: 2013-02-12 15:13:59 +0100 (Tue, 12 Feb 2013) $
7 * $Revision: 907 $
8 *
9 *
10 * \author M. Selvaggi - UCL, Louvain-la-Neuve
11 *
12 */
13
14#include <unordered_map>
15#include "modules/VertexSorter.h"
16
17//#include "CLHEP/Units/GlobalSystemOfUnits.h"
18//#include "CLHEP/Units/GlobalPhysicalConstants.h"
19
20#include "classes/DelphesClasses.h"
21#include "classes/DelphesFactory.h"
22#include "classes/DelphesFormula.h"
23#include "classes/DelphesPileUpReader.h"
24
25#include "ExRootAnalysis/ExRootResult.h"
26#include "ExRootAnalysis/ExRootFilter.h"
27#include "ExRootAnalysis/ExRootClassifier.h"
28
29#include "TMath.h"
30#include "TString.h"
31#include "TFormula.h"
32#include "TRandom3.h"
33#include "TObjArray.h"
34#include "TDatabasePDG.h"
35#include "TLorentzVector.h"
36#include "TMatrixT.h"
37#include "TVector3.h"
38
39static const double mm = 1.;
40static const double m = 1000.*mm;
41static const double ns = 1.;
42static const double s = 1.e+9 *ns;
43static const double c_light = 2.99792458e+8 * m/s;
44
45//------------------------------------------------------------------------------
46
47VertexSorter::VertexSorter() :
48 fInputArray(NULL), fTrackInputArray(NULL), fItTrackInputArray(NULL), fJetInputArray(NULL), fItJetInputArray(NULL), fOutputArray(NULL)
49{
50}
51
52//------------------------------------------------------------------------------
53
54VertexSorter::~VertexSorter()
55{
56}
57
58//------------------------------------------------------------------------------
59
60void VertexSorter::Init()
61{
62 fInputArray = ImportArray(GetString("InputArray", "VertexFinder/clusters"));
63
64 fTrackInputArray = ImportArray(GetString("TrackInputArray", "VertexFinder/tracks"));
65 fItTrackInputArray = fTrackInputArray->MakeIterator();
66
67 if (string (GetString("JetInputArray", "")) != "")
68 {
69 fJetInputArray = ImportArray(GetString("JetInputArray", ""));
70 fItJetInputArray = fJetInputArray->MakeIterator();
71 }
72
73 // import beamspot
74 try
75 {
76 fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
77 }
78 catch(runtime_error &e)
79 {
80 fBeamSpotInputArray = 0;
81 }
82
83 fOutputArray = ExportArray(GetString("OutputArray", "clusters"));
84
85 fMethod = GetString ("Method", "BTV");
86}
87
88//------------------------------------------------------------------------------
89
90void VertexSorter::Finish()
91{
92 if(fItTrackInputArray) delete fItTrackInputArray;
93 if(fItJetInputArray) delete fItJetInputArray;
94}
95
96//------------------------------------------------------------------------------
97//
98bool VertexSorter::secondDescending (pair<unsigned, double> pair0, pair<unsigned, double> pair1)
99{
100 return (pair0.second > pair1.second);
101}
102bool VertexSorter::secondAscending (pair<unsigned, double> pair0, pair<unsigned, double> pair1)
103{
104 return (pair0.second < pair1.second);
105}
106
107void VertexSorter::Process()
108{
109 Candidate *candidate, *jetCandidate, *beamSpotCandidate;
110 unordered_map<int, unsigned> clusterIDToIndex;
111 unordered_map<int, double> clusterIDToSumPT2;
112 vector<pair<int, double> > sortedClusterIDs;
113
114 for (int iCluster = 0; iCluster < fInputArray->GetEntries (); iCluster++)
115 {
116 const Candidate &cluster = *((Candidate *) fInputArray->At (iCluster));
117 clusterIDToIndex[cluster.ClusterIndex] = iCluster;
118 clusterIDToSumPT2[cluster.ClusterIndex] = 0.0;
119 }
120
121 if (fMethod == "BTV")
122 {
123 if (!fJetInputArray)
124 {
125 cout << "BTV PV sorting selected, but no jet collection given!" << endl;
126 throw 0;
127 }
128
129 fItTrackInputArray->Reset();
130 while((candidate = static_cast<Candidate*>(fItTrackInputArray->Next())))
131 {
132 if (candidate->Momentum.Pt () < 1.0)
133 continue;
134 if (candidate->ClusterIndex < 0)
135 continue;
136 TLorentzVector p (candidate->Momentum.Px (), candidate->Momentum.Py (), candidate->Momentum.Pz (), candidate->Momentum.E ());
137 bool isInJet = false;
138
139 fItJetInputArray->Reset();
140 while((jetCandidate = static_cast<Candidate*>(fItJetInputArray->Next())))
141 {
142 if (jetCandidate->Momentum.Pt () < 30.0)
143 continue;
144 TLorentzVector q (jetCandidate->Momentum.Px (), jetCandidate->Momentum.Py (), jetCandidate->Momentum.Pz (), jetCandidate->Momentum.E ());
145
146 if (p.DeltaR (q) > 0.4)
147 continue;
148 isInJet = true;
149 break;
150 }
151 if (!isInJet)
152 continue;
153
154 clusterIDToSumPT2.at (candidate->ClusterIndex) += candidate->Momentum.Pt () * candidate->Momentum.Pt ();
155 }
156
157 for (const auto &clusterID : clusterIDToSumPT2)
158 sortedClusterIDs.push_back (make_pair (clusterID.first, clusterID.second));
159 sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondDescending);
160 }
161 else if (fMethod == "GenClosest")
162 {
163 if (!fBeamSpotInputArray)
164 {
165 cout << "GenClosest PV sorting selected, but no beamspot collection given!" << endl;
166 throw 0;
167 }
168 if (!fBeamSpotInputArray->GetSize ())
169 {
170 cout << "Beamspot collection is empty!" << endl;
171 throw 0;
172 }
173
174 beamSpotCandidate = (Candidate *) fBeamSpotInputArray->At (0);
175 for (int iCluster = 0; iCluster < fInputArray->GetEntries (); iCluster++)
176 {
177 const Candidate &cluster = *((Candidate *) fInputArray->At (iCluster));
178 sortedClusterIDs.push_back (make_pair (cluster.ClusterIndex, fabs (cluster.Position.Z () - beamSpotCandidate->Position.Z ())));
179 }
180 sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondAscending);
181 }
182 else if (fMethod == "GenBest")
183 {
184 fItTrackInputArray->Reset();
185 while((candidate = static_cast<Candidate*>(fItTrackInputArray->Next())))
186 {
187 if (candidate->IsPU)
188 continue;
189 for (const auto &clusterID : clusterIDToIndex)
190 {
191 if (candidate->ClusterIndex != clusterID.first)
192 continue;
193 clusterIDToSumPT2.at (clusterID.first) += candidate->Momentum.Pt () * candidate->Momentum.Pt ();
194 }
195 }
196
197 for (const auto &clusterID : clusterIDToSumPT2)
198 sortedClusterIDs.push_back (make_pair (clusterID.first, clusterID.second));
199 sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondDescending);
200 }
201 else
202 {
203 cout << "\"" << fMethod << "\" is not a valid sorting method!" << endl;
204 cout << "Valid methods are:" << endl;
205 cout << " BTV" << endl;
206 cout << " GenClosest" << endl;
207 cout << " GenBest" << endl;
208 throw 0;
209 }
210
211 for (const auto &clusterID : sortedClusterIDs)
212 {
213 Candidate *cluster = (Candidate *) fInputArray->At (clusterIDToIndex.at (clusterID.first));
214 if (fMethod == "BTV")
215 cluster->BTVSumPT2 = clusterID.second;
216 else if (fMethod == "GenClosest")
217 cluster->GenDeltaZ = clusterID.second;
218 else if (fMethod == "GenBest")
219 cluster->GenSumPT2 = clusterID.second;
220 fOutputArray->Add (cluster);
221 }
222}
223
224//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.