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 
17 
20 
26 
29 
32 
35 
40 
41 #include <unordered_map>
42 
43 #include "TMath.h"
44 #include "TMatrixD.h"
45 #include "TVectorD.h"
46 
47 //----------------------------------------------------------------------------------------------------
48 
50 public:
53 
54  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
55 
56 private:
57  void produce(edm::Event &, const edm::EventSetup &) override;
58 
59  void processProton(const HepMC::GenVertex *in_vtx,
60  const HepMC::GenParticle *in_trk,
61  const CTPPSGeometry &geometry,
62  const LHCInfo &lhcInfo,
63  const CTPPSBeamParameters &beamParameters,
65  std::vector<CTPPSLocalTrackLite> &out_tracks,
66 
67  edm::DetSetVector<TotemRPRecHit> &out_strip_hits,
69  edm::DetSetVector<CTPPSDiamondRecHit> &out_diamond_hits,
70 
71  std::map<int, edm::DetSetVector<TotemRPRecHit>> &out_strip_hits_per_particle,
72  std::map<int, edm::DetSetVector<CTPPSPixelRecHit>> &out_pixel_hits_per_particle,
73  std::map<int, edm::DetSetVector<CTPPSDiamondRecHit>> &out_diamond_hits_per_particle) const;
74 
75  // ------------ config file parameters ------------
76 
80 
82 
86 
89 
95 
99 
100  double pitchStrips_;
102 
105 
106  unsigned int verbosity_;
107 
108  // ------------ internal parameters ------------
109 
112 };
113 
114 //----------------------------------------------------------------------------------------------------
115 
117  : lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
118  opticsLabel_(iConfig.getParameter<std::string>("opticsLabel")),
119  hepMCToken_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("hepMCTag"))),
120 
121  produceScoringPlaneHits_(iConfig.getParameter<bool>("produceScoringPlaneHits")),
122  produceRecHits_(iConfig.getParameter<bool>("produceRecHits")),
123 
124  useEmpiricalApertures_(iConfig.getParameter<bool>("useEmpiricalApertures")),
125  empiricalAperture45_xi0_int_(iConfig.getParameter<double>("empiricalAperture45_xi0_int")),
126  empiricalAperture45_xi0_slp_(iConfig.getParameter<double>("empiricalAperture45_xi0_slp")),
127  empiricalAperture45_a_int_(iConfig.getParameter<double>("empiricalAperture45_a_int")),
128  empiricalAperture45_a_slp_(iConfig.getParameter<double>("empiricalAperture45_a_slp")),
129  empiricalAperture56_xi0_int_(iConfig.getParameter<double>("empiricalAperture56_xi0_int")),
130  empiricalAperture56_xi0_slp_(iConfig.getParameter<double>("empiricalAperture56_xi0_slp")),
131  empiricalAperture56_a_int_(iConfig.getParameter<double>("empiricalAperture56_a_int")),
132  empiricalAperture56_a_slp_(iConfig.getParameter<double>("empiricalAperture56_a_slp")),
133 
134  produceHitsRelativeToBeam_(iConfig.getParameter<bool>("produceHitsRelativeToBeam")),
135  roundToPitch_(iConfig.getParameter<bool>("roundToPitch")),
136  checkIsHit_(iConfig.getParameter<bool>("checkIsHit")),
137 
138  pitchStrips_(iConfig.getParameter<double>("pitchStrips")),
139  insensitiveMarginStrips_(iConfig.getParameter<double>("insensitiveMarginStrips")),
140 
141  pitchPixelsHor_(iConfig.getParameter<double>("pitchPixelsHor")),
142  pitchPixelsVer_(iConfig.getParameter<double>("pitchPixelsVer")),
143 
144  verbosity_(iConfig.getUntrackedParameter<unsigned int>("verbosity", 0)) {
146  produces<std::vector<CTPPSLocalTrackLite>>();
147 
148  if (produceRecHits_) {
149  produces<edm::DetSetVector<TotemRPRecHit>>();
150  produces<edm::DetSetVector<CTPPSDiamondRecHit>>();
151  produces<edm::DetSetVector<CTPPSPixelRecHit>>();
152 
153  produces<std::map<int, edm::DetSetVector<TotemRPRecHit>>>();
154  produces<std::map<int, edm::DetSetVector<CTPPSDiamondRecHit>>>();
155  produces<std::map<int, edm::DetSetVector<CTPPSPixelRecHit>>>();
156  }
157 
158  // v position of strip 0
161 }
162 
163 //----------------------------------------------------------------------------------------------------
164 
166  // get input
168  iEvent.getByToken(hepMCToken_, hepmc_prod);
169 
170  // get conditions
171  edm::ESHandle<LHCInfo> hLHCInfo;
172  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
173 
174  edm::ESHandle<CTPPSBeamParameters> hBeamParameters;
175  iSetup.get<CTPPSBeamParametersRcd>().get(hBeamParameters);
176 
178  iSetup.get<CTPPSInterpolatedOpticsRcd>().get(opticsLabel_, hOpticalFunctions);
179 
181  iSetup.get<VeryForwardMisalignedGeometryRecord>().get(geometry);
182 
183  // prepare outputs
184  std::unique_ptr<edm::DetSetVector<TotemRPRecHit>> pStripRecHits(new edm::DetSetVector<TotemRPRecHit>());
185  std::unique_ptr<edm::DetSetVector<CTPPSDiamondRecHit>> pDiamondRecHits(new edm::DetSetVector<CTPPSDiamondRecHit>());
186  std::unique_ptr<edm::DetSetVector<CTPPSPixelRecHit>> pPixelRecHits(new edm::DetSetVector<CTPPSPixelRecHit>());
187 
188  auto pStripRecHitsPerParticle = std::make_unique<std::map<int, edm::DetSetVector<TotemRPRecHit>>>();
189  auto pDiamondRecHitsPerParticle = std::make_unique<std::map<int, edm::DetSetVector<CTPPSDiamondRecHit>>>();
190  auto pPixelRecHitsPerParticle = std::make_unique<std::map<int, edm::DetSetVector<CTPPSPixelRecHit>>>();
191 
192  std::unique_ptr<std::vector<CTPPSLocalTrackLite>> pTracks(new std::vector<CTPPSLocalTrackLite>());
193 
194  // loop over event vertices
195  auto evt = hepmc_prod->GetEvent();
196  for (auto it_vtx = evt->vertices_begin(); it_vtx != evt->vertices_end(); ++it_vtx) {
197  auto vtx = *(it_vtx);
198 
199  // loop over outgoing particles
200  for (auto it_part = vtx->particles_out_const_begin(); it_part != vtx->particles_out_const_end(); ++it_part) {
201  auto part = *(it_part);
202 
203  // accept only stable protons
204  if (part->pdg_id() != 2212)
205  continue;
206 
207  if (part->status() != 1 && part->status() < 83)
208  continue;
209 
211  part,
212  *geometry,
213  *hLHCInfo,
214  *hBeamParameters,
215  *hOpticalFunctions,
216  *pTracks,
217  *pStripRecHits,
218  *pPixelRecHits,
219  *pDiamondRecHits,
220  *pStripRecHitsPerParticle,
221  *pPixelRecHitsPerParticle,
222  *pDiamondRecHitsPerParticle);
223  }
224  }
225 
227  iEvent.put(std::move(pTracks));
228 
229  if (produceRecHits_) {
230  iEvent.put(std::move(pStripRecHits));
231  iEvent.put(std::move(pPixelRecHits));
232  iEvent.put(std::move(pDiamondRecHits));
233 
234  iEvent.put(std::move(pStripRecHitsPerParticle));
235  iEvent.put(std::move(pPixelRecHitsPerParticle));
236  iEvent.put(std::move(pDiamondRecHitsPerParticle));
237  }
238 }
239 
240 //----------------------------------------------------------------------------------------------------
241 
243  const HepMC::GenVertex *in_vtx,
244  const HepMC::GenParticle *in_trk,
245  const CTPPSGeometry &geometry,
246  const LHCInfo &lhcInfo,
247  const CTPPSBeamParameters &beamParameters,
249  std::vector<CTPPSLocalTrackLite> &out_tracks,
250  edm::DetSetVector<TotemRPRecHit> &out_strip_hits,
251  edm::DetSetVector<CTPPSPixelRecHit> &out_pixel_hits,
252  edm::DetSetVector<CTPPSDiamondRecHit> &out_diamond_hits,
253  std::map<int, edm::DetSetVector<TotemRPRecHit>> &out_strip_hits_per_particle,
254  std::map<int, edm::DetSetVector<CTPPSPixelRecHit>> &out_pixel_hits_per_particle,
255  std::map<int, edm::DetSetVector<CTPPSDiamondRecHit>> &out_diamond_hits_per_particle) const {
258 
259  std::stringstream ssLog;
260 
261  // vectors in CMS convention
262  const HepMC::FourVector &vtx_cms = in_vtx->position(); // in mm
263  const HepMC::FourVector &mom_cms = in_trk->momentum();
264 
265  // transformation to LHC/TOTEM convention
266  HepMC::FourVector vtx_lhc(-vtx_cms.x(), vtx_cms.y(), -vtx_cms.z(), vtx_cms.t());
267  HepMC::FourVector mom_lhc(-mom_cms.x(), mom_cms.y(), -mom_cms.z(), mom_cms.t());
268 
269  // determine the LHC arm and related parameters
270  unsigned int arm = 3;
271  double z_sign;
272  double beamMomentum = 0.;
273  double xangle = 0.;
274  double empiricalAperture_xi0_int, empiricalAperture_xi0_slp;
275  double empiricalAperture_a_int, empiricalAperture_a_slp;
276 
277  if (mom_lhc.z() < 0) // sector 45
278  {
279  arm = 0;
280  z_sign = -1;
281  beamMomentum = beamParameters.getBeamMom45();
282  xangle = beamParameters.getHalfXangleX45();
283  empiricalAperture_xi0_int = empiricalAperture45_xi0_int_;
284  empiricalAperture_xi0_slp = empiricalAperture45_xi0_slp_;
285  empiricalAperture_a_int = empiricalAperture45_a_int_;
286  empiricalAperture_a_slp = empiricalAperture45_a_slp_;
287  } else { // sector 56
288  arm = 1;
289  z_sign = +1;
290  beamMomentum = beamParameters.getBeamMom56();
291  xangle = beamParameters.getHalfXangleX56();
292  empiricalAperture_xi0_int = empiricalAperture56_xi0_int_;
293  empiricalAperture_xi0_slp = empiricalAperture56_xi0_slp_;
294  empiricalAperture_a_int = empiricalAperture56_a_int_;
295  empiricalAperture_a_slp = empiricalAperture56_a_slp_;
296  }
297 
298  // calculate kinematics for optics parametrisation
299  const double p = mom_lhc.rho();
300  const double xi = 1. - p / beamMomentum;
301  const double th_x_phys = mom_lhc.x() / p;
302  const double th_y_phys = mom_lhc.y() / p;
303  const double vtx_lhc_eff_x = vtx_lhc.x() - vtx_lhc.z() * (mom_lhc.x() / mom_lhc.z() + xangle);
304  const double vtx_lhc_eff_y = vtx_lhc.y() - vtx_lhc.z() * (mom_lhc.y() / mom_lhc.z());
305 
306  if (verbosity_) {
307  ssLog << "simu: xi = " << xi << ", th_x_phys = " << th_x_phys << ", th_y_phys = " << th_y_phys
308  << ", vtx_lhc_eff_x = " << vtx_lhc_eff_x << ", vtx_lhc_eff_y = " << vtx_lhc_eff_y << std::endl;
309  }
310 
311  // check empirical aperture
313  const auto &xangle = lhcInfo.crossingAngle();
314  const double xi_th = (empiricalAperture_xi0_int + xangle * empiricalAperture_xi0_slp) +
315  (empiricalAperture_a_int + xangle * empiricalAperture_a_slp) * th_x_phys;
316 
317  if (xi > xi_th) {
318  if (verbosity_) {
319  ssLog << "stop because of empirical appertures";
320  edm::LogInfo("CTPPSDirectProtonSimulation") << ssLog.str();
321  }
322 
323  return;
324  }
325  }
326 
327  // transport the proton into each pot/scoring plane
328  for (const auto &ofp : opticalFunctions) {
329  CTPPSDetId rpId(ofp.first);
330  const unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
331 
332  // first check the arm
333  if (rpId.arm() != arm)
334  continue;
335 
336  if (verbosity_)
337  ssLog << " RP " << rpDecId << std::endl;
338 
339  // transport proton
341  vtx_lhc_eff_x * 1E-1, th_x_phys, vtx_lhc_eff_y * 1E-1, th_y_phys, xi}; // conversions: mm -> cm
343  ofp.second.transport(k_in, k_out, true);
344 
345  double b_x = k_out.x * 1E1, b_y = k_out.y * 1E1; // conversions: cm -> mm
346  double a_x = k_out.th_x, a_y = k_out.th_y;
347 
348  // if needed, subtract beam position and angle
350  // determine beam position
351  LHCInterpolatedOpticalFunctionsSet::Kinematics k_be_in = {0., 0., 0., 0., 0.};
353  ofp.second.transport(k_be_in, k_be_out, true);
354 
355  a_x -= k_be_out.th_x;
356  a_y -= k_be_out.th_y;
357  b_x -= k_be_out.x * 1E1;
358  b_y -= k_be_out.y * 1E1; // conversions: cm -> mm
359  }
360 
361  const double z_scoringPlane = ofp.second.getScoringPlaneZ() * 1E1; // conversion: cm --> mm
362 
363  if (verbosity_) {
364  ssLog << " proton transported: a_x = " << a_x << " rad, a_y = " << a_y << " rad, b_x = " << b_x
365  << " mm, b_y = " << b_y << " mm, z = " << z_scoringPlane << " mm" << std::endl;
366  }
367 
368  // save scoring plane hit
370  out_tracks.emplace_back(
371  rpId.rawId(), b_x, 0., b_y, 0., 0., 0., 0., 0., 0., CTPPSpixelLocalTrackReconstructionInfo::invalid, 0, 0., 0.);
372 
373  // stop if rec hits are not to be produced
374  if (!produceRecHits_)
375  continue;
376 
377  // loop over all sensors in the RP
378  for (const auto &detIdInt : geometry.getSensorsInRP(rpId)) {
379  CTPPSDetId detId(detIdInt);
380 
381  // determine the track impact point (in global coordinates)
382  // !! this assumes that local axes (1, 0, 0) and (0, 1, 0) describe the sensor surface
383  CLHEP::Hep3Vector gl_o = geometry.localToGlobal(detId, CLHEP::Hep3Vector(0, 0, 0));
384  CLHEP::Hep3Vector gl_a1 = geometry.localToGlobal(detId, CLHEP::Hep3Vector(1, 0, 0)) - gl_o;
385  CLHEP::Hep3Vector gl_a2 = geometry.localToGlobal(detId, CLHEP::Hep3Vector(0, 1, 0)) - gl_o;
386 
387  TMatrixD A(3, 3);
388  TVectorD B(3);
389  A(0, 0) = a_x;
390  A(0, 1) = -gl_a1.x();
391  A(0, 2) = -gl_a2.x();
392  B(0) = gl_o.x() - b_x;
393  A(1, 0) = a_y;
394  A(1, 1) = -gl_a1.y();
395  A(1, 2) = -gl_a2.y();
396  B(1) = gl_o.y() - b_y;
397  A(2, 0) = z_sign;
398  A(2, 1) = -gl_a1.z();
399  A(2, 2) = -gl_a2.z();
400  B(2) = gl_o.z() - z_scoringPlane;
401  TMatrixD Ai(3, 3);
402  Ai = A.Invert();
403  TVectorD P(3);
404  P = Ai * B;
405 
406  double ze = P(0);
407  CLHEP::Hep3Vector h_glo(a_x * ze + b_x, a_y * ze + b_y, z_sign * ze + z_scoringPlane);
408 
409  // hit in local coordinates
410  CLHEP::Hep3Vector h_loc = geometry.globalToLocal(detId, h_glo);
411 
412  if (verbosity_) {
413  ssLog << std::endl
414  << " de z = " << P(0) << " mm, p1 = " << P(1) << " mm, p2 = " << P(2) << " mm" << std::endl
415  << " h_glo: x = " << h_glo.x() << " mm, y = " << h_glo.y() << " mm, z = " << h_glo.z() << " mm"
416  << std::endl
417  << " h_loc: c1 = " << h_loc.x() << " mm, c2 = " << h_loc.y() << " mm, c3 = " << h_loc.z() << " mm"
418  << std::endl;
419  }
420 
421  // strips
422  if (detId.subdetId() == CTPPSDetId::sdTrackingStrip) {
423  double u = h_loc.x();
424  double v = h_loc.y();
425 
426  if (verbosity_ > 5)
427  ssLog << " u=" << u << ", v=" << v;
428 
429  // is it within detector?
431  if (verbosity_ > 5)
432  ssLog << " | no hit" << std::endl;
433  continue;
434  }
435 
436  // round the measurement
437  if (roundToPitch_) {
438  double m = stripZeroPosition_ - v;
439  signed int strip = (int)floor(m / pitchStrips_ + 0.5);
440 
442 
443  if (verbosity_ > 5)
444  ssLog << " | strip=" << strip;
445  }
446 
447  double sigma = pitchStrips_ / sqrt(12.);
448 
449  if (verbosity_ > 5)
450  ssLog << " | m=" << v << ", sigma=" << sigma << std::endl;
451 
452  edm::DetSet<TotemRPRecHit> &hits = out_strip_hits.find_or_insert(detId);
453  hits.push_back(TotemRPRecHit(v, sigma));
454 
455  edm::DetSet<TotemRPRecHit> &hits_per_particle =
456  out_strip_hits_per_particle[in_trk->barcode()].find_or_insert(detId);
457  hits_per_particle.push_back(TotemRPRecHit(v, sigma));
458  }
459 
460  // diamonds
461  if (detId.subdetId() == CTPPSDetId::sdTimingDiamond) {
462  throw cms::Exception("CTPPSDirectProtonSimulation") << "Diamonds are not yet supported.";
463  }
464 
465  // pixels
466  if (detId.subdetId() == CTPPSDetId::sdTrackingPixel) {
467  if (verbosity_) {
468  CTPPSPixelDetId pixelDetId(detIdInt);
469  ssLog << " pixel plane " << pixelDetId.plane() << ": local hit x = " << h_loc.x()
470  << " mm, y = " << h_loc.y() << " mm, z = " << h_loc.z() << " mm" << std::endl;
471  }
472 
473  bool module3By2 = (geometry.getSensor(detIdInt)->sensorType() != DDD_CTPPS_PIXELS_SENSOR_TYPE_2x2);
474  if (checkIsHit_ && !CTPPSPixelTopology::isPixelHit(h_loc.x(), h_loc.y(), module3By2))
475  continue;
476 
477  if (roundToPitch_) {
478  h_loc.setX(pitchPixelsHor_ * floor(h_loc.x() / pitchPixelsHor_ + 0.5));
479  h_loc.setY(pitchPixelsVer_ * floor(h_loc.y() / pitchPixelsVer_ + 0.5));
480  }
481 
482  if (verbosity_ > 5)
483  ssLog << " hit accepted: m1 = " << h_loc.x() << " mm, m2 = " << h_loc.y() << " mm" << std::endl;
484 
485  const double sigmaHor = pitchPixelsHor_ / sqrt(12.);
486  const double sigmaVer = pitchPixelsVer_ / sqrt(12.);
487 
488  const LocalPoint lp(h_loc.x(), h_loc.y(), h_loc.z());
489  const LocalError le(sigmaHor, 0., sigmaVer);
490 
491  edm::DetSet<CTPPSPixelRecHit> &hits = out_pixel_hits.find_or_insert(detId);
492  hits.push_back(CTPPSPixelRecHit(lp, le));
493 
494  edm::DetSet<CTPPSPixelRecHit> &hits_per_particle =
495  out_pixel_hits_per_particle[in_trk->barcode()].find_or_insert(detId);
496  hits_per_particle.push_back(CTPPSPixelRecHit(lp, le));
497  }
498  }
499  }
500 
501  if (verbosity_)
502  edm::LogInfo("CTPPSDirectProtonSimulation") << ssLog.str();
503 }
504 
505 //----------------------------------------------------------------------------------------------------
506 
509  desc.addUntracked<unsigned int>("verbosity", 0);
510 
511  desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
512  desc.add<std::string>("opticsLabel", "")->setComment("label of the optics records");
513  desc.add<edm::InputTag>("hepMCTag", edm::InputTag("generator", "unsmeared"));
514 
515  desc.add<bool>("produceScoringPlaneHits", true);
516  desc.add<bool>("produceRecHits", true);
517 
518  desc.add<bool>("useEmpiricalApertures", false);
519  desc.add<double>("empiricalAperture45_xi0_int", 0.);
520  desc.add<double>("empiricalAperture45_xi0_slp", 0.);
521  desc.add<double>("empiricalAperture45_a_int", 0.);
522  desc.add<double>("empiricalAperture45_a_slp", 0.);
523  desc.add<double>("empiricalAperture56_xi0_int", 0.);
524  desc.add<double>("empiricalAperture56_xi0_slp", 0.);
525  desc.add<double>("empiricalAperture56_a_int", 0.);
526  desc.add<double>("empiricalAperture56_a_slp", 0.);
527 
528  desc.add<bool>("produceHitsRelativeToBeam", false);
529  desc.add<bool>("roundToPitch", true);
530  desc.add<bool>("checkIsHit", true);
531  desc.add<double>("pitchStrips", 66.e-3); // in mm
532  desc.add<double>("insensitiveMarginStrips", 34.e-3); // in mm
533 
534  desc.add<double>("pitchPixelsHor", 100.e-3)->setComment("x in local coordinates, in mm");
535  desc.add<double>("pitchPixelsVer", 150.e-3)->setComment("y in local coordinates, in mm");
536 
537  descriptions.add("ctppsDirectProtonSimulation", desc);
538 }
539 
540 //----------------------------------------------------------------------------------------------------
541 
uint32_t station() const
Definition: CTPPSDetId.h:58
uint32_t plane() const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
void push_back(const T &t)
Definition: DetSet.h:68
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const std::string & sensorType() const
Definition: DetGeomDesc.h:70
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool produceScoringPlaneHits_
flags what output to be produced
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
static bool isPixelHit(float xLocalCoordinate, float yLocalCoordinate, bool is3x2=true)
static bool IsHit(double u, double v, double insensitiveMargin=0)
Definition: RPTopology.cc:36
static const unsigned short no_of_strips_
Definition: RPTopology.h:60
void produce(edm::Event &, const edm::EventSetup &) override
const std::string DDD_CTPPS_PIXELS_SENSOR_TYPE_2x2
Definition: CTPPSDDDNames.h:16
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:254
edm::EDGetTokenT< edm::HepMCProduct > hepMCToken_
float const crossingAngle() const
Definition: LHCInfo.cc:176
static const double y_width_
Definition: RPTopology.h:62
const DetGeomDesc * getSensor(unsigned int id) const
returns geometry of a detector performs necessary checks, returns NULL if fails
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
void processProton(const HepMC::GenVertex *in_vtx, const HepMC::GenParticle *in_trk, const CTPPSGeometry &geometry, const LHCInfo &lhcInfo, const CTPPSBeamParameters &beamParameters, const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions, 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
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
Reconstructed hit in TOTEM RP.
Definition: TotemRPRecHit.h:18
double getBeamMom56() const
double insensitiveMarginStrips_
size of insensitive margin at sensor&#39;s edge facing the beam, in mm
double getHalfXangleX56() const
static const std::string B
uint32_t arm() const
Definition: CTPPSDetId.h:51
CLHEP::Hep3Vector localToGlobal(const DetGeomDesc *, const CLHEP::Hep3Vector &) const
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool checkApertures_
simulation parameters
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
The manager class for TOTEM RP geometry.
Definition: CTPPSGeometry.h:33
static const double last_strip_to_border_dist_
Definition: RPTopology.h:64
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
part
Definition: HCALResponse.h:20
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::pair< OmniClusterRef, TrackingParticleRef > P
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
double stripZeroPosition_
internal variable: v position of strip 0, in mm
ESHandle< TrackerGeometry > geometry
HLT enums.
Event setup record containing the misaligned geometry information. It is used for alignment studies o...
static const double pitch_
Definition: RPTopology.h:58
T get() const
Definition: EventSetup.h:71
double getBeamMom45() const
const std::set< unsigned int > & getSensorsInRP(unsigned int) const
after checks returns set of sensor ids corresponding to the given RP id
double getHalfXangleX45() const
CLHEP::Hep3Vector globalToLocal(const DetGeomDesc *, const CLHEP::Hep3Vector &) const
def move(src, dest)
Definition: eostools.py:511
CTPPSDirectProtonSimulation(const edm::ParameterSet &)
uint32_t rp() const
Definition: CTPPSDetId.h:65