ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/dasu/UltraFastSim/UltraFastSim.cc
Revision: 1.1.1.1 (vendor branch)
Committed: Wed Feb 9 05:39:29 2011 UTC (14 years, 3 months ago) by dasu
Content type: text/plain
Branch: v0
CVS Tags: r0
Changes since 1.1: +0 -0 lines
Log Message:
First version 

File Contents

# Content
1 #include "UltraFastSim.h"
2
3 #include "Event.h"
4
5 #include "Analysis.h"
6
7 #include "fastjet/PseudoJet.hh"
8 #include "fastjet/JetDefinition.hh"
9 #include "fastjet/ClusterSequence.hh"
10
11 using namespace Pythia8;
12 using namespace fastjet;
13 using namespace std;
14
15 UltraFastSim::UltraFastSim(Rndm* r) : rndmPtr(r),
16 jetDefPtr(0),
17 cs(0),
18 trackerResolution(0.001),
19 ecalResolution(0.01),
20 ecalConstantTerm(0.01),
21 hcalResolution(1.),
22 hcalConstantTerm(0.1)
23 {
24 jetDefPtr = new JetDefinition(antikt_algorithm, 0.5);
25 return;
26 }
27
28 bool UltraFastSim::run(Event &event) {
29
30 // Clear the previous event
31
32 clear();
33
34 // Select particles of interest for later analysis
35
36 for (int i = 0; i < event.size(); ++i) {
37 Particle& particle = event[i];
38
39 // Select generated b quarks in detector acceptance
40
41 if(abs(particle.id()) == 5 && particle.pT() > 10. && abs(particle.eta()) < 2.5) {
42 selectedBQuarks.push_back(PseudoJet(particle.px(), particle.py(), particle.pz(), particle.e()));
43 }
44
45 // Select generated taus in detector acceptance
46
47 if(abs(particle.id()) == 15 && particle.pT() > 10. && abs(particle.eta()) < 2.5) {
48 selectedTaus.push_back(PseudoJet(particle.px(), particle.py(), particle.pz(), particle.e()));
49 }
50
51 // Consider only particles with good status
52 if(particle.status() > 0) {
53 // Consider only visible particles
54 if(particle.isVisible()) {
55 // Ignore soft particles and those outside the detector
56 if(particle.pT() > 1.0 && abs(particle.eta()) < 5.0) {
57
58 // Select electrons within detector acceptance
59 if(abs(particle.id()) == 11 && particle.pT() > 10. && abs(particle.eta()) < 2.5) {
60 selectedElectrons.push_back(PseudoJet(particle.px(), particle.py(), particle.pz(), particle.e()));
61 }
62
63 // Select muons within detector acceptance
64 if(abs(particle.id()) == 13 && particle.pT() > 10. && abs(particle.eta()) < 2.5) {
65 selectedMuons.push_back(PseudoJet(particle.px(), particle.py(), particle.pz(), particle.e()));
66 }
67
68 // Select other charged tracks and smear them
69 // Select photons and smear them
70 // Select other neutral particles and smear them
71
72 // Select for making jets using fastjet
73 selectedParticles.push_back(PseudoJet(particle.px(), particle.py(), particle.pz(), particle.e()));
74
75 }
76 }
77 }
78 }
79
80 makeTaus();
81 makeJets();
82 makeBJets();
83
84 cout << "Number of Particles = " << selectedParticles.size() << endl;
85 cout << "Number of Gen Elecs = " << selectedElectrons.size() << endl;
86 cout << "Number of Gen Muons = " << selectedMuons.size() << endl;
87 cout << "Number of Gen Taus = " << selectedTaus.size() << endl;
88 cout << "Number of b Quarks = " << selectedBQuarks.size() << endl;
89 cout << "Number of Jets = " << sortedJets.size() << endl;
90 cout << "Number of bJets = " << bJets.size() << endl;
91 cout << "Number of taus = " << taus.size() << endl;
92 cout << endl;
93
94 return true;
95
96 }
97
98 void UltraFastSim::clear() {
99 delete cs;
100 selectedParticles.clear();
101 selectedElectrons.clear();
102 selectedMuons.clear();
103 selectedTaus.clear();
104 selectedBQuarks.clear();
105 jets.clear();
106 sortedJets.clear();
107 bJets.clear();
108 taus.clear();
109 }
110
111 void UltraFastSim::makeJets() {
112 cs = new ClusterSequence(selectedParticles, *jetDefPtr);
113 double ptmin = 15.0;
114 jets = cs->inclusive_jets(ptmin);
115 sortedJets = sorted_by_pt(jets);
116 }
117
118 void UltraFastSim::makeBJets() {
119 for (unsigned int i = 0; i < sortedJets.size(); i++) {
120 if(abs(sortedJets[i].rap()) < 2.5) {
121 for(unsigned int j = 0; j < selectedBQuarks.size(); j++) {
122 float dRap = fabs(sortedJets[i].rap() - selectedBQuarks[j].rap());
123 float dPhi = fabs(sortedJets[i].phi() - selectedBQuarks[j].phi());
124 float dR = sqrt(dRap*dRap + dPhi*dPhi);
125 if(dR < 0.3) {
126 bJets.push_back(sortedJets[i]);
127 break;
128 }
129 }
130 }
131 }
132 }
133
134 void UltraFastSim::makeTaus() {
135 for (unsigned int i = 0; i < selectedParticles.size(); i++) {
136 // Make sure that there is high PT track seed
137 if(selectedParticles[i].Et() > 5. && abs(selectedParticles[i].rap()) < 2.5) {
138 int nObjectsInInnerCone = 0;
139 float isolationEnergy = 0.;
140 PseudoJet tau(selectedParticles[i]);
141 for(unsigned int j = 0; j < selectedParticles.size(); j++) {
142 if(i != j) {
143 if(selectedParticles[i].Et() >= selectedParticles[j].Et()) {
144 float dRap = fabs(selectedParticles[i].rap() - selectedParticles[j].rap());
145 float dPhi = fabs(selectedParticles[i].phi() - selectedParticles[j].phi());
146 float dR = sqrt(dRap*dRap + dPhi*dPhi);
147 // Make sure that there are no more than five objects in the inner 0.3 cone
148 if(dR < 0.3) {
149 nObjectsInInnerCone++;
150 if(nObjectsInInnerCone > 5) break;
151 tau += selectedParticles[j];
152 }
153 else if(dR < 0.5) {
154 isolationEnergy += selectedParticles[j].e();
155 }
156 }
157 else {
158 break;
159 }
160 }
161 }
162 // Make sure that there is no more than 5% energy in the 0.5
163 if(tau.Et() > 15. && (isolationEnergy / tau.e()) < 0.05) {
164 taus.push_back(tau);
165 }
166 }
167 }
168 }