CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 edm::EventSetup &iSetup)
13  : det_id_(det_id) {
14  verbosity_ = params.getParameter<int>("RPVerbosity");
16  theNoiseInElectrons = params.getParameter<double>("RPEquivalentNoiseCharge300um");
17  theStripThresholdInE = params.getParameter<double>("RPVFATThreshold");
18  noNoise_ = params.getParameter<bool>("RPNoNoise");
19  misalignment_simulation_on_ = params.getParameter<bool>("RPDisplacementOn");
20  links_persistence_ = params.getParameter<bool>("RPDigiSimHitRelationsPresistence");
21 
22  theRPGaussianTailNoiseAdder = std::make_unique<RPGaussianTailNoiseAdder>(
24  theRPPileUpSignals = std::make_unique<RPPileUpSignals>(params, det_id_);
25  theRPVFATSimulator = std::make_unique<RPVFATSimulator>(params, det_id_);
26  theRPHitChargeConverter = std::make_unique<RPHitChargeConverter>(params, eng, det_id_);
27  theRPDisplacementGenerator = std::make_unique<RPDisplacementGenerator>(params, det_id_, iSetup);
28 }
29 
30 void RPDetDigitizer::run(const std::vector<PSimHit> &input,
31  const std::vector<int> &input_links,
32  std::vector<TotemRPDigi> &output_digi,
33  simromanpot::DigiPrimaryMapType &output_digi_links) {
34  if (verbosity_)
35  LogDebug("RPDetDigitizer ") << det_id_ << " received input.size()=" << input.size() << "\n";
36  theRPPileUpSignals->reset();
37 
38  bool links_persistence_checked = links_persistence_ && input_links.size() == input.size();
39 
40  int input_size = input.size();
41  for (int i = 0; i < input_size; ++i) {
42  simromanpot::strip_charge_map the_strip_charge_map;
44  the_strip_charge_map = theRPHitChargeConverter->processHit(theRPDisplacementGenerator->displace(input[i]));
45  else
46  the_strip_charge_map = theRPHitChargeConverter->processHit(input[i]);
47 
48  if (verbosity_)
49  LogDebug("RPHitChargeConverter ") << det_id_ << " returned hits=" << the_strip_charge_map.size() << "\n";
50  if (links_persistence_checked)
51  theRPPileUpSignals->add(the_strip_charge_map, input_links[i]);
52  else
53  theRPPileUpSignals->add(the_strip_charge_map, 0);
54  }
55 
56  const simromanpot::strip_charge_map &theSignal = theRPPileUpSignals->dumpSignal();
57  simromanpot::strip_charge_map_links_type &theSignalProvenance = theRPPileUpSignals->dumpLinks();
59  if (noNoise_)
60  afterNoise = theSignal;
61  else
62  afterNoise = theRPGaussianTailNoiseAdder->addNoise(theSignal);
63 
64  theRPVFATSimulator->ConvertChargeToHits(afterNoise, theSignalProvenance, output_digi, output_digi_links);
65 }
double theNoiseInElectrons
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:22
std::map< unsigned short, std::vector< std::pair< int, double > > > strip_charge_map_links_type
Definition: RPSimTypes.h:28
static std::string const input
Definition: EdmProvDump.cc:47
double theStripThresholdInE
std::unique_ptr< RPGaussianTailNoiseAdder > theRPGaussianTailNoiseAdder
std::unique_ptr< RPHitChargeConverter > theRPHitChargeConverter
std::unique_ptr< RPVFATSimulator > theRPVFATSimulator
uint32_t RPDetId
Definition: RPSimTypes.h:11
RPDetDigitizer(const edm::ParameterSet &params, CLHEP::HepRandomEngine &eng, RPDetId det_id, const edm::EventSetup &iSetup)
bool misalignment_simulation_on_
std::map< unsigned short, double > strip_charge_map
Definition: RPSimTypes.h:14
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Geometrical and topological information on RP silicon detector. Uses coordinate a frame with origin i...
Definition: RPTopology.h:19
unsigned short DetStripNo() const
Definition: RPTopology.h:42
std::unique_ptr< RPDisplacementGenerator > theRPDisplacementGenerator
#define LogDebug(id)