CMS 3D CMS Logo

CTPPSDirectProtonSimulation.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  * Laurent Forthomme
5  ****************************************************************************/
6 
15 
19 
23 
29 
32 
35 
38 
41 
47 
50 
51 #include "CLHEP/Random/RandGauss.h"
52 #include "CLHEP/Units/GlobalPhysicalConstants.h"
53 
54 #include <unordered_map>
55 
56 #include "TMath.h"
57 #include "TMatrixD.h"
58 #include "TVectorD.h"
59 #include "TF1.h"
60 #include "TF2.h"
61 #include "TFile.h"
62 #include "CLHEP/Random/RandFlat.h"
63 
64 //----------------------------------------------------------------------------------------------------
65 
67 public:
69  ~CTPPSDirectProtonSimulation() override {}
70 
72 
73 private:
74  void produce(edm::Event &, const edm::EventSetup &) override;
75 
76  void processProton(const HepMC::GenVertex *in_vtx,
77  const HepMC::GenParticle *in_trk,
78  const CTPPSGeometry &geometry,
79  const LHCInfo &lhcInfo,
80  const CTPPSBeamParameters &beamParameters,
81  const PPSPixelTopology &ppt,
83  CLHEP::HepRandomEngine *rndEngine,
84 
85  std::vector<CTPPSLocalTrackLite> &out_tracks,
86 
87  edm::DetSetVector<TotemRPRecHit> &out_strip_hits,
89  edm::DetSetVector<CTPPSDiamondRecHit> &out_diamond_hits,
90 
91  std::map<int, edm::DetSetVector<TotemRPRecHit>> &out_strip_hits_per_particle,
92  std::map<int, edm::DetSetVector<CTPPSPixelRecHit>> &out_pixel_hits_per_particle,
93  std::map<int, edm::DetSetVector<CTPPSDiamondRecHit>> &out_diamond_hits_per_particle) const;
94 
95  // ------------ config file parameters ------------
96 
97  // conditions
104 
106 
107  // input
109 
110  // flags what output to be produced
112  bool produceRecHits_;
113 
114  // settings of LHC aperture limitations (high xi)
116  std::unique_ptr<TF2> empiricalAperture45_;
117  std::unique_ptr<TF2> empiricalAperture56_;
118 
119  // efficiency flags
124 
125  // efficiency maps
126  std::map<unsigned int, std::unique_ptr<TH2F>> efficiencyMapsPerRP_;
127  std::map<unsigned int, std::unique_ptr<TH2F>> efficiencyMapsPerPlane_;
128 
129  // other parameters
131  bool roundToPitch_;
133 
134  double pitchStrips_;
135  double insensitiveMarginStrips_;
136 
138  double pitchPixelsVer_;
139 
140  unsigned int verbosity_;
141 
143 
144  // ------------ internal parameters ------------
145 
147  double stripZeroPosition_;
148 };
149 
150 //----------------------------------------------------------------------------------------------------
151 
153  : tokenLHCInfo_(esConsumes(edm::ESInputTag{"", iConfig.getParameter<std::string>("lhcInfoLabel")})),
154  tokenBeamParameters_(esConsumes()),
155  pixelTopologyToken_(esConsumes()),
156  tokenOpticalFunctions_(esConsumes(edm::ESInputTag{"", iConfig.getParameter<std::string>("opticsLabel")})),
157  tokenGeometry_(esConsumes()),
158  tokenDirectSimuData_(esConsumes()),
159 
160  hepMCToken_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("hepMCTag"))),
161 
162  produceScoringPlaneHits_(iConfig.getParameter<bool>("produceScoringPlaneHits")),
163  produceRecHits_(iConfig.getParameter<bool>("produceRecHits")),
164 
165  useEmpiricalApertures_(iConfig.getParameter<bool>("useEmpiricalApertures")),
166 
167  useTrackingEfficiencyPerRP_(iConfig.getParameter<bool>("useTrackingEfficiencyPerRP")),
168  useTimingEfficiencyPerRP_(iConfig.getParameter<bool>("useTimingEfficiencyPerRP")),
169  useTrackingEfficiencyPerPlane_(iConfig.getParameter<bool>("useTrackingEfficiencyPerPlane")),
170  useTimingEfficiencyPerPlane_(iConfig.getParameter<bool>("useTimingEfficiencyPerPlane")),
171 
172  produceHitsRelativeToBeam_(iConfig.getParameter<bool>("produceHitsRelativeToBeam")),
173  roundToPitch_(iConfig.getParameter<bool>("roundToPitch")),
174  checkIsHit_(iConfig.getParameter<bool>("checkIsHit")),
175 
176  pitchStrips_(iConfig.getParameter<double>("pitchStrips")),
177  insensitiveMarginStrips_(iConfig.getParameter<double>("insensitiveMarginStrips")),
178 
179  pitchPixelsHor_(iConfig.getParameter<double>("pitchPixelsHor")),
180  pitchPixelsVer_(iConfig.getParameter<double>("pitchPixelsVer")),
181 
182  verbosity_(iConfig.getUntrackedParameter<unsigned int>("verbosity", 0)) {
183  if (produceScoringPlaneHits_)
184  produces<std::vector<CTPPSLocalTrackLite>>();
185 
186  if (produceRecHits_) {
187  produces<edm::DetSetVector<TotemRPRecHit>>();
188  produces<edm::DetSetVector<CTPPSDiamondRecHit>>();
189  produces<edm::DetSetVector<CTPPSPixelRecHit>>();
190 
191  produces<std::map<int, edm::DetSetVector<TotemRPRecHit>>>();
192  produces<std::map<int, edm::DetSetVector<CTPPSDiamondRecHit>>>();
193  produces<std::map<int, edm::DetSetVector<CTPPSPixelRecHit>>>();
194  }
195 
196  // check user input
197  if (useTrackingEfficiencyPerRP_ && useTrackingEfficiencyPerPlane_)
198  throw cms::Exception("PPS")
199  << "useTrackingEfficiencyPerRP and useTrackingEfficiencyPerPlane should not be simultaneously set true.";
200 
201  if (useTimingEfficiencyPerRP_ && useTimingEfficiencyPerPlane_)
202  throw cms::Exception("PPS")
203  << "useTimingEfficiencyPerRP and useTimingEfficiencyPerPlane should not be simultaneously set true.";
204 
205  // v position of strip 0
208 }
209 
210 //----------------------------------------------------------------------------------------------------
211 
214  desc.addUntracked<unsigned int>("verbosity", 0);
215 
216  desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
217  desc.add<std::string>("opticsLabel", "")->setComment("label of the optics records");
218 
219  desc.add<edm::InputTag>("hepMCTag", edm::InputTag("generator", "unsmeared"));
220 
221  desc.add<bool>("produceScoringPlaneHits", true);
222  desc.add<bool>("produceRecHits", true);
223 
224  desc.add<bool>("useEmpiricalApertures", true);
225 
226  desc.add<bool>("useTrackingEfficiencyPerRP", false);
227  desc.add<bool>("useTimingEfficiencyPerRP", false);
228  desc.add<bool>("useTrackingEfficiencyPerPlane", false);
229  desc.add<bool>("useTimingEfficiencyPerPlane", false);
230 
231  desc.add<bool>("produceHitsRelativeToBeam", true);
232  desc.add<bool>("roundToPitch", true);
233  desc.add<bool>("checkIsHit", true);
234  desc.add<double>("pitchStrips", 66.e-3); // in mm
235  desc.add<double>("insensitiveMarginStrips", 34.e-3); // in mm
236 
237  desc.add<double>("pitchPixelsHor", 100.e-3);
238  desc.add<double>("pitchPixelsVer", 150.e-3);
239 
240  descriptions.add("ctppsDirectProtonSimulation", desc);
241 }
242 
243 //----------------------------------------------------------------------------------------------------
244 
246  // get input
248  iEvent.getByToken(hepMCToken_, hepmc_prod);
249 
250  // get conditions
251  auto const &lhcInfo = iSetup.getData(tokenLHCInfo_);
252  auto const &beamParameters = iSetup.getData(tokenBeamParameters_);
253  auto const &ppt = iSetup.getData(pixelTopologyToken_);
254  auto const &opticalFunctions = iSetup.getData(tokenOpticalFunctions_);
255  auto const &geometry = iSetup.getData(tokenGeometry_);
256  auto const &directSimuData = iSetup.getData(tokenDirectSimuData_);
257 
258  if (directSimuDataRcdWatcher_.check(iSetup)) {
260  std::make_unique<TF1>(TF1("timeResolutionDiamonds45", directSimuData.getTimeResolutionDiamonds45().c_str()));
262  std::make_unique<TF1>(TF1("timeResolutionDiamonds56", directSimuData.getTimeResolutionDiamonds56().c_str()));
263 
265  std::make_unique<TF2>(TF2("empiricalAperture45", directSimuData.getEmpiricalAperture45().c_str()));
267  std::make_unique<TF2>(TF2("empiricalAperture56", directSimuData.getEmpiricalAperture56().c_str()));
268 
269  // load the efficiency maps
271  efficiencyMapsPerRP_ = directSimuData.loadEffeciencyHistogramsPerRP();
273  efficiencyMapsPerPlane_ = directSimuData.loadEffeciencyHistogramsPerPlane();
274  }
275 
276  // prepare outputs
277  std::unique_ptr<edm::DetSetVector<TotemRPRecHit>> pStripRecHits(new edm::DetSetVector<TotemRPRecHit>());
278  std::unique_ptr<edm::DetSetVector<CTPPSDiamondRecHit>> pDiamondRecHits(new edm::DetSetVector<CTPPSDiamondRecHit>());
279  std::unique_ptr<edm::DetSetVector<CTPPSPixelRecHit>> pPixelRecHits(new edm::DetSetVector<CTPPSPixelRecHit>());
280 
281  auto pStripRecHitsPerParticle = std::make_unique<std::map<int, edm::DetSetVector<TotemRPRecHit>>>();
282  auto pDiamondRecHitsPerParticle = std::make_unique<std::map<int, edm::DetSetVector<CTPPSDiamondRecHit>>>();
283  auto pPixelRecHitsPerParticle = std::make_unique<std::map<int, edm::DetSetVector<CTPPSPixelRecHit>>>();
284 
285  std::unique_ptr<std::vector<CTPPSLocalTrackLite>> pTracks(new std::vector<CTPPSLocalTrackLite>());
286 
287  // get random engine
289  CLHEP::HepRandomEngine *engine = &rng->getEngine(iEvent.streamID());
290 
291  // loop over event vertices
292  auto evt = hepmc_prod->GetEvent();
293  for (auto it_vtx = evt->vertices_begin(); it_vtx != evt->vertices_end(); ++it_vtx) {
294  auto vtx = *(it_vtx);
295 
296  // loop over outgoing particles
297  for (auto it_part = vtx->particles_out_const_begin(); it_part != vtx->particles_out_const_end(); ++it_part) {
298  auto part = *(it_part);
299 
300  // accept only stable protons
301  if (part->pdg_id() != 2212)
302  continue;
303 
304  if (part->status() != 1 && part->status() < 83)
305  continue;
306 
308  part,
309  geometry,
310  lhcInfo,
311  beamParameters,
312  ppt,
314  engine,
315  *pTracks,
316  *pStripRecHits,
317  *pPixelRecHits,
318  *pDiamondRecHits,
319  *pStripRecHitsPerParticle,
320  *pPixelRecHitsPerParticle,
321  *pDiamondRecHitsPerParticle);
322  }
323  }
324 
326  iEvent.put(std::move(pTracks));
327 
328  if (produceRecHits_) {
329  iEvent.put(std::move(pStripRecHits));
330  iEvent.put(std::move(pPixelRecHits));
331  iEvent.put(std::move(pDiamondRecHits));
332 
333  iEvent.put(std::move(pStripRecHitsPerParticle));
334  iEvent.put(std::move(pPixelRecHitsPerParticle));
335  iEvent.put(std::move(pDiamondRecHitsPerParticle));
336  }
337 }
338 
339 //----------------------------------------------------------------------------------------------------
340 
342  const HepMC::GenVertex *in_vtx,
343  const HepMC::GenParticle *in_trk,
344  const CTPPSGeometry &geometry,
345  const LHCInfo &lhcInfo,
346  const CTPPSBeamParameters &beamParameters,
347  const PPSPixelTopology &ppt,
349  CLHEP::HepRandomEngine *rndEngine,
350  std::vector<CTPPSLocalTrackLite> &out_tracks,
351  edm::DetSetVector<TotemRPRecHit> &out_strip_hits,
352  edm::DetSetVector<CTPPSPixelRecHit> &out_pixel_hits,
353  edm::DetSetVector<CTPPSDiamondRecHit> &out_diamond_hits,
354  std::map<int, edm::DetSetVector<TotemRPRecHit>> &out_strip_hits_per_particle,
355  std::map<int, edm::DetSetVector<CTPPSPixelRecHit>> &out_pixel_hits_per_particle,
356  std::map<int, edm::DetSetVector<CTPPSDiamondRecHit>> &out_diamond_hits_per_particle) const {
359 
360  std::stringstream ssLog;
361 
362  // vectors in CMS convention
363  const HepMC::FourVector &vtx_cms = in_vtx->position(); // in mm
364  const HepMC::FourVector &mom_cms = in_trk->momentum();
365 
366  // transformation to LHC/TOTEM convention
367  HepMC::FourVector vtx_lhc(-vtx_cms.x(), vtx_cms.y(), -vtx_cms.z(), vtx_cms.t());
368  HepMC::FourVector mom_lhc(-mom_cms.x(), mom_cms.y(), -mom_cms.z(), mom_cms.t());
369 
370  // determine the LHC arm and related parameters
371  unsigned int arm = 3;
372  double z_sign;
373  double beamMomentum = 0.;
374  double xangle = 0.;
375  const std::unique_ptr<TF2> *empiricalAperture;
376  if (mom_lhc.z() < 0) // sector 45
377  {
378  arm = 0;
379  z_sign = -1;
380  beamMomentum = beamParameters.getBeamMom45();
381  xangle = beamParameters.getHalfXangleX45();
382  empiricalAperture = &empiricalAperture45_;
383  } else { // sector 56
384  arm = 1;
385  z_sign = +1;
386  beamMomentum = beamParameters.getBeamMom56();
387  xangle = beamParameters.getHalfXangleX56();
388  empiricalAperture = &empiricalAperture56_;
389  }
390 
391  // calculate effective RP arrival time
392  // effective time mimics the timing calibration -> effective times are distributed about 0
393  // units:
394  // vertex: all components in mm
395  // c_light: in mm/ns
396  // time_eff: in ns
397  const double time_eff = (vtx_lhc.t() - z_sign * vtx_lhc.z()) / CLHEP::c_light;
398 
399  // calculate kinematics for optics parametrisation
400  const double p = mom_lhc.rho();
401  const double xi = 1. - p / beamMomentum;
402  const double th_x_phys = mom_lhc.x() / p;
403  const double th_y_phys = mom_lhc.y() / p;
404  const double vtx_lhc_eff_x = vtx_lhc.x() - vtx_lhc.z() * (mom_lhc.x() / mom_lhc.z() + xangle);
405  const double vtx_lhc_eff_y = vtx_lhc.y() - vtx_lhc.z() * (mom_lhc.y() / mom_lhc.z());
406 
407  if (verbosity_) {
408  ssLog << "simu: xi = " << xi << ", th_x_phys = " << th_x_phys << ", th_y_phys = " << th_y_phys
409  << ", vtx_lhc_eff_x = " << vtx_lhc_eff_x << ", vtx_lhc_eff_y = " << vtx_lhc_eff_y << std::endl;
410  }
411 
412  // check empirical aperture
414  const auto &xangle = lhcInfo.crossingAngle();
415  (*empiricalAperture)->SetParameter("xi", xi);
416  (*empiricalAperture)->SetParameter("xangle", xangle);
417  const double th_x_th = (*empiricalAperture)->EvalPar(nullptr);
418 
419  if (th_x_th > th_x_phys) {
420  if (verbosity_) {
421  ssLog << "stop because of empirical appertures";
422  edm::LogInfo("PPS") << ssLog.str();
423  }
424 
425  return;
426  }
427  }
428 
429  // transport the proton into each pot/scoring plane
430  for (const auto &ofp : opticalFunctions) {
431  CTPPSDetId rpId(ofp.first);
432  const unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
433 
434  // first check the arm
435  if (rpId.arm() != arm)
436  continue;
437 
438  if (verbosity_)
439  ssLog << " RP " << rpDecId << std::endl;
440 
441  // transport proton
443  vtx_lhc_eff_x * 1E-1, th_x_phys, vtx_lhc_eff_y * 1E-1, th_y_phys, xi}; // conversions: mm -> cm
445  ofp.second.transport(k_in, k_out, true);
446 
447  double b_x = k_out.x * 1E1, b_y = k_out.y * 1E1; // conversions: cm -> mm
448  double a_x = k_out.th_x, a_y = k_out.th_y;
449 
450  // if needed, subtract beam position and angle
452  // determine beam position
453  LHCInterpolatedOpticalFunctionsSet::Kinematics k_be_in = {0., 0., 0., 0., 0.};
455  ofp.second.transport(k_be_in, k_be_out, true);
456 
457  a_x -= k_be_out.th_x;
458  a_y -= k_be_out.th_y;
459  b_x -= k_be_out.x * 1E1;
460  b_y -= k_be_out.y * 1E1; // conversions: cm -> mm
461  }
462 
463  const double z_scoringPlane = ofp.second.getScoringPlaneZ() * 1E1; // conversion: cm --> mm
464 
465  if (verbosity_) {
466  ssLog << " proton transported: a_x = " << a_x << " rad, a_y = " << a_y << " rad, b_x = " << b_x
467  << " mm, b_y = " << b_y << " mm, z = " << z_scoringPlane << " mm" << std::endl;
468  }
469 
470  // RP type
471  const bool isTrackingRP =
472  (rpId.subdetId() == CTPPSDetId::sdTrackingStrip || rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
473  const bool isTimingRP = (rpId.subdetId() == CTPPSDetId::sdTimingDiamond);
474 
475  // apply per-RP efficiency
476  if ((useTimingEfficiencyPerRP_ && isTimingRP) || (useTrackingEfficiencyPerRP_ && isTrackingRP)) {
477  const auto it = efficiencyMapsPerRP_.find(rpId);
478 
479  if (it != efficiencyMapsPerRP_.end()) {
480  const double r = CLHEP::RandFlat::shoot(rndEngine, 0., 1.);
481  auto *effMap = it->second.get();
482  const double eff = effMap->GetBinContent(effMap->FindBin(b_x, b_y));
483  if (r > eff) {
484  if (verbosity_)
485  ssLog << " stop due to per-RP efficiency" << std::endl;
486  continue;
487  }
488  }
489  }
490 
491  // save scoring plane hit
493  out_tracks.emplace_back(
494  rpId.rawId(), b_x, 0., b_y, 0., 0., 0., 0., 0., 0., CTPPSpixelLocalTrackReconstructionInfo::invalid, 0, 0., 0.);
495 
496  // stop if rec hits are not to be produced
497  if (!produceRecHits_)
498  continue;
499 
500  // loop over all sensors in the RP
501  for (const auto &detIdInt : geometry.sensorsInRP(rpId)) {
502  CTPPSDetId detId(detIdInt);
503 
504  // determine the track impact point (in global coordinates)
505  // !! this assumes that local axes (1, 0, 0) and (0, 1, 0) describe the sensor surface
506  const auto &gl_o = geometry.localToGlobal(detId, CTPPSGeometry::Vector(0, 0, 0));
507  const auto &gl_a1 = geometry.localToGlobal(detId, CTPPSGeometry::Vector(1, 0, 0)) - gl_o;
508  const auto &gl_a2 = geometry.localToGlobal(detId, CTPPSGeometry::Vector(0, 1, 0)) - gl_o;
509 
510  const double gl_o_z = gl_o.z();
511 
512  TMatrixD A(3, 3);
513  TVectorD B(3);
514  A(0, 0) = a_x;
515  A(0, 1) = -gl_a1.x();
516  A(0, 2) = -gl_a2.x();
517  B(0) = gl_o.x() - b_x;
518  A(1, 0) = a_y;
519  A(1, 1) = -gl_a1.y();
520  A(1, 2) = -gl_a2.y();
521  B(1) = gl_o.y() - b_y;
522  A(2, 0) = z_sign;
523  A(2, 1) = -gl_a1.z();
524  A(2, 2) = -gl_a2.z();
525  B(2) = gl_o_z - z_scoringPlane;
526  TMatrixD Ai(3, 3);
527  Ai = A.Invert();
528  TVectorD P(3);
529  P = Ai * B;
530 
531  double ze = P(0);
532  const CTPPSGeometry::Vector h_glo(a_x * ze + b_x, a_y * ze + b_y, z_sign * ze + z_scoringPlane);
533 
534  // hit in local coordinates
535  auto h_loc = geometry.globalToLocal(detId, h_glo);
536 
537  if (verbosity_) {
538  ssLog << std::endl
539  << " de z = " << P(0) << " mm, p1 = " << P(1) << " mm, p2 = " << P(2) << " mm" << std::endl
540  << " h_glo: x = " << h_glo.x() << " mm, y = " << h_glo.y() << " mm, z = " << h_glo.z() << " mm"
541  << std::endl
542  << " h_loc: c1 = " << h_loc.x() << " mm, c2 = " << h_loc.y() << " mm, c3 = " << h_loc.z() << " mm"
543  << std::endl;
544  }
545 
546  // apply per-plane efficiency
547  if ((useTimingEfficiencyPerPlane_ && isTimingRP) || (useTrackingEfficiencyPerPlane_ && isTrackingRP)) {
548  const auto it = efficiencyMapsPerPlane_.find(detId);
549 
550  if (it != efficiencyMapsPerPlane_.end()) {
551  const double r = CLHEP::RandFlat::shoot(rndEngine, 0., 1.);
552  auto *effMap = it->second.get();
553  const double eff = effMap->GetBinContent(effMap->FindBin(h_glo.x(), h_glo.y()));
554  if (r > eff) {
555  if (verbosity_)
556  ssLog << " stop due to per-plane efficiency" << std::endl;
557  continue;
558  }
559  }
560  }
561 
562  // strips
563  if (detId.subdetId() == CTPPSDetId::sdTrackingStrip) {
564  double u = h_loc.x();
565  double v = h_loc.y();
566 
567  if (verbosity_ > 5)
568  ssLog << " u=" << u << ", v=" << v;
569 
570  // is it within detector?
572  if (verbosity_ > 5)
573  ssLog << " | no hit" << std::endl;
574  continue;
575  }
576 
577  // round the measurement
578  if (roundToPitch_) {
579  double m = stripZeroPosition_ - v;
580  signed int strip = (int)floor(m / pitchStrips_ + 0.5);
581 
583 
584  if (verbosity_ > 5)
585  ssLog << " | strip=" << strip;
586  }
587 
588  double sigma = pitchStrips_ / sqrt(12.);
589 
590  if (verbosity_ > 5)
591  ssLog << " | m=" << v << ", sigma=" << sigma << std::endl;
592 
593  edm::DetSet<TotemRPRecHit> &hits = out_strip_hits.find_or_insert(detId);
594  hits.emplace_back(v, sigma);
595 
596  edm::DetSet<TotemRPRecHit> &hits_per_particle =
597  out_strip_hits_per_particle[in_trk->barcode()].find_or_insert(detId);
598  hits_per_particle.emplace_back(v, sigma);
599  }
600 
601  // diamonds
602  if (detId.subdetId() == CTPPSDetId::sdTimingDiamond) {
603  CTPPSDiamondDetId diamondDetId(detIdInt);
604 
605  // check acceptance
606  const auto *dg = geometry.sensor(detIdInt);
607  const auto &diamondDimensions = dg->getDiamondDimensions();
608  const auto x_half_width = diamondDimensions.xHalfWidth;
609  const auto y_half_width = diamondDimensions.yHalfWidth;
610  const auto z_half_width = diamondDimensions.zHalfWidth;
611 
612  if (h_loc.x() < -x_half_width || h_loc.x() > +x_half_width || h_loc.y() < -y_half_width ||
613  h_loc.y() > +y_half_width)
614  continue;
615 
616  // timing information
617  const double time_resolution = (diamondDetId.arm() == 0) ? timeResolutionDiamonds45_->Eval(h_glo.x())
618  : timeResolutionDiamonds56_->Eval(h_glo.x());
619 
620  const double t0 = time_eff + CLHEP::RandGauss::shoot(rndEngine, 0., time_resolution);
621  const double tot = 1.23456;
622  const double ch_t_precis = time_resolution;
623  const int time_slice = 0;
624 
625  // build rec hit
626  const bool multiHit = false;
627 
628  CTPPSDiamondRecHit rc(gl_o.x(),
629  2. * x_half_width,
630  gl_o.y(),
631  2. * y_half_width,
632  gl_o_z,
633  2. * z_half_width,
634  t0,
635  tot,
636  ch_t_precis,
637  time_slice,
638  HPTDCErrorFlags(),
639  multiHit);
640 
641  edm::DetSet<CTPPSDiamondRecHit> &hits = out_diamond_hits.find_or_insert(detId);
642  hits.push_back(rc);
643 
644  edm::DetSet<CTPPSDiamondRecHit> &hits_per_particle =
645  out_diamond_hits_per_particle[in_trk->barcode()].find_or_insert(detId);
646  hits_per_particle.push_back(rc);
647  }
648 
649  // pixels
650  if (detId.subdetId() == CTPPSDetId::sdTrackingPixel) {
651  if (verbosity_) {
652  CTPPSPixelDetId pixelDetId(detIdInt);
653  ssLog << " pixel plane " << pixelDetId.plane() << ": local hit x = " << h_loc.x()
654  << " mm, y = " << h_loc.y() << " mm, z = " << h_loc.z() << " mm" << std::endl;
655  }
656 
657  bool module3By2 = (geometry.sensor(detIdInt)->sensorType() != DDD_CTPPS_PIXELS_SENSOR_TYPE_2x2);
658  if (checkIsHit_ && !ppt.isPixelHit(h_loc.x(), h_loc.y(), module3By2))
659  continue;
660 
661  if (roundToPitch_) {
662  h_loc.SetX(pitchPixelsHor_ * floor(h_loc.x() / pitchPixelsHor_ + 0.5));
663  h_loc.SetY(pitchPixelsVer_ * floor(h_loc.y() / pitchPixelsVer_ + 0.5));
664  }
665 
666  if (verbosity_ > 5)
667  ssLog << " hit accepted: m1 = " << h_loc.x() << " mm, m2 = " << h_loc.y() << " mm" << std::endl;
668 
669  const double sigmaHor = pitchPixelsHor_ / sqrt(12.);
670  const double sigmaVer = pitchPixelsVer_ / sqrt(12.);
671 
672  const LocalPoint lp(h_loc.x(), h_loc.y(), h_loc.z());
673  const LocalError le(sigmaHor, 0., sigmaVer);
674 
675  edm::DetSet<CTPPSPixelRecHit> &hits = out_pixel_hits.find_or_insert(detId);
676  hits.emplace_back(lp, le);
677 
678  edm::DetSet<CTPPSPixelRecHit> &hits_per_particle =
679  out_pixel_hits_per_particle[in_trk->barcode()].find_or_insert(detId);
680  hits_per_particle.emplace_back(lp, le);
681  }
682  }
683  }
684 
685  if (verbosity_)
686  edm::LogInfo("PPS") << ssLog.str();
687 }
688 
689 //----------------------------------------------------------------------------------------------------
690 
CTPPSBeamParameters
Definition: CTPPSBeamParameters.h:22
edm::ESWatcher::check
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
edm::DetSet::push_back
void push_back(const T &t)
Definition: DetSet.h:66
edm::DetSetVector< TotemRPRecHit >
CTPPSDirectProtonSimulation::~CTPPSDirectProtonSimulation
~CTPPSDirectProtonSimulation() override
Definition: CTPPSDirectProtonSimulation.cc:71
PPSPixelTopology.h
RPTopology::y_width_
static const double y_width_
Definition: RPTopology.h:60
edm::RandomNumberGenerator::getEngine
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
CTPPSBeamParameters.h
CTPPSDirectProtonSimulation::useTrackingEfficiencyPerRP_
bool useTrackingEfficiencyPerRP_
Definition: CTPPSDirectProtonSimulation.cc:122
edm::ESInputTag
Definition: ESInputTag.h:87
CTPPSBeamParametersRcd.h
CTPPSGeometry
The manager class for TOTEM RP geometry.
Definition: CTPPSGeometry.h:29
LHCInterpolatedOpticalFunctionsSet::Kinematics
proton kinematics description
Definition: LHCInterpolatedOpticalFunctionsSet.h:28
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
edm::ESWatcher< PPSDirectSimulationDataRcd >
ESHandle.h
ESInputTag
edm::DetSet
Definition: DetSet.h:23
CTPPSDirectProtonSimulation::useTrackingEfficiencyPerPlane_
bool useTrackingEfficiencyPerPlane_
Definition: CTPPSDirectProtonSimulation.cc:124
CTPPSDirectProtonSimulation::timeResolutionDiamonds56_
std::unique_ptr< TF1 > timeResolutionDiamonds56_
Definition: CTPPSDirectProtonSimulation.cc:144
edm::EDGetTokenT< edm::HepMCProduct >
CTPPSDirectProtonSimulation::insensitiveMarginStrips_
double insensitiveMarginStrips_
size of insensitive margin at sensor's edge facing the beam, in mm
Definition: CTPPSDirectProtonSimulation.cc:137
CTPPSDiamondRecHit
Reconstructed hit in diamond detectors.
Definition: CTPPSDiamondRecHit.h:16
edm
HLT enums.
Definition: AlignableModifier.h:19
RandomNumberGenerator.h
digitizers_cfi.strip
strip
Definition: digitizers_cfi.py:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
CTPPSBeamParameters::getHalfXangleX56
double getHalfXangleX56() const
Definition: CTPPSBeamParameters.cc:70
geometry
Definition: geometry.py:1
CTPPSPixelDetId.h
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89285
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
LHCInfo.h
CTPPSDirectProtonSimulation
Definition: CTPPSDirectProtonSimulation.cc:65
CTPPSBeamParameters::getBeamMom56
double getBeamMom56() const
Definition: CTPPSBeamParameters.cc:56
LHCInfo
Definition: LHCInfo.h:12
EDProducer.h
CTPPSBeamParameters::getHalfXangleX45
double getHalfXangleX45() const
Definition: CTPPSBeamParameters.cc:68
CTPPSDirectProtonSimulation::efficiencyMapsPerRP_
std::map< unsigned int, std::unique_ptr< TH2F > > efficiencyMapsPerRP_
Definition: CTPPSDirectProtonSimulation.cc:128
CTPPSLocalTrackLite.h
TotemRPRecHit.h
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
findQualityFiles.v
v
Definition: findQualityFiles.py:179
CTPPSpixelLocalTrackReconstructionInfo::invalid
DDD_CTPPS_PIXELS_SENSOR_TYPE_2x2
const std::string DDD_CTPPS_PIXELS_SENSOR_TYPE_2x2
Definition: CTPPSDDDNames.h:16
edm::Handle< edm::HepMCProduct >
CTPPSPixelRecHit.h
VeryForwardMisalignedGeometryRecord.h
LHCInterpolatedOpticalFunctionsSetCollection
Definition: LHCInterpolatedOpticalFunctionsSetCollection.h:10
CTPPSDetId::sdTrackingStrip
Definition: CTPPSDetId.h:44
CTPPSGeometry.h
MakerMacros.h
HPTDCErrorFlags
Definition: HPTDCErrorFlags.h:15
CTPPSDirectProtonSimulation::roundToPitch_
bool roundToPitch_
Definition: CTPPSDirectProtonSimulation.cc:133
part
part
Definition: HCALResponse.h:20
CTPPSDirectProtonSimulation::stripZeroPosition_
double stripZeroPosition_
internal variable: v position of strip 0, in mm
Definition: CTPPSDirectProtonSimulation.cc:149
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
LHCInfo::crossingAngle
const float crossingAngle() const
Definition: LHCInfo.cc:182
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
RPTopology::last_strip_to_border_dist_
static const double last_strip_to_border_dist_
Definition: RPTopology.h:62
CTPPSDirectProtonSimulation::tokenGeometry_
edm::ESGetToken< CTPPSGeometry, VeryForwardMisalignedGeometryRecord > tokenGeometry_
Definition: CTPPSDirectProtonSimulation.cc:104
LHCInterpolatedOpticalFunctionsSet::Kinematics::y
double y
Definition: LHCInterpolatedOpticalFunctionsSet.h:31
FrontierCondition_GT_autoExpress_cfi.t0
t0
Definition: FrontierCondition_GT_autoExpress_cfi.py:149
CTPPSDirectProtonSimulation::CTPPSDirectProtonSimulation
CTPPSDirectProtonSimulation(const edm::ParameterSet &)
Definition: CTPPSDirectProtonSimulation.cc:151
Service.h
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
CTPPSDirectProtonSimulation::hepMCToken_
edm::EDGetTokenT< edm::HepMCProduct > hepMCToken_
Definition: CTPPSDirectProtonSimulation.cc:110
CTPPSDirectProtonSimulation::produceRecHits_
bool produceRecHits_
Definition: CTPPSDirectProtonSimulation.cc:114
CTPPSDiamondRecHit.h
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
CTPPSDirectProtonSimulation::useTimingEfficiencyPerPlane_
bool useTimingEfficiencyPerPlane_
Definition: CTPPSDirectProtonSimulation.cc:125
CTPPSDirectProtonSimulation::pitchStrips_
double pitchStrips_
strip pitch in mm
Definition: CTPPSDirectProtonSimulation.cc:136
LHCInterpolatedOpticalFunctionsSetCollection.h
CTPPSDetId::sdTimingDiamond
Definition: CTPPSDetId.h:44
CTPPSDirectProtonSimulation::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: CTPPSDirectProtonSimulation.cc:211
LHCInterpolatedOpticalFunctionsSet::Kinematics::x
double x
Definition: LHCInterpolatedOpticalFunctionsSet.h:29
Point3DBase< float, LocalTag >
CTPPSDetId::sdTrackingPixel
Definition: CTPPSDetId.h:44
CTPPSDirectProtonSimulation::useTimingEfficiencyPerRP_
bool useTimingEfficiencyPerRP_
Definition: CTPPSDirectProtonSimulation.cc:123
CTPPSDirectProtonSimulation::tokenDirectSimuData_
edm::ESGetToken< PPSDirectSimulationData, PPSDirectSimulationDataRcd > tokenDirectSimuData_
Definition: CTPPSDirectProtonSimulation.cc:105
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
CTPPSDirectProtonSimulation::pitchPixelsHor_
double pitchPixelsHor_
Definition: CTPPSDirectProtonSimulation.cc:139
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
CTPPSDirectProtonSimulation::verbosity_
unsigned int verbosity_
Definition: CTPPSDirectProtonSimulation.cc:142
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
OpticalFunctionsConfig_cfi.xangle
xangle
Definition: OpticalFunctionsConfig_cfi.py:17
CTPPSDirectProtonSimulation::checkIsHit_
bool checkIsHit_
Definition: CTPPSDirectProtonSimulation.cc:134
CTPPSDirectProtonSimulation::efficiencyMapsPerPlane_
std::map< unsigned int, std::unique_ptr< TH2F > > efficiencyMapsPerPlane_
Definition: CTPPSDirectProtonSimulation.cc:129
CTPPSDetId::arm
uint32_t arm() const
Definition: CTPPSDetId.h:55
CTPPSDirectProtonSimulation::produceHitsRelativeToBeam_
bool produceHitsRelativeToBeam_
Definition: CTPPSDirectProtonSimulation.cc:132
LHCInterpolatedOpticalFunctionsSet::Kinematics::th_y
double th_y
Definition: LHCInterpolatedOpticalFunctionsSet.h:32
CTPPSDiamondDetId.h
CTPPSDiamondDetId
Detector ID class for CTPPS Timing Diamond detectors. Bits [19:31] : Assigend in CTPPSDetId Calss Bit...
Definition: CTPPSDiamondDetId.h:24
edm::ParameterSet
Definition: ParameterSet.h:47
edm::DetSet::emplace_back
decltype(auto) emplace_back(Args &&... args)
Definition: DetSet.h:68
Event.h
RPTopology::pitch_
static const double pitch_
Definition: RPTopology.h:56
CTPPSDirectProtonSimulation::processProton
void processProton(const HepMC::GenVertex *in_vtx, const HepMC::GenParticle *in_trk, const CTPPSGeometry &geometry, const LHCInfo &lhcInfo, const CTPPSBeamParameters &beamParameters, const PPSPixelTopology &ppt, const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions, CLHEP::HepRandomEngine *rndEngine, std::vector< CTPPSLocalTrackLite > &out_tracks, edm::DetSetVector< TotemRPRecHit > &out_strip_hits, edm::DetSetVector< CTPPSPixelRecHit > &out_pixel_hits, edm::DetSetVector< CTPPSDiamondRecHit > &out_diamond_hits, std::map< int, edm::DetSetVector< TotemRPRecHit >> &out_strip_hits_per_particle, std::map< int, edm::DetSetVector< CTPPSPixelRecHit >> &out_pixel_hits_per_particle, std::map< int, edm::DetSetVector< CTPPSDiamondRecHit >> &out_diamond_hits_per_particle) const
Definition: CTPPSDirectProtonSimulation.cc:340
LocalError
Definition: LocalError.h:12
CTPPSDirectProtonSimulation::tokenBeamParameters_
edm::ESGetToken< CTPPSBeamParameters, CTPPSBeamParametersRcd > tokenBeamParameters_
Definition: CTPPSDirectProtonSimulation.cc:101
CTPPSDirectProtonSimulation::produceScoringPlaneHits_
bool produceScoringPlaneHits_
Definition: CTPPSDirectProtonSimulation.cc:113
A
CTPPSDetId
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:31
edm::Service< edm::RandomNumberGenerator >
createfilelist.int
int
Definition: createfilelist.py:10
CTPPSDirectProtonSimulation::tokenLHCInfo_
edm::ESGetToken< LHCInfo, LHCInfoRcd > tokenLHCInfo_
Definition: CTPPSDirectProtonSimulation.cc:100
iEvent
int iEvent
Definition: GenABIO.cc:224
CTPPSInterpolatedOpticsRcd.h
RPTopology::IsHit
static bool IsHit(double u, double v, double insensitiveMargin=0)
Definition: RPTopology.cc:28
edm::stream::EDProducer
Definition: EDProducer.h:38
CTPPSPixelDetId::plane
uint32_t plane() const
Definition: CTPPSPixelDetId.h:37
RPTopology::no_of_strips_
static const unsigned short no_of_strips_
Definition: RPTopology.h:58
edm::EventSetup
Definition: EventSetup.h:58
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
CTPPSBeamParameters::getBeamMom45
double getBeamMom45() const
Definition: CTPPSBeamParameters.cc:55
DetSetVector.h
CTPPSPixelDetId
Definition: CTPPSPixelDetId.h:16
CTPPSDirectProtonSimulation::empiricalAperture56_
std::unique_ptr< TF2 > empiricalAperture56_
Definition: CTPPSDirectProtonSimulation.cc:119
PPSPixelTopology
Definition: PPSPixelTopology.h:22
CTPPSDirectProtonSimulation::pitchPixelsVer_
double pitchPixelsVer_
Definition: CTPPSDirectProtonSimulation.cc:140
edm::ESGetToken< LHCInfo, LHCInfoRcd >
PPSDirectSimulationData.h
CTPPSDirectProtonSimulation::directSimuDataRcdWatcher_
edm::ESWatcher< PPSDirectSimulationDataRcd > directSimuDataRcdWatcher_
Definition: CTPPSDirectProtonSimulation.cc:107
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
alignCSCRings.r
r
Definition: alignCSCRings.py:93
CTPPSDirectProtonSimulation::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: CTPPSDirectProtonSimulation.cc:244
protons_cff.arm
arm
Definition: protons_cff.py:43
profile_2016_postTS2_cff.rpId
rpId
Definition: profile_2016_postTS2_cff.py:21
protons_cff.xi
xi
Definition: protons_cff.py:35
TtFullHadDaughter::B
static const std::string B
Definition: TtFullHadronicEvent.h:9
MaterialEffects_cfi.A
A
Definition: MaterialEffects_cfi.py:11
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
ESInputTag.h
GenParticle.GenParticle
GenParticle
Definition: GenParticle.py:18
eostools.move
def move(src, dest)
Definition: eostools.py:511
extraflags_cff.vtx
vtx
Definition: extraflags_cff.py:18
Frameworkfwd.h
ESWatcher.h
Exception
Definition: hltDiff.cc:245
CTPPSDirectProtonSimulation::pixelTopologyToken_
edm::ESGetToken< PPSPixelTopology, PPSPixelTopologyRcd > pixelTopologyToken_
Definition: CTPPSDirectProtonSimulation.cc:102
EventSetup.h
RPTopology.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
CTPPSDetId.h
Exception.h
edm::DetSetVector::find_or_insert
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:234
genParticles_cff.map
map
Definition: genParticles_cff.py:11
LHCInterpolatedOpticalFunctionsSet::Kinematics::th_x
double th_x
Definition: LHCInterpolatedOpticalFunctionsSet.h:30
ParameterSet.h
LHCInfoRcd.h
HepMCProduct.h
PPSDirectSimulationDataRcd.h
P
std::pair< OmniClusterRef, TrackingParticleRef > P
Definition: BDHadronTrackMonitoringAnalyzer.cc:202
DeDxTools::esConsumes
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::Event
Definition: Event.h:73
CTPPSDirectProtonSimulation::empiricalAperture45_
std::unique_ptr< TF2 > empiricalAperture45_
Definition: CTPPSDirectProtonSimulation.cc:118
PPSPixelTopology::isPixelHit
bool isPixelHit(float xLocalCoordinate, float yLocalCoordinate, bool is3x2) const
Definition: PPSPixelTopology.cc:29
PPSPixelTopologyRcd.h
edm::InputTag
Definition: InputTag.h:15
OpticalFunctionsConfig_cfi.opticalFunctions
opticalFunctions
Definition: OpticalFunctionsConfig_cfi.py:16
CTPPSDirectProtonSimulation::tokenOpticalFunctions_
edm::ESGetToken< LHCInterpolatedOpticalFunctionsSetCollection, CTPPSInterpolatedOpticsRcd > tokenOpticalFunctions_
Definition: CTPPSDirectProtonSimulation.cc:103
CTPPSDirectProtonSimulation::useEmpiricalApertures_
bool useEmpiricalApertures_
Definition: CTPPSDirectProtonSimulation.cc:117
CTPPSDirectProtonSimulation::timeResolutionDiamonds45_
std::unique_ptr< TF1 > timeResolutionDiamonds45_
Definition: CTPPSDirectProtonSimulation.cc:144
CTPPSGeometry::Vector
DetGeomDesc::Translation Vector
Definition: CTPPSGeometry.h:39