CMS 3D CMS Logo

RPDetDigitizer.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <iostream>
3 
8 
10  CLHEP::HepRandomEngine &eng,
11  RPDetId det_id,
12  const CTPPSRPAlignmentCorrectionsData *alignments,
13  const CTPPSGeometry &geom)
14 
15  : det_id_(det_id) {
16  verbosity_ = params.getParameter<int>("RPVerbosity");
18  theNoiseInElectrons = params.getParameter<double>("RPEquivalentNoiseCharge300um");
19  theStripThresholdInE = params.getParameter<double>("RPVFATThreshold");
20  noNoise_ = params.getParameter<bool>("RPNoNoise");
21  misalignment_simulation_on_ = params.getParameter<bool>("RPDisplacementOn");
22  links_persistence_ = params.getParameter<bool>("RPDigiSimHitRelationsPresistence");
23 
24  theRPGaussianTailNoiseAdder = std::make_unique<RPGaussianTailNoiseAdder>(
26  theRPPileUpSignals = std::make_unique<RPPileUpSignals>(params, det_id_);
27  theRPVFATSimulator = std::make_unique<RPVFATSimulator>(params, det_id_);
28  theRPHitChargeConverter = std::make_unique<RPHitChargeConverter>(params, eng, det_id_);
29  theRPDisplacementGenerator = std::make_unique<RPDisplacementGenerator>(
30  params.getParameter<bool>("RPDisplacementOn"), det_id_, alignments, geom);
31 }
32 
33 void RPDetDigitizer::run(const std::vector<PSimHit> &input,
34  const std::vector<int> &input_links,
35  std::vector<TotemRPDigi> &output_digi,
36  simromanpot::DigiPrimaryMapType &output_digi_links) {
37  if (verbosity_)
38  LogDebug("RPDetDigitizer ") << det_id_ << " received input.size()=" << input.size() << "\n";
39  theRPPileUpSignals->reset();
40 
41  bool links_persistence_checked = links_persistence_ && input_links.size() == input.size();
42 
43  int input_size = input.size();
44  for (int i = 0; i < input_size; ++i) {
45  simromanpot::strip_charge_map the_strip_charge_map;
47  the_strip_charge_map = theRPHitChargeConverter->processHit(theRPDisplacementGenerator->displace(input[i]));
48  else
49  the_strip_charge_map = theRPHitChargeConverter->processHit(input[i]);
50 
51  if (verbosity_)
52  LogDebug("RPHitChargeConverter ") << det_id_ << " returned hits=" << the_strip_charge_map.size() << "\n";
53  if (links_persistence_checked)
54  theRPPileUpSignals->add(the_strip_charge_map, input_links[i]);
55  else
56  theRPPileUpSignals->add(the_strip_charge_map, 0);
57  }
58 
59  const simromanpot::strip_charge_map &theSignal = theRPPileUpSignals->dumpSignal();
60  simromanpot::strip_charge_map_links_type &theSignalProvenance = theRPPileUpSignals->dumpLinks();
62  if (noNoise_)
63  afterNoise = theSignal;
64  else
65  afterNoise = theRPGaussianTailNoiseAdder->addNoise(theSignal);
66 
67  theRPVFATSimulator->ConvertChargeToHits(afterNoise, theSignalProvenance, output_digi, output_digi_links);
68 }
double theNoiseInElectrons
RPDetDigitizer(const edm::ParameterSet &params, CLHEP::HepRandomEngine &eng, RPDetId det_id, const CTPPSRPAlignmentCorrectionsData *alignments, const CTPPSGeometry &geom)
std::unique_ptr< RPPileUpSignals > theRPPileUpSignals
void run(const std::vector< PSimHit > &input, const std::vector< int > &input_links, std::vector< TotemRPDigi > &output_digi, simromanpot::DigiPrimaryMapType &output_digi_links)
std::vector< std::vector< std::pair< int, double > > > DigiPrimaryMapType
Definition: RPSimTypes.h:23
std::map< unsigned short, std::vector< std::pair< int, double > > > strip_charge_map_links_type
Definition: RPSimTypes.h:29
static std::string const input
Definition: EdmProvDump.cc:50
double theStripThresholdInE
std::unique_ptr< RPGaussianTailNoiseAdder > theRPGaussianTailNoiseAdder
std::unique_ptr< RPHitChargeConverter > theRPHitChargeConverter
std::unique_ptr< RPVFATSimulator > theRPVFATSimulator
uint32_t RPDetId
Definition: RPSimTypes.h:12
The manager class for TOTEM RP geometry.
Definition: CTPPSGeometry.h:30
bool misalignment_simulation_on_
std::map< unsigned short, double > strip_charge_map
Definition: RPSimTypes.h:15
Container for CTPPS RP alignment corrections. The corrections are stored on two levels - RP and senso...
Geometrical and topological information on RP silicon detector. Uses coordinate a frame with origin i...
Definition: RPTopology.h:19
std::unique_ptr< RPDisplacementGenerator > theRPDisplacementGenerator
unsigned short DetStripNo() const
Definition: RPTopology.h:42
#define LogDebug(id)