CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ProtonTaggerFilter.cc
Go to the documentation of this file.
1 
11 // Version: $Id: ProtonTaggerFilter.cc,v 1.3 2009/03/03 14:02:39 abdullin Exp $
12 
14 
22 
23 
24 //#include "CLHEP/Random/RandGaussQ.h"
25 
26 
27 #include <iostream>
28 #include <list>
29 
30 
34 {
35  std::cout << "ProtonTaggerFilter: Initializing ..." << std::endl;
36 
37  // ... get parameters
38 
39  beam1mode = p.getParameter<unsigned int>("beam1mode");
40  beam2mode = p.getParameter<unsigned int>("beam2mode");
41 
42  beamCombiningMode = p.getParameter<unsigned int>("beamCombiningMode");
43 
44  switch (beam1mode)
45  {
46  case 0: std::cout << "Option chosen for beam 1: ingnore" << std::endl; break;
47  case 1: std::cout << "Option chosen for beam 1: 420" << std::endl; break;
48  case 2: std::cout << "Option chosen for beam 1: 220" << std::endl; break;
49  case 3: std::cout << "Option chosen for beam 1: 420 and 220" << std::endl; break;
50  case 4: std::cout << "Option chosen for beam 1: 420 or 220" << std::endl; break;
51  default: throw cms::Exception("FastSimulation/ProtonTaggers")
52  << "Error: beam1mode cannot be " << beam1mode;
53  }
54 
55  switch (beam2mode)
56  {
57  case 0: std::cout << "Option chosen for beam 2: ingnore" << std::endl; break;
58  case 1: std::cout << "Option chosen for beam 2: 420" << std::endl; break;
59  case 2: std::cout << "Option chosen for beam 2: 220" << std::endl; break;
60  case 3: std::cout << "Option chosen for beam 2: 420 and 220" << std::endl; break;
61  case 4: std::cout << "Option chosen for beam 2: 420 or 220" << std::endl; break;
62  default: throw cms::Exception("FastSimulation/ProtonTaggers")
63  << "Error: beam2mode cannot be " << beam2mode;
64  }
65 
66  switch (beamCombiningMode)
67  {
68  case 1: std::cout << "Option chosen: one proton is sufficient" << std::endl; break;
69  case 2: std::cout << "Option chosen: two protons should be tagged" << std::endl; break;
70  case 3: std::cout << "Option chosen: two protons should be tagged as 220+220 or 420+420" << std::endl; break;
71  case 4: std::cout << "Option chosen: two protons should be tagged as 220+420 or 420+220" << std::endl; break;
72  default: throw cms::Exception("FastSimulation/ProtonTaggers")
73  << "Error: beamCombiningMode cannot be " << beamCombiningMode;
74  }
75 
76  if (((beam1mode != 4) || (beam2mode != 4)) && (beamCombiningMode > 2)) {
77  std::cerr << "Warning: beamCombiningMode = " << beamCombiningMode
78  << " only makes sence with beam1mode = beam2mode = 4" << std::endl;
79  }
80 
81  if (((beam1mode == 0) || (beam2mode == 0)) && (beamCombiningMode > 1)) {
82  std::cerr << "Warning: You ask for 2 protons while one of the beams is set to ignore" << std::endl;
83  }
84 
85  if ((beam1mode == 0) && (beam2mode == 0)) {
86  std::cerr << "Warning: Both beams are set to ignore." << std::endl;
87  }
88 
89  std::cout << "ProtonTaggerFilter: Initialized" << std::endl;
90 }
91 
92 
96 
97 
101 {
102  std::cout << "ProtonTaggerFilter: Getting ready ..." << std::endl;
103 
104  edm::FileInPath myDataFile("FastSimulation/ForwardDetectors/data/acceptance_420_220.root");
105  std::string fullPath = myDataFile.fullPath();
106 
107  std::cout << "Opening " << fullPath << std::endl;
108  TFile f(fullPath.c_str());
109 
110  if (f.Get("description") != NULL)
111  std::cout << "Description found: " << f.Get("description")->GetTitle() << std::endl;
112 
113  std::cout << "Reading acceptance tables @#@#%@$%@$#%@%" << std::endl;
114 
115  helper420beam1.Init(f, "a420");
116  helper420beam2.Init(f, "a420_b2");
117 
118  helper220beam1.Init(f, "a220");
119  helper220beam2.Init(f, "a220_b2");
120 
121  helper420a220beam1.Init(f, "a420a220");
122  helper420a220beam2.Init(f, "a420a220_b2");
123 
124  f.Close();
125 
126  std::cout << "ProtonTaggerFilter: Ready" << std::endl;
127 }
128 
129 
133 
134 
138 {
139  // ... get generated event
140 
142  iEvent.getByLabel("generatorSmeared",evtSource);
143  const HepMC::GenEvent* genEvent = evtSource->GetEvent();
144 
145  //std::cout << "event contains " << genEvent->particles_size() << " particles " << std::endl;
146  if (genEvent->particles_empty()) {
147  std::cout << "empty source event" << std::endl;
148  return false;
149  }
150 
151  // ... get pileup event
152 
153  edm::Handle<edm::HepMCProduct> pileUpSource;
154  const HepMC::GenEvent* pileUpEvent = 0;
155  bool isPileUp = true;
156 
157  bool isProduct = iEvent.getByLabel("famosPileUp", "PileUpEvents", pileUpSource);
158 
159  if (isProduct) {
160  pileUpEvent = pileUpSource->GetEvent();
161  //std::cout << "got pileup" << std::endl;
162  } else {
163  isPileUp = false;
164  //std::cout << "no pileup in the event" << std::endl;
165  }
166 
167  //if (isPileUp) std::cout << "event contains " << pileUpEvent->particles_size() << " pileup particles " << std::endl;
168 
169  // ... some constants
170 
171  const double mp = 0.938272029; // just a proton mass
172  const double Ebeam = 7000.0; // beam energy - would be better to read from parameters
173  const float pzCut = 2500; // ignore particles with less than |Pz| < pzCut
174  const float acceptThreshold = 0.5; // proton will be "accepted" if the value of acceptance is above this threshold
175 
176  // ... loop over particles, find the most energetic proton in either direction
177 
178  std::list<HepMC::GenParticle*> veryForwardParicles;
179 
180  for ( HepMC::GenEvent::particle_const_iterator piter = genEvent->particles_begin();
181  piter != genEvent->particles_end();
182  ++piter )
183  {
184 
185  HepMC::GenParticle* p = *piter;
186 
187  float pz = p->momentum().pz();
188  if (((pz > pzCut) || (pz < -pzCut)) && ((p->status() == 0) || (p->status() == 1))) {
189  veryForwardParicles.push_back(p);
190  //std::cout << "pdgid: " << p->pdg_id() << " status: " << p->status() << std::endl;
191  }
192  }
193 
194  //std::cout << "# generated forward particles : " << veryForwardParicles.size() << std::endl;
195 
196  if ( isPileUp ) {
197 
198  //std::cout << "Adding pileup " << std::endl;
199 
200  for ( HepMC::GenEvent::particle_const_iterator piter = pileUpEvent->particles_begin();
201  piter != pileUpEvent->particles_end();
202  ++piter )
203  {
204 
205  HepMC::GenParticle* p = *piter;
206 
207  float pz = p->momentum().pz();
208  if (((pz > pzCut) || (pz < -pzCut)) && ((p->status() == 0) || (p->status() == 1))) {
209  veryForwardParicles.push_back(p);
210  //std::cout << "pdgid: " << p->pdg_id() << " status: " << p->status() << std::endl;
211  }
212  }
213  }
214 
215  //std::cout << "# forward particles to be tried: " << veryForwardParicles.size() << std::endl;
216 
217 
218  // ... return false if no forward protons found
219 
220  if (veryForwardParicles.empty()) return false;
221 
222 
223  // ... set all acceptances to zero
224 
225  float acc420b1, acc220b1, acc420and220b1, acc420or220b1; // beam 1 (clockwise)
226  float acc420b2, acc220b2, acc420and220b2, acc420or220b2; // beam 2 (anti-clockwise)
227 
228  acc420b1 = acc220b1 = acc420and220b1 = acc420or220b1 = 0;
229  acc420b2 = acc220b2 = acc420and220b2 = acc420or220b2 = 0;
230 
231  int nP1at220m = 0;
232  int nP1at420m = 0;
233 
234  int nP2at220m = 0;
235  int nP2at420m = 0;
236 
237  // ... loop over (pre-selected) forward particles
238 
239  for (std::list<HepMC::GenParticle*>::const_iterator part = veryForwardParicles.begin();
240  part != veryForwardParicles.end();
241  part++)
242  {
243 
245 
246  float pz = p->momentum().pz();
247  float pt = p->momentum().perp();
248  float phi = p->momentum().phi();;
249 
250  if ((pz > Ebeam) || (pz < -Ebeam)) continue;
251 
252  // ... compute kinimatical variable
253 
254  float xi = 1.0; // fractional momentum loss
255  if (pz>0) xi -= pz/Ebeam;
256  else xi += pz/Ebeam;
257 
258  double t = (-pt*pt - mp*mp*xi*xi) / (1-xi); // "t"
259 
260  //std::cout << " pdg_id: " << p->pdg_id() << " eta: " << p->momentum().eta() << " e: "
261  // << p->momentum().e() << std::endl;
262  //std::cout << "pz: " << pz << " pt: " << pt << " xi: " << xi
263  // << " t: " << t << " phi: " << phi << std::endl;
264 
265  if (xi<0.0) xi=-10.0;
266  if (xi>1.0) xi=10.0;
267 
268  //float rnd1 = RandGauss::shoot(0.,1.);
269  //float rnd2 = RandGauss::shoot(0.,1.);
270 
271  // ... get acceptance from tables: beam 1 (if it is not ignored)
272 
273  if ((pz > 0) && (beam1mode != 0)) {
274 
275  acc420b1 = helper420beam1.GetAcceptance(t, xi, phi);
276  acc220b1 = helper220beam1.GetAcceptance(t, xi, phi);
277  acc420and220b1 = helper420a220beam1.GetAcceptance(t, xi, phi);
278 
279  acc420or220b1 = acc420b1 + acc220b1 - acc420and220b1;
280 
281  //std::cout << "+acc420b1: " << acc420b1 << " acc220b1: " << acc220b1 << " acc420and220b1: " << acc420and220b1 << " acc420or220b1: " << acc420or220b1 << std::endl;
282 
283  bool res420and220 = (acc420and220b1 > acceptThreshold);
284  bool res420 = (acc420b1 > acceptThreshold);
285  bool res220 = (acc220b1 > acceptThreshold);
286 
287  if (res420and220) {nP1at220m++; nP1at420m++;}
288  else if (res420) nP1at420m++;
289  else if (res220) nP1at220m++;
290 
291  if ((p->pdg_id() != 2212) && (res220 || res420 || res420and220)) {
292  std::cout << " !!! P got proton 1 at 420 m: pz = " << pz << std::endl;
293  if (res220) std::cout << "got a particle with pid" << p->pdg_id() << " at 220 m along beam 1, pz = " << pz << std::endl;
294  if (res420) std::cout << "got a particle with pid" << p->pdg_id() << " at 420 m along beam 1, pz = " << pz << std::endl;
295  if (res420and220) std::cout << "got a particle with pid" << p->pdg_id() << " at 220 m & 420 m along beam 1, pz = " << pz << std::endl;
296  }
297  }
298 
299  // ... the same for beam 2
300 
301  if ((pz < 0) && (beam2mode != 0)) {
302 
303  acc420b2 = helper420beam2.GetAcceptance(t, xi, phi);
304  acc220b2 = helper220beam2.GetAcceptance(t, xi, phi);
305  acc420and220b2 = helper420a220beam2.GetAcceptance(t, xi, phi);
306 
307  acc420or220b2 = acc420b2 + acc220b2 - acc420and220b2;
308 
309  //std::cout << "+acc420b2: " << acc420b2 << " acc220b2: " << acc220b2 << " acc420and220b2: " << acc420and220b2 << " acc420or220b2: " << acc420or220b2 << std::endl;
310 
311  bool res420and220 = (acc420and220b2 > acceptThreshold);
312  bool res420 = (acc420b2 > acceptThreshold);
313  bool res220 = (acc220b2 > acceptThreshold);
314 
315  if (res420and220) {nP2at220m++; nP2at420m++;}
316  else if (res420) nP2at420m++;
317  else if (res220) nP2at220m++;
318 
319  if ((p->pdg_id() != 2212) && (res220 || res420 || res420and220)) {
320  std::cout << " !!! P got proton 1 at 420 m: pz = " << pz << std::endl;
321  if (res220) std::cout << "got a particle with pid" << p->pdg_id() << " at 220 m along beam 2, pz = " << pz << std::endl;
322  if (res420) std::cout << "got a particle with pid" << p->pdg_id() << " at 420 m along beam 2, pz = " << pz << std::endl;
323  if (res420and220) std::cout << "got a particle with pid" << p->pdg_id() << " at 220 m & 420 m along beam 2, pz = " << pz << std::endl;
324  }
325  }
326  }
327 
328  // ... boolean result for each detector
329 
330  bool p1at220m = (nP1at220m>0) ? true : false;
331  bool p1at420m = (nP1at420m>0) ? true : false;
332  bool p2at220m = (nP2at220m>0) ? true : false;
333  bool p2at420m = (nP2at420m>0) ? true : false;
334 
335  if ((nP1at220m > 1) && (beam1mode!=1)) std::cout << " !!! " << nP1at220m << " proton(s) from beam 1 at 220 m" << std::endl;
336  if ((nP1at420m > 1) && (beam1mode!=2)) std::cout << " !!! " << nP1at420m << " proton(s) from beam 1 at 420 m" << std::endl;
337  if ((nP2at220m > 1) && (beam2mode!=1)) std::cout << " !!! " << nP2at220m << " proton(s) from beam 2 at 220 m" << std::endl;
338  if ((nP2at420m > 1) && (beam2mode!=2)) std::cout << " !!! " << nP2at420m << " proton(s) from beam 2 at 420 m" << std::endl;
339 
340 
341  // ... make a decision based on requested filter configuration
342 
343  bool p1accepted = false;
344  bool p2accepted = false;
345 
346  if ((beam1mode==1) && p1at420m ) p1accepted = true;
347  if ((beam1mode==2) && p1at220m ) p1accepted = true;
348  if ((beam1mode==3) && p1at220m && p1at420m) p1accepted = true;
349  if ((beam1mode==4) && (p1at220m || p1at420m)) p1accepted = true;
350 
351  if ((beam2mode==1) && p2at420m ) p2accepted = true;
352  if ((beam2mode==2) && p2at220m ) p2accepted = true;
353  if ((beam2mode==3) && p2at220m && p2at420m ) p2accepted = true;
354  if ((beam2mode==4) && (p2at220m || p2at420m)) p2accepted = true;
355 
356  //if (p1accepted) std::cout << "proton 1 accepted" << std::endl;
357  //if (p2accepted) std::cout << "proton 2 accepted" << std::endl;
358 
359  switch (beamCombiningMode) {
360 
361  case 1: // ... either of two protons
362  if (p1accepted || p2accepted) return true; else return false;
363 
364  case 2: // ... both protons
365  if (p1accepted && p2accepted) return true; else return false;
366 
367  case 3: // ... 220+220 or 420+420
368  if ((p1at220m && p2at220m) || (p1at420m && p2at420m)) return true; else return false;
369 
370  case 4: // ... 220+420 or 420+220
371  if ((p1at220m && p2at420m) || (p1at420m && p2at220m)) return true; else return false;
372 
373  }
374 
375  return false;
376 }
377 
378 // ... define CMSSW module
379 
T getParameter(std::string const &) const
void Init(TFile &, const std::string)
Get acceptance tables from root file.
virtual bool filter(edm::Event &e, const edm::EventSetup &c)
decide if the event is accepted by the proton taggers
virtual ~ProtonTaggerFilter()
empty destructor
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
AcceptanceTableHelper helper420beam1
Objects which actually compute the acceptance (one per detector or combination of detectors) ...
#define NULL
Definition: scimark2.h:8
virtual void endJob()
endjob function of the EDFilter
AcceptanceTableHelper helper420a220beam1
AcceptanceTableHelper helper420a220beam2
int iEvent
Definition: GenABIO.cc:230
unsigned int beamCombiningMode
choose how to combine data from the two beams (ask for 1/2 proton)
AcceptanceTableHelper helper220beam2
double f[11][100]
tuple genEvent
Definition: MCTruth.py:33
ProtonTaggerFilter(edm::ParameterSet const &p)
default constructor
unsigned int beam1mode
choose which of the detectors (FP420/TOTEM/both) will be used for beam 1
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
AcceptanceTableHelper helper220beam1
unsigned int beam2mode
choose which of the detectors (FP420/TOTEM/both) will be used for beam 2
part
Definition: HCALResponse.h:20
AcceptanceTableHelper helper420beam2
virtual void beginJob()
startup function of the EDFilter
tuple cout
Definition: gather_cfg.py:121
std::string fullPath() const
Definition: FileInPath.cc:165
float GetAcceptance(float, float, float)
Acceptance as a function of t, xi and phi.
Fast simulation of near-beam detector acceptance.