CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Protected Member Functions | Protected Attributes
PPSFastLocalSimulation Class Reference

Fast (no G4) proton simulation in within one station. Uses misaligned geometry. More...

Inheritance diagram for PPSFastLocalSimulation:
edm::stream::EDProducer<>

Classes

struct  Distribution
 

Public Member Functions

 PPSFastLocalSimulation (const edm::ParameterSet &)
 
 ~PPSFastLocalSimulation () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Protected Member Functions

void GenerateTrack (unsigned int pi, CLHEP::HepRandomEngine &rndEng, HepMC::GenEvent *gEv, std::unique_ptr< edm::DetSetVector< TotemRPRecHit >> &stripHitColl, std::unique_ptr< edm::DetSetVector< CTPPSDiamondRecHit >> &diamondHitColl, std::unique_ptr< edm::DetSetVector< CTPPSPixelRecHit >> &pixelHitColl, const CTPPSGeometry &geometry)
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Protected Attributes

Distribution angular_dist_
 angular parameters in rad More...
 
edm::ESGetToken< CTPPSGeometry, VeryForwardMisalignedGeometryRecordesTokenGeometry_
 
double insensitiveMarginStrips_
 size of insensitive margin at sensor's edge facing the beam, in mm More...
 
bool makeHepMC_
 whether a HepMC description of the proton shall be saved in the event More...
 
bool makeHits_
 whether the hits of the proton shall be calculated and saved More...
 
double particle_E_
 particle energy and momentum More...
 
double particle_p_
 
unsigned int particlesPerEvent_
 number of particles to generate per event More...
 
double pitchDiamonds_
 
double pitchPixels_
 
double pitchStrips_
 in mm More...
 
Distribution position_dist_
 position parameters in mm More...
 
bool roundToPitch_
 whether measurement values shall be rounded to the nearest strip More...
 
std::vector< unsigned int > RPs_
 the list of RPs to simulate More...
 
double stripZeroPosition_
 v position of strip 0, in mm More...
 
unsigned int verbosity_
 verbosity level More...
 
double z0_
 the "origin" of tracks, in mm More...
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Fast (no G4) proton simulation in within one station. Uses misaligned geometry.

Definition at line 44 of file PPSFastLocalSimulation.cc.

Constructor & Destructor Documentation

◆ PPSFastLocalSimulation()

PPSFastLocalSimulation::PPSFastLocalSimulation ( const edm::ParameterSet ps)

Definition at line 186 of file PPSFastLocalSimulation.cc.

References RPTopology::last_strip_to_border_dist_, makeHepMC_, makeHits_, RPTopology::no_of_strips_, RPTopology::pitch_, stripZeroPosition_, and RPTopology::y_width_.

187  : verbosity_(ps.getUntrackedParameter<unsigned int>("verbosity", 0)),
188 
189  makeHepMC_(ps.getParameter<bool>("makeHepMC")),
190  makeHits_(ps.getParameter<bool>("makeHits")),
191 
192  RPs_(ps.getParameter<vector<unsigned int>>("RPs")),
193 
194  particlesPerEvent_(ps.getParameter<unsigned int>("particlesPerEvent")),
195  particle_E_(ps.getParameter<double>("particle_E")),
196  particle_p_(ps.getParameter<double>("particle_p")),
197  z0_(ps.getParameter<double>("z0")),
198 
199  roundToPitch_(ps.getParameter<bool>("roundToPitch")),
200  pitchStrips_(ps.getParameter<double>("pitchStrips")),
201  pitchDiamonds_(ps.getParameter<double>("pitchDiamonds")),
202  pitchPixels_(ps.getParameter<double>("pitchPixels")),
203 
204  insensitiveMarginStrips_(ps.getParameter<double>("insensitiveMarginStrips")),
205 
206  position_dist_(ps.getParameterSet("position_distribution")),
207  angular_dist_(ps.getParameterSet("angular_distribution")),
208 
210  // v position of strip 0
213 
214  // register the output
215  if (makeHepMC_)
216  produces<HepMCProduct>();
217 
218  if (makeHits_) {
219  produces<DetSetVector<TotemRPRecHit>>();
220  produces<DetSetVector<CTPPSDiamondRecHit>>();
221  produces<DetSetVector<CTPPSPixelRecHit>>();
222  }
223 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
unsigned int verbosity_
verbosity level
Distribution position_dist_
position parameters in mm
ParameterSet const & getParameterSet(std::string const &) const
bool makeHits_
whether the hits of the proton shall be calculated and saved
static const unsigned short no_of_strips_
Definition: RPTopology.h:53
double z0_
the "origin" of tracks, in mm
T getUntrackedParameter(std::string const &, T const &) const
static const double y_width_
Definition: RPTopology.h:55
bool makeHepMC_
whether a HepMC description of the proton shall be saved in the event
Distribution angular_dist_
angular parameters in rad
bool roundToPitch_
whether measurement values shall be rounded to the nearest strip
static const double last_strip_to_border_dist_
Definition: RPTopology.h:57
unsigned int particlesPerEvent_
number of particles to generate per event
double particle_E_
particle energy and momentum
std::vector< unsigned int > RPs_
the list of RPs to simulate
double stripZeroPosition_
v position of strip 0, in mm
static const double pitch_
Definition: RPTopology.h:51
double insensitiveMarginStrips_
size of insensitive margin at sensor&#39;s edge facing the beam, in mm
edm::ESGetToken< CTPPSGeometry, VeryForwardMisalignedGeometryRecord > esTokenGeometry_

◆ ~PPSFastLocalSimulation()

PPSFastLocalSimulation::~PPSFastLocalSimulation ( )
override

Definition at line 227 of file PPSFastLocalSimulation.cc.

227 {}

Member Function Documentation

◆ GenerateTrack()

void PPSFastLocalSimulation::GenerateTrack ( unsigned int  pi,
CLHEP::HepRandomEngine &  rndEng,
HepMC::GenEvent gEv,
std::unique_ptr< edm::DetSetVector< TotemRPRecHit >> &  stripHitColl,
std::unique_ptr< edm::DetSetVector< CTPPSDiamondRecHit >> &  diamondHitColl,
std::unique_ptr< edm::DetSetVector< CTPPSPixelRecHit >> &  pixelHitColl,
const CTPPSGeometry geometry 
)
protected

Definition at line 231 of file PPSFastLocalSimulation.cc.

References A, angular_dist_, B, nano_mu_digi_cff::bx, CTPPSDiamondDetId::channel(), protons_cff::decRPId, hcalRecHitTable_cff::detId, spr::find(), PPSFastLocalSimulation::Distribution::Generate(), GenParticle::GenParticle, hfClusterShapes_cfi::hits, heavyIonCSV_trainingSettings::idx, insensitiveMarginStrips_, createfilelist::int, RPTopology::IsHit(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, visualization-live-secondInstance_cfg::m, makeHepMC_, makeHits_, particle_E_, particle_p_, pitchDiamonds_, pitchPixels_, pitchStrips_, position_dist_, printId(), roundToPitch_, RPs_, CTPPSDetId::sdTimingDiamond, CTPPSDetId::sdTrackingPixel, CTPPSDetId::sdTrackingStrip, mathSSE::sqrt(), nano_mu_digi_cff::strip, stripZeroPosition_, findQualityFiles::v, verbosity_, ApeEstimator_cff::width, and z0_.

Referenced by produce().

237  {
238  // generate track
239  double bx = 0., by = 0., ax = 0., ay = 0.;
240  position_dist_.Generate(rndEng, bx, by);
241  angular_dist_.Generate(rndEng, ax, ay);
242 
243  if (verbosity_ > 5)
244  printf("\tax = %.3f mrad, bx = %.3f mm, ay = %.3f mrad, by = %.3f mm, z0 = %.3f m\n",
245  ax * 1E3,
246  bx,
247  ay * 1E3,
248  by,
249  z0_ * 1E-3);
250 
251  // add HepMC track description
252  if (makeHepMC_) {
253  GenVertex *gVx = new GenVertex(HepMC::FourVector(bx, by, z0_, 0.));
254  gEv->add_vertex(gVx);
255 
256  GenParticle *gPe;
257  double az = sqrt(1. - ax * ax - ay * ay);
258  gPe = new GenParticle(HepMC::FourVector(particle_p_ * ax, particle_p_ * ay, particle_p_ * az, particle_E_),
259  2212,
260  1); // add a proton in final state
261  gPe->suggest_barcode(idx + 1);
262  gVx->add_particle_out(gPe);
263  }
264 
265  if (makeHits_) {
266  // check all sensors known to geometry
267  for (CTPPSGeometry::mapType::const_iterator it = geometry.beginSensor(); it != geometry.endSensor(); ++it) {
268  // get RP decimal id
269  CTPPSDetId detId(it->first);
270  unsigned int decRPId = detId.arm() * 100 + detId.station() * 10 + detId.rp();
271 
272  // stop if the RP is not selected
273  if (find(RPs_.begin(), RPs_.end(), decRPId) == RPs_.end())
274  continue;
275 
276  // keep only 1 diamond channel to represent 1 plane
277  if (detId.subdetId() == CTPPSDetId::sdTimingDiamond) {
278  CTPPSDiamondDetId channelId(it->first);
279  if (channelId.channel() != 0)
280  continue;
281  }
282 
283  if (verbosity_ > 5) {
284  printf(" ");
285  printId(it->first);
286  printf(": ");
287  }
288 
289  // determine the track impact point (in global coordinates)
290  // !! this assumes that local axes (1, 0, 0) and (0, 1, 0) describe the sensor surface
291  const auto gl_o = geometry.localToGlobal(detId, CTPPSGeometry::Vector(0, 0, 0));
292  const auto gl_a1 = geometry.localToGlobal(detId, CTPPSGeometry::Vector(1, 0, 0)) - gl_o;
293  const auto gl_a2 = geometry.localToGlobal(detId, CTPPSGeometry::Vector(0, 1, 0)) - gl_o;
294 
295  TMatrixD A(3, 3);
296  TVectorD B(3);
297  A(0, 0) = ax;
298  A(0, 1) = -gl_a1.x();
299  A(0, 2) = -gl_a2.x();
300  B(0) = gl_o.x() - bx;
301  A(1, 0) = ay;
302  A(1, 1) = -gl_a1.y();
303  A(1, 2) = -gl_a2.y();
304  B(1) = gl_o.y() - by;
305  A(2, 0) = 1.;
306  A(2, 1) = -gl_a1.z();
307  A(2, 2) = -gl_a2.z();
308  B(2) = gl_o.z() - z0_;
309  TMatrixD Ai(3, 3);
310  Ai = A.Invert();
311  TVectorD P(3);
312  P = Ai * B;
313 
314  double de_z = P(0);
315  CTPPSGeometry::Vector h_glo(ax * de_z + bx, ay * de_z + by, de_z + z0_);
316 
317  // hit in local coordinates
318  CTPPSGeometry::Vector h_loc = geometry.globalToLocal(detId, h_glo);
319 
320  // strips
321  if (detId.subdetId() == CTPPSDetId::sdTrackingStrip) {
322  double u = h_loc.x();
323  double v = h_loc.y();
324 
325  if (verbosity_ > 5)
326  printf(" u=%+8.4f, v=%+8.4f", u, v);
327 
328  // is it within detector?
330  if (verbosity_ > 5)
331  printf(" | no hit\n");
332  continue;
333  }
334 
335  // round the measurement
336  if (roundToPitch_) {
337  double m = stripZeroPosition_ - v;
338  signed int strip = (int)floor(m / pitchStrips_ + 0.5);
339 
341 
342  if (verbosity_ > 5)
343  printf(" | strip=%+4i", strip);
344  }
345 
346  double sigma = pitchStrips_ / sqrt(12.);
347 
348  if (verbosity_ > 5)
349  printf(" | m=%+8.4f, sigma=%+8.4f\n", v, sigma);
350 
351  DetSet<TotemRPRecHit> &hits = stripHitColl->find_or_insert(detId);
352  hits.emplace_back(v, sigma);
353  }
354 
355  // diamonds
356  if (detId.subdetId() == CTPPSDetId::sdTimingDiamond) {
357  if (roundToPitch_) {
358  h_loc.SetX(pitchDiamonds_ * floor(h_loc.x() / pitchDiamonds_ + 0.5));
359  }
360 
361  if (verbosity_ > 5)
362  printf(" m = %.3f\n", h_loc.x());
363 
364  const double width = pitchDiamonds_;
365 
366  DetSet<CTPPSDiamondRecHit> &hits = diamondHitColl->find_or_insert(detId);
367  hits.emplace_back(h_loc.x(), width, 0., 0., 0., 0., 0., 0., 0., 0, HPTDCErrorFlags(), false);
368  }
369 
370  // pixels
371  if (detId.subdetId() == CTPPSDetId::sdTrackingPixel) {
372  if (roundToPitch_) {
373  h_loc.SetX(pitchPixels_ * floor(h_loc.x() / pitchPixels_ + 0.5));
374  h_loc.SetY(pitchPixels_ * floor(h_loc.y() / pitchPixels_ + 0.5));
375  }
376 
377  if (verbosity_ > 5)
378  printf(" m1 = %.3f, m2 = %.3f\n", h_loc.x(), h_loc.y());
379 
380  const double sigma = pitchPixels_ / sqrt(12.);
381 
382  const LocalPoint lp(h_loc.x(), h_loc.y(), h_loc.z());
383  const LocalError le(sigma, 0., sigma);
384 
385  DetSet<CTPPSPixelRecHit> &hits = pixelHitColl->find_or_insert(detId);
386  hits.emplace_back(lp, le);
387  }
388  }
389  }
390 }
Definition: APVGainStruct.h:7
unsigned int verbosity_
verbosity level
Distribution position_dist_
position parameters in mm
bool makeHits_
whether the hits of the proton shall be calculated and saved
static bool IsHit(double u, double v, double insensitiveMargin=0)
Definition: RPTopology.cc:29
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
double z0_
the "origin" of tracks, in mm
DetGeomDesc::Translation Vector
Definition: CTPPSGeometry.h:36
bool makeHepMC_
whether a HepMC description of the proton shall be saved in the event
T sqrt(T t)
Definition: SSEVec.h:19
Distribution angular_dist_
angular parameters in rad
void Generate(CLHEP::HepRandomEngine &rndEng, double &x, double &y)
bool roundToPitch_
whether measurement values shall be rounded to the nearest strip
void printId(unsigned int id)
Definition: Utilities.cc:26
double particle_E_
particle energy and momentum
std::vector< unsigned int > RPs_
the list of RPs to simulate
std::pair< OmniClusterRef, TrackingParticleRef > P
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
double stripZeroPosition_
v position of strip 0, in mm
Detector ID class for CTPPS Timing Diamond detectors. Bits [19:31] : Assigend in CTPPSDetId Calss Bit...
Definition: APVGainStruct.h:7
double insensitiveMarginStrips_
size of insensitive margin at sensor&#39;s edge facing the beam, in mm

◆ produce()

void PPSFastLocalSimulation::produce ( edm::Event event,
const edm::EventSetup es 
)
overrideprotected

Definition at line 394 of file PPSFastLocalSimulation.cc.

References esTokenGeometry_, GenerateTrack(), edm::EventSetup::getData(), edm::RandomNumberGenerator::getEngine(), ZgammaFilter_cfi::HepMCProduct, makeHepMC_, makeHits_, eostools::move(), particlesPerEvent_, pi, and verbosity_.

394  {
395  if (verbosity_ > 2)
396  printf(">> PPSFastLocalSimulation::produce > event %llu\n", event.id().event());
397 
399  HepRandomEngine &rndEng = rng->getEngine(event.streamID());
400 
401  if (verbosity_ > 5)
402  printf("\tseed = %li\n", rndEng.getSeed());
403 
404  // get geometry
405  auto const &geometry = es.getData(esTokenGeometry_);
406 
407  // initialize products
408  GenEvent *gEv = new GenEvent();
409  gEv->set_event_number(event.id().event());
410 
411  unique_ptr<DetSetVector<TotemRPRecHit>> stripHitColl(new DetSetVector<TotemRPRecHit>());
412  unique_ptr<DetSetVector<CTPPSDiamondRecHit>> diamondHitColl(new DetSetVector<CTPPSDiamondRecHit>());
413  unique_ptr<DetSetVector<CTPPSPixelRecHit>> pixelHitColl(new DetSetVector<CTPPSPixelRecHit>());
414 
415  // run particle loop
416  for (unsigned int pi = 0; pi < particlesPerEvent_; pi++) {
417  if (verbosity_ > 5)
418  printf(" generating track %u\n", pi);
419 
420  GenerateTrack(pi, rndEng, gEv, stripHitColl, diamondHitColl, pixelHitColl, geometry);
421  }
422 
423  // save products
424  if (makeHepMC_) {
425  unique_ptr<HepMCProduct> hepMCoutput(new HepMCProduct());
426  hepMCoutput->addHepMCData(gEv);
427  event.put(std::move(hepMCoutput));
428  }
429 
430  if (makeHits_) {
431  event.put(std::move(stripHitColl));
432  event.put(std::move(diamondHitColl));
433  event.put(std::move(pixelHitColl));
434  }
435 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
unsigned int verbosity_
verbosity level
bool makeHits_
whether the hits of the proton shall be calculated and saved
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
const Double_t pi
bool makeHepMC_
whether a HepMC description of the proton shall be saved in the event
unsigned int particlesPerEvent_
number of particles to generate per event
void GenerateTrack(unsigned int pi, CLHEP::HepRandomEngine &rndEng, HepMC::GenEvent *gEv, std::unique_ptr< edm::DetSetVector< TotemRPRecHit >> &stripHitColl, std::unique_ptr< edm::DetSetVector< CTPPSDiamondRecHit >> &diamondHitColl, std::unique_ptr< edm::DetSetVector< CTPPSPixelRecHit >> &pixelHitColl, const CTPPSGeometry &geometry)
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
edm::ESGetToken< CTPPSGeometry, VeryForwardMisalignedGeometryRecord > esTokenGeometry_

Member Data Documentation

◆ angular_dist_

Distribution PPSFastLocalSimulation::angular_dist_
protected

angular parameters in rad

Definition at line 94 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ esTokenGeometry_

edm::ESGetToken<CTPPSGeometry, VeryForwardMisalignedGeometryRecord> PPSFastLocalSimulation::esTokenGeometry_
protected

Definition at line 101 of file PPSFastLocalSimulation.cc.

Referenced by produce().

◆ insensitiveMarginStrips_

double PPSFastLocalSimulation::insensitiveMarginStrips_
protected

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

Definition at line 78 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ makeHepMC_

bool PPSFastLocalSimulation::makeHepMC_
protected

whether a HepMC description of the proton shall be saved in the event

Definition at line 54 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack(), PPSFastLocalSimulation(), and produce().

◆ makeHits_

bool PPSFastLocalSimulation::makeHits_
protected

whether the hits of the proton shall be calculated and saved

Definition at line 57 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack(), PPSFastLocalSimulation(), and produce().

◆ particle_E_

double PPSFastLocalSimulation::particle_E_
protected

particle energy and momentum

Definition at line 66 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ particle_p_

double PPSFastLocalSimulation::particle_p_
protected

Definition at line 66 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ particlesPerEvent_

unsigned int PPSFastLocalSimulation::particlesPerEvent_
protected

number of particles to generate per event

Definition at line 63 of file PPSFastLocalSimulation.cc.

Referenced by produce().

◆ pitchDiamonds_

double PPSFastLocalSimulation::pitchDiamonds_
protected

Definition at line 75 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ pitchPixels_

double PPSFastLocalSimulation::pitchPixels_
protected

Definition at line 75 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ pitchStrips_

double PPSFastLocalSimulation::pitchStrips_
protected

in mm

Definition at line 75 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ position_dist_

Distribution PPSFastLocalSimulation::position_dist_
protected

position parameters in mm

Definition at line 91 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ roundToPitch_

bool PPSFastLocalSimulation::roundToPitch_
protected

whether measurement values shall be rounded to the nearest strip

Definition at line 72 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ RPs_

std::vector<unsigned int> PPSFastLocalSimulation::RPs_
protected

the list of RPs to simulate

Definition at line 60 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().

◆ stripZeroPosition_

double PPSFastLocalSimulation::stripZeroPosition_
protected

v position of strip 0, in mm

Definition at line 99 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack(), and PPSFastLocalSimulation().

◆ verbosity_

unsigned int PPSFastLocalSimulation::verbosity_
protected

verbosity level

Definition at line 51 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack(), and produce().

◆ z0_

double PPSFastLocalSimulation::z0_
protected

the "origin" of tracks, in mm

Definition at line 69 of file PPSFastLocalSimulation.cc.

Referenced by GenerateTrack().