CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
CTPPSDirectProtonSimulation Class Reference
Inheritance diagram for CTPPSDirectProtonSimulation:
edm::stream::EDProducer<>

Public Member Functions

 CTPPSDirectProtonSimulation (const edm::ParameterSet &)
 
 ~CTPPSDirectProtonSimulation () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void processProton (const HepMC::GenVertex *in_vtx, const HepMC::GenParticle *in_trk, const CTPPSGeometry &geometry, 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) const
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Static Private Member Functions

static bool isPixelHit (float xLocalCoordinate, float yLocalCoordinate, bool is3x2=true)
 

Private Attributes

bool checkApertures_
 simulation parameters More...
 
bool checkIsHit_
 
double empiricalAperture45_a_
 
double empiricalAperture45_xi0_
 
double empiricalAperture56_a_
 
double empiricalAperture56_xi0_
 
edm::EDGetTokenT< edm::HepMCProducthepMCToken_
 input tag More...
 
double insensitiveMarginStrips_
 size of insensitive margin at sensor's edge facing the beam, in mm More...
 
double pitchPixelsHor_
 
double pitchPixelsVer_
 
double pitchStrips_
 strip pitch in mm More...
 
bool produceHitsRelativeToBeam_
 
bool produceRecHits_
 
bool produceScoringPlaneHits_
 flags what output to be produced More...
 
bool roundToPitch_
 
double stripZeroPosition_
 internal variable: v position of strip 0, in mm More...
 
bool useEmpiricalApertures_
 
unsigned int verbosity_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Definition at line 46 of file CTPPSDirectProtonSimulation.cc.

Constructor & Destructor Documentation

CTPPSDirectProtonSimulation::CTPPSDirectProtonSimulation ( const edm::ParameterSet iConfig)
explicit

Definition at line 117 of file CTPPSDirectProtonSimulation.cc.

References RPTopology::last_strip_to_border_dist_, RPTopology::no_of_strips_, RPTopology::pitch_, produceRecHits_, produceScoringPlaneHits_, stripZeroPosition_, and RPTopology::y_width_.

117  :
118  hepMCToken_( consumes<edm::HepMCProduct>( iConfig.getParameter<edm::InputTag>( "hepMCTag" ) ) ),
119 
120  produceScoringPlaneHits_( iConfig.getParameter<bool>( "produceScoringPlaneHits" ) ),
121  produceRecHits_( iConfig.getParameter<bool>( "produceRecHits" ) ),
122 
123  useEmpiricalApertures_( iConfig.getParameter<bool>( "useEmpiricalApertures" ) ),
124  empiricalAperture45_xi0_( iConfig.getParameter<double>( "empiricalAperture45_xi0" ) ),
125  empiricalAperture45_a_( iConfig.getParameter<double>( "empiricalAperture45_a" ) ),
126  empiricalAperture56_xi0_( iConfig.getParameter<double>( "empiricalAperture56_xi0" ) ),
127  empiricalAperture56_a_( iConfig.getParameter<double>( "empiricalAperture56_a" ) ),
128 
129  produceHitsRelativeToBeam_( iConfig.getParameter<bool>( "produceHitsRelativeToBeam" ) ),
130  roundToPitch_( iConfig.getParameter<bool>( "roundToPitch" ) ),
131  checkIsHit_( iConfig.getParameter<bool>( "checkIsHit" ) ),
132 
133  pitchStrips_( iConfig.getParameter<double>( "pitchStrips" ) ),
134  insensitiveMarginStrips_( iConfig.getParameter<double>( "insensitiveMarginStrips" ) ),
135 
136  pitchPixelsHor_( iConfig.getParameter<double>( "pitchPixelsHor" ) ),
137  pitchPixelsVer_( iConfig.getParameter<double>( "pitchPixelsVer" ) ),
138 
139  verbosity_( iConfig.getUntrackedParameter<unsigned int>( "verbosity", 0 ) )
140 {
142  produces<std::vector<CTPPSLocalTrackLite>>();
143 
144  if (produceRecHits_)
145  {
146  produces<edm::DetSetVector<TotemRPRecHit>>();
147  produces<edm::DetSetVector<CTPPSDiamondRecHit>>();
148  produces<edm::DetSetVector<CTPPSPixelRecHit>>();
149  }
150 
151  // v position of strip 0
153 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool produceScoringPlaneHits_
flags what output to be produced
static const unsigned short no_of_strips_
Definition: RPTopology.h:60
edm::EDGetTokenT< edm::HepMCProduct > hepMCToken_
input tag
static const double y_width_
Definition: RPTopology.h:62
double insensitiveMarginStrips_
size of insensitive margin at sensor&#39;s edge facing the beam, in mm
static const double last_strip_to_border_dist_
Definition: RPTopology.h:64
double stripZeroPosition_
internal variable: v position of strip 0, in mm
static const double pitch_
Definition: RPTopology.h:58
CTPPSDirectProtonSimulation::~CTPPSDirectProtonSimulation ( )
inlineoverride

Member Function Documentation

void CTPPSDirectProtonSimulation::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 456 of file CTPPSDirectProtonSimulation.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addUntracked(), and DEFINE_FWK_MODULE.

Referenced by ~CTPPSDirectProtonSimulation().

457 {
459  desc.addUntracked<unsigned int>("verbosity", 0);
460  desc.add<edm::InputTag>("hepMCTag", edm::InputTag("generator", "unsmeared"));
461  desc.add<bool>("produceScoringPlaneHits", true);
462  desc.add<bool>("produceRecHits", true);
463  desc.add<bool>("useEmpiricalApertures", false);
464  desc.add<double>("empiricalAperture45_xi0", 0.);
465  desc.add<double>("empiricalAperture45_a", 0.);
466  desc.add<double>("empiricalAperture56_xi0", 0.);
467  desc.add<double>("empiricalAperture56_a", 0.);
468  desc.add<bool>("produceHitsRelativeToBeam", false);
469  desc.add<bool>("roundToPitch", true);
470  desc.add<bool>("checkIsHit", true);
471  desc.add<double>("pitchStrips", 66.e-3); // in mm
472  desc.add<double>("insensitiveMarginStrips", 34.e-3); // in mm
473  desc.add<double>("pitchPixelsHor", 100.e-3)->setComment("x in local coordinates, in mm");
474  desc.add<double>("pitchPixelsVer", 150.e-3)->setComment("y in local coordinates, in mm");
475 
476  descriptions.add("ctppsDirectProtonSimulation", desc);
477 }
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static bool CTPPSDirectProtonSimulation::isPixelHit ( float  xLocalCoordinate,
float  yLocalCoordinate,
bool  is3x2 = true 
)
inlinestaticprivate

Definition at line 65 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

66  {
67  float tmpXlocalCoordinate = xLocalCoordinate + (79*0.1 + 0.2);
68  float tmpYlocalCoordinate = yLocalCoordinate + (0.15*51 + 0.3*2 + 0.15*25);
69 
70  if(tmpXlocalCoordinate<0) return false;
71  if(tmpYlocalCoordinate<0) return false;
72  int xModuleSize = 0.1*79 + 0.2*2 + 0.1*79; // mm - 100 um pitch direction
73  int yModuleSize; // mm - 150 um pitch direction
74  if (is3x2) yModuleSize = 0.15*51 + 0.3*2 + 0.15*50 + 0.3*2 + 0.15*51;
75  else yModuleSize = 0.15*51 + 0.3*2 + 0.15*51;
76  if(tmpXlocalCoordinate>xModuleSize) return false;
77  if(tmpYlocalCoordinate>yModuleSize) return false;
78  return true;
79  }
void CTPPSDirectProtonSimulation::processProton ( const HepMC::GenVertex *  in_vtx,
const HepMC::GenParticle *  in_trk,
const CTPPSGeometry geometry,
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 
) const
private

xi is positive for diffractive protons, thus proton momentum p = (1-xi) * p_nom horizontal component of proton momentum: p_x = th_x * (1-xi) * p_nom

Definition at line 215 of file CTPPSDirectProtonSimulation.cc.

References patCaloMETCorrections_cff::A, CTPPSDetId::arm(), TtFullHadDaughter::B, checkIsHit_, empiricalAperture45_a_, empiricalAperture45_xi0_, empiricalAperture56_a_, empiricalAperture56_xi0_, Exception, edm::DetSetVector< T >::find_or_insert(), CTPPSBeamParameters::getBeamMom45(), CTPPSBeamParameters::getBeamMom56(), CTPPSBeamParameters::getHalfXangleX45(), CTPPSBeamParameters::getHalfXangleX56(), CTPPSGeometry::getSensorsInRP(), CTPPSGeometry::globalToLocal(), hfClusterShapes_cfi::hits, insensitiveMarginStrips_, createfilelist::int, invalid, RPTopology::IsHit(), isPixelHit(), CTPPSGeometry::localToGlobal(), funct::m, AlCaHLTBitMon_ParallelJobs::p, pitchPixelsHor_, pitchPixelsVer_, pitchStrips_, CTPPSPixelDetId::plane(), produceHitsRelativeToBeam_, produceRecHits_, produceScoringPlaneHits_, edm::DetSet< T >::push_back(), DetId::rawId(), roundToPitch_, CTPPSDetId::rp(), CTPPSDetId::sdTimingDiamond, CTPPSDetId::sdTrackingPixel, CTPPSDetId::sdTrackingStrip, mathSSE::sqrt(), CTPPSDetId::station(), digitizers_cfi::strip, stripZeroPosition_, DetId::subdetId(), LHCInterpolatedOpticalFunctionsSet::Kinematics::th_x, LHCInterpolatedOpticalFunctionsSet::Kinematics::th_y, useEmpiricalApertures_, findQualityFiles::v, verbosity_, LHCInterpolatedOpticalFunctionsSet::Kinematics::x, hybridSuperClusters_cfi::xi, and LHCInterpolatedOpticalFunctionsSet::Kinematics::y.

Referenced by produce(), and ~CTPPSDirectProtonSimulation().

220 {
223 
224  std::stringstream ssLog;
225 
226  // vectors in CMS convention
227  const HepMC::FourVector& vtx_cms = in_vtx->position(); // in mm
228  const HepMC::FourVector& mom_cms = in_trk->momentum();
229 
230  // transformation to LHC/TOTEM convention
231  HepMC::FourVector vtx_lhc(-vtx_cms.x(), vtx_cms.y(), -vtx_cms.z(), vtx_cms.t());
232  HepMC::FourVector mom_lhc(-mom_cms.x(), mom_cms.y(), -mom_cms.z(), mom_cms.t());
233 
234  // determine the LHC arm and related parameters
235  unsigned int arm = 3;
236  double z_sign;
237  double beamMomentum = 0.;
238  double xangle = 0.;
239  double empiricalAperture_xi0, empiricalAperture_a;
240 
241  if (mom_lhc.z() < 0) // sector 45
242  {
243  arm = 0;
244  z_sign = -1;
245  beamMomentum = beamParameters.getBeamMom45();
246  xangle = beamParameters.getHalfXangleX45();
247  empiricalAperture_xi0 = empiricalAperture45_xi0_;
248  empiricalAperture_a = empiricalAperture45_a_;
249  } else { // sector 56
250  arm = 1;
251  z_sign = +1;
252  beamMomentum = beamParameters.getBeamMom56();
253  xangle = beamParameters.getHalfXangleX56();
254  empiricalAperture_xi0 = empiricalAperture56_xi0_;
255  empiricalAperture_a = empiricalAperture56_a_;
256  }
257 
258  // calculate kinematics for optics parametrisation
259  const double p = mom_lhc.rho();
260  const double xi = 1. - p / beamMomentum;
261  const double th_x_phys = mom_lhc.x() / p;
262  const double th_y_phys = mom_lhc.y() / p;
263  const double vtx_lhc_eff_x = vtx_lhc.x() - vtx_lhc.z() * (mom_lhc.x() / mom_lhc.z() + xangle);
264  const double vtx_lhc_eff_y = vtx_lhc.y() - vtx_lhc.z() * (mom_lhc.y() / mom_lhc.z());
265 
266  if (verbosity_)
267  {
268  ssLog << "simu: xi = " << xi << ", th_x_phys = " << th_x_phys << ", th_y_phys = " << th_y_phys
269  << ", vtx_lhc_eff_x = " << vtx_lhc_eff_x << ", vtx_lhc_eff_y = " << vtx_lhc_eff_y << std::endl;
270  }
271 
272  // check empirical aperture
274  {
275  const double xi_th = empiricalAperture_xi0 + th_x_phys * empiricalAperture_a;
276  if (xi > xi_th)
277  {
278  if (verbosity_)
279  {
280  ssLog << "stop because of empirical appertures";
281  edm::LogInfo("CTPPSDirectProtonSimulation") << ssLog.str();
282  }
283 
284  return;
285  }
286  }
287 
288  // transport the proton into each pot/scoring plane
289  for (const auto &ofp : opticalFunctions)
290  {
291  CTPPSDetId rpId(ofp.first);
292  const unsigned int rpDecId = rpId.arm()*100 + rpId.station()*10 + rpId.rp();
293 
294  // first check the arm
295  if (rpId.arm() != arm)
296  continue;
297 
298  if (verbosity_)
299  ssLog << " RP " << rpDecId << std::endl;
300 
301  // transport proton
302  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in = { vtx_lhc_eff_x*1E-1, th_x_phys, vtx_lhc_eff_y*1E-1, th_y_phys, xi }; // conversions: mm -> cm
304  ofp.second.transport(k_in, k_out, true);
305 
306  double b_x = k_out.x * 1E1, b_y = k_out.y * 1E1; // conversions: cm -> mm
307  double a_x = k_out.th_x, a_y = k_out.th_y;
308 
309  // if needed, subtract beam position and angle
311  {
312  // determine beam position
313  LHCInterpolatedOpticalFunctionsSet::Kinematics k_be_in = { 0., 0., 0., 0., 0. };
315  ofp.second.transport(k_be_in, k_be_out, true);
316 
317  a_x -= k_be_out.th_x; a_y -= k_be_out.th_y;
318  b_x -= k_be_out.x * 1E1; b_y -= k_be_out.y * 1E1; // conversions: cm -> mm
319  }
320 
321  const double z_scoringPlane = ofp.second.getScoringPlaneZ() * 1E1; // conversion: cm --> mm
322 
323  if (verbosity_)
324  {
325  ssLog << " proton transported: a_x = " << a_x << " rad, a_y = " << a_y << " rad, b_x = " << b_x <<
326  " mm, b_y = " << b_y << " mm, z = " << z_scoringPlane << " mm" << std::endl;
327  }
328 
329  // save scoring plane hit
331  out_tracks.emplace_back(rpId.rawId(), b_x, 0., b_y, 0., 0., 0., 0., 0., 0., CTPPSpixelLocalTrackReconstructionInfo::invalid, 0, 0., 0.);
332 
333  // stop if rec hits are not to be produced
334  if (!produceRecHits_)
335  continue;
336 
337  // loop over all sensors in the RP
338  for (const auto& detIdInt : geometry.getSensorsInRP(rpId))
339  {
340  CTPPSDetId detId(detIdInt);
341 
342  // determine the track impact point (in global coordinates)
343  // !! this assumes that local axes (1, 0, 0) and (0, 1, 0) describe the sensor surface
344  CLHEP::Hep3Vector gl_o = geometry.localToGlobal(detId, CLHEP::Hep3Vector(0, 0, 0));
345  CLHEP::Hep3Vector gl_a1 = geometry.localToGlobal(detId, CLHEP::Hep3Vector(1, 0, 0)) - gl_o;
346  CLHEP::Hep3Vector gl_a2 = geometry.localToGlobal(detId, CLHEP::Hep3Vector(0, 1, 0)) - gl_o;
347 
348  TMatrixD A(3, 3);
349  TVectorD B(3);
350  A(0, 0) = a_x; A(0, 1) = -gl_a1.x(); A(0, 2) = -gl_a2.x(); B(0) = gl_o.x() - b_x;
351  A(1, 0) = a_y; A(1, 1) = -gl_a1.y(); A(1, 2) = -gl_a2.y(); B(1) = gl_o.y() - b_y;
352  A(2, 0) = z_sign; A(2, 1) = -gl_a1.z(); A(2, 2) = -gl_a2.z(); B(2) = gl_o.z() - z_scoringPlane;
353  TMatrixD Ai(3, 3);
354  Ai = A.Invert();
355  TVectorD P(3);
356  P = Ai * B;
357 
358  double ze = P(0);
359  CLHEP::Hep3Vector h_glo(a_x * ze + b_x, a_y * ze + b_y, z_sign*ze + z_scoringPlane);
360 
361  // hit in local coordinates
362  CLHEP::Hep3Vector h_loc = geometry.globalToLocal(detId, h_glo);
363 
364  if (verbosity_)
365  {
366  ssLog << std::endl
367  << " de z = " << P(0) << " mm, p1 = " << P(1) << " mm, p2 = " << P(2) << " mm" << std::endl
368  << " h_glo: x = " << h_glo.x() << " mm, y = " << h_glo.y() << " mm, z = " << h_glo.z() << " mm" << std::endl
369  << " h_loc: c1 = " << h_loc.x() << " mm, c2 = " << h_loc.y() << " mm, c3 = " << h_loc.z() << " mm" << std::endl;
370  }
371 
372  // strips
373  if (detId.subdetId() == CTPPSDetId::sdTrackingStrip)
374  {
375  double u = h_loc.x();
376  double v = h_loc.y();
377 
378  if (verbosity_ > 5)
379  ssLog << " u=" << u << ", v=" << v;
380 
381  // is it within detector?
383  {
384  if (verbosity_ > 5)
385  ssLog << " | no hit" << std::endl;
386  continue;
387  }
388 
389  // round the measurement
390  if (roundToPitch_)
391  {
392  double m = stripZeroPosition_ - v;
393  signed int strip = (int) floor(m / pitchStrips_ + 0.5);
394 
396 
397  if (verbosity_ > 5)
398  ssLog << " | strip=" << strip;
399  }
400 
401  double sigma = pitchStrips_ / sqrt(12.);
402 
403  if (verbosity_ > 5)
404  ssLog << " | m=" << v << ", sigma=" << sigma << std::endl;
405 
406  edm::DetSet<TotemRPRecHit> &hits = out_strip_hits.find_or_insert(detId);
407  hits.push_back(TotemRPRecHit(v, sigma));
408  }
409 
410  // diamonds
411  if (detId.subdetId() == CTPPSDetId::sdTimingDiamond)
412  {
413  throw cms::Exception("CTPPSDirectProtonSimulation") << "Diamonds are not yet supported.";
414  }
415 
416  // pixels
417  if (detId.subdetId() == CTPPSDetId::sdTrackingPixel)
418  {
419  if (verbosity_)
420  {
421  CTPPSPixelDetId pixelDetId(detIdInt);
422  ssLog << " pixel plane " << pixelDetId.plane() << ": local hit x = " << h_loc.x()
423  << " mm, y = " << h_loc.y() << " mm, z = " << h_loc.z() << " mm" << std::endl;
424  }
425 
426  if (checkIsHit_ && !isPixelHit(h_loc.x(), h_loc.y()))
427  continue;
428 
429  if (roundToPitch_)
430  {
431  h_loc.setX( pitchPixelsHor_ * floor(h_loc.x()/pitchPixelsHor_ + 0.5) );
432  h_loc.setY( pitchPixelsVer_ * floor(h_loc.y()/pitchPixelsVer_ + 0.5) );
433  }
434 
435  if (verbosity_ > 5)
436  ssLog << " hit accepted: m1 = " << h_loc.x() << " mm, m2 = " << h_loc.y() << " mm" << std::endl;
437 
438  const double sigmaHor = pitchPixelsHor_ / sqrt(12.);
439  const double sigmaVer = pitchPixelsVer_ / sqrt(12.);
440 
441  const LocalPoint lp(h_loc.x(), h_loc.y(), h_loc.z());
442  const LocalError le(sigmaHor, 0., sigmaVer);
443 
444  edm::DetSet<CTPPSPixelRecHit> &hits = out_pixel_hits.find_or_insert(detId);
445  hits.push_back(CTPPSPixelRecHit(lp, le));
446  }
447  }
448  }
449 
450  if (verbosity_)
451  edm::LogInfo("CTPPSDirectProtonSimulation") << ssLog.str();
452 }
void push_back(const T &t)
Definition: DetSet.h:68
bool produceScoringPlaneHits_
flags what output to be produced
static bool IsHit(double u, double v, double insensitiveMargin=0)
Definition: RPTopology.cc:36
static bool isPixelHit(float xLocalCoordinate, float yLocalCoordinate, bool is3x2=true)
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:254
T sqrt(T t)
Definition: SSEVec.h:18
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
CLHEP::Hep3Vector localToGlobal(const DetGeomDesc *, const CLHEP::Hep3Vector &) const
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
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
void CTPPSDirectProtonSimulation::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 157 of file CTPPSDirectProtonSimulation.cc.

References geometry, edm::EventSetup::get(), edm::Event::getByToken(), edm::HepMCProduct::GetEvent(), hepMCToken_, eostools::move(), processProton(), produceRecHits_, produceScoringPlaneHits_, edm::Event::put(), and badGlobalMuonTaggersAOD_cff::vtx.

Referenced by ~CTPPSDirectProtonSimulation().

158 {
159  // get input
161  iEvent.getByToken(hepMCToken_, hepmc_prod);
162 
163  // get conditions
164  edm::ESHandle<CTPPSBeamParameters> hBeamParameters;
165  iSetup.get<CTPPSBeamParametersRcd>().get(hBeamParameters);
166 
168  iSetup.get<CTPPSInterpolatedOpticsRcd>().get(hOpticalFunctions);
169 
171  iSetup.get<VeryForwardMisalignedGeometryRecord>().get(geometry);
172 
173  // prepare outputs
174  std::unique_ptr<edm::DetSetVector<TotemRPRecHit>> pStripRecHits(new edm::DetSetVector<TotemRPRecHit>());
175  std::unique_ptr<edm::DetSetVector<CTPPSDiamondRecHit>> pDiamondRecHits(new edm::DetSetVector<CTPPSDiamondRecHit>());
176  std::unique_ptr<edm::DetSetVector<CTPPSPixelRecHit>> pPixelRecHits(new edm::DetSetVector<CTPPSPixelRecHit>());
177 
178  std::unique_ptr<std::vector<CTPPSLocalTrackLite>> pTracks(new std::vector<CTPPSLocalTrackLite>());
179 
180  // loop over event vertices
181  auto evt = new HepMC::GenEvent( *hepmc_prod->GetEvent() );
182  for ( auto it_vtx = evt->vertices_begin(); it_vtx != evt->vertices_end(); ++it_vtx )
183  {
184  auto vtx = *( it_vtx );
185 
186  // loop over outgoing particles
187  for ( auto it_part = vtx->particles_out_const_begin(); it_part != vtx->particles_out_const_end(); ++it_part )
188  {
189  auto part = *( it_part );
190 
191  // accept only stable protons
192  if ( part->pdg_id() != 2212 )
193  continue;
194 
195  if ( part->status() != 1 && part->status() < 83 )
196  continue;
197 
198  processProton(vtx, part, *geometry, *hBeamParameters, *hOpticalFunctions, *pTracks, *pStripRecHits, *pPixelRecHits, *pDiamondRecHits);
199  }
200  }
201 
203  iEvent.put(std::move(pTracks));
204 
205  if (produceRecHits_)
206  {
207  iEvent.put(std::move(pStripRecHits));
208  iEvent.put(std::move(pPixelRecHits));
209  iEvent.put(std::move(pDiamondRecHits));
210  }
211 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool produceScoringPlaneHits_
flags what output to be produced
edm::EDGetTokenT< edm::HepMCProduct > hepMCToken_
input tag
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
part
Definition: HCALResponse.h:20
ESHandle< TrackerGeometry > geometry
Event setup record containing the misaligned geometry information. It is used for alignment studies o...
T get() const
Definition: EventSetup.h:71
void processProton(const HepMC::GenVertex *in_vtx, const HepMC::GenParticle *in_trk, const CTPPSGeometry &geometry, 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) const
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

bool CTPPSDirectProtonSimulation::checkApertures_
private

simulation parameters

Definition at line 91 of file CTPPSDirectProtonSimulation.cc.

bool CTPPSDirectProtonSimulation::checkIsHit_
private

Definition at line 99 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

double CTPPSDirectProtonSimulation::empiricalAperture45_a_
private

Definition at line 94 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

double CTPPSDirectProtonSimulation::empiricalAperture45_xi0_
private

Definition at line 94 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

double CTPPSDirectProtonSimulation::empiricalAperture56_a_
private

Definition at line 95 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

double CTPPSDirectProtonSimulation::empiricalAperture56_xi0_
private

Definition at line 95 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

edm::EDGetTokenT<edm::HepMCProduct> CTPPSDirectProtonSimulation::hepMCToken_
private

input tag

Definition at line 84 of file CTPPSDirectProtonSimulation.cc.

Referenced by produce().

double CTPPSDirectProtonSimulation::insensitiveMarginStrips_
private

size of insensitive margin at sensor's edge facing the beam, in mm

Definition at line 102 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

double CTPPSDirectProtonSimulation::pitchPixelsHor_
private

Definition at line 104 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

double CTPPSDirectProtonSimulation::pitchPixelsVer_
private

Definition at line 105 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

double CTPPSDirectProtonSimulation::pitchStrips_
private

strip pitch in mm

Definition at line 101 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

bool CTPPSDirectProtonSimulation::produceHitsRelativeToBeam_
private

Definition at line 97 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

bool CTPPSDirectProtonSimulation::produceRecHits_
private
bool CTPPSDirectProtonSimulation::produceScoringPlaneHits_
private

flags what output to be produced

Definition at line 87 of file CTPPSDirectProtonSimulation.cc.

Referenced by CTPPSDirectProtonSimulation(), processProton(), and produce().

bool CTPPSDirectProtonSimulation::roundToPitch_
private

Definition at line 98 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

double CTPPSDirectProtonSimulation::stripZeroPosition_
private

internal variable: v position of strip 0, in mm

Definition at line 112 of file CTPPSDirectProtonSimulation.cc.

Referenced by CTPPSDirectProtonSimulation(), and processProton().

bool CTPPSDirectProtonSimulation::useEmpiricalApertures_
private

Definition at line 93 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().

unsigned int CTPPSDirectProtonSimulation::verbosity_
private

Definition at line 107 of file CTPPSDirectProtonSimulation.cc.

Referenced by processProton().