Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:25:34

0001 #include <vector>
0002 #include <iostream>
0003 
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "SimPPS/RPDigiProducer/plugins/RPDetDigitizer.h"
0006 #include "Geometry/VeryForwardRPTopology/interface/RPTopology.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 
0009 RPDetDigitizer::RPDetDigitizer(const edm::ParameterSet &params,
0010                                CLHEP::HepRandomEngine &eng,
0011                                RPDetId det_id,
0012                                const CTPPSRPAlignmentCorrectionsData *alignments,
0013                                const CTPPSGeometry &geom)
0014 
0015     : det_id_(det_id) {
0016   verbosity_ = params.getParameter<int>("RPVerbosity");
0017   numStrips_ = RPTopology().DetStripNo();
0018   theNoiseInElectrons = params.getParameter<double>("RPEquivalentNoiseCharge300um");
0019   theStripThresholdInE = params.getParameter<double>("RPVFATThreshold");
0020   noNoise_ = params.getParameter<bool>("RPNoNoise");
0021   misalignment_simulation_on_ = params.getParameter<bool>("RPDisplacementOn");
0022   links_persistence_ = params.getParameter<bool>("RPDigiSimHitRelationsPresistence");
0023 
0024   theRPGaussianTailNoiseAdder = std::make_unique<RPGaussianTailNoiseAdder>(
0025       numStrips_, theNoiseInElectrons, theStripThresholdInE, eng, verbosity_);
0026   theRPPileUpSignals = std::make_unique<RPPileUpSignals>(params, det_id_);
0027   theRPVFATSimulator = std::make_unique<RPVFATSimulator>(params, det_id_);
0028   theRPHitChargeConverter = std::make_unique<RPHitChargeConverter>(params, eng, det_id_);
0029   theRPDisplacementGenerator = std::make_unique<RPDisplacementGenerator>(
0030       params.getParameter<bool>("RPDisplacementOn"), det_id_, alignments, geom);
0031 }
0032 
0033 void RPDetDigitizer::run(const std::vector<PSimHit> &input,
0034                          const std::vector<int> &input_links,
0035                          std::vector<TotemRPDigi> &output_digi,
0036                          simromanpot::DigiPrimaryMapType &output_digi_links) {
0037   if (verbosity_)
0038     LogDebug("RPDetDigitizer ") << det_id_ << " received input.size()=" << input.size() << "\n";
0039   theRPPileUpSignals->reset();
0040 
0041   bool links_persistence_checked = links_persistence_ && input_links.size() == input.size();
0042 
0043   int input_size = input.size();
0044   for (int i = 0; i < input_size; ++i) {
0045     simromanpot::strip_charge_map the_strip_charge_map;
0046     if (misalignment_simulation_on_)
0047       the_strip_charge_map = theRPHitChargeConverter->processHit(theRPDisplacementGenerator->displace(input[i]));
0048     else
0049       the_strip_charge_map = theRPHitChargeConverter->processHit(input[i]);
0050 
0051     if (verbosity_)
0052       LogDebug("RPHitChargeConverter ") << det_id_ << " returned hits=" << the_strip_charge_map.size() << "\n";
0053     if (links_persistence_checked)
0054       theRPPileUpSignals->add(the_strip_charge_map, input_links[i]);
0055     else
0056       theRPPileUpSignals->add(the_strip_charge_map, 0);
0057   }
0058 
0059   const simromanpot::strip_charge_map &theSignal = theRPPileUpSignals->dumpSignal();
0060   simromanpot::strip_charge_map_links_type &theSignalProvenance = theRPPileUpSignals->dumpLinks();
0061   simromanpot::strip_charge_map afterNoise;
0062   if (noNoise_)
0063     afterNoise = theSignal;
0064   else
0065     afterNoise = theRPGaussianTailNoiseAdder->addNoise(theSignal);
0066 
0067   theRPVFATSimulator->ConvertChargeToHits(afterNoise, theSignalProvenance, output_digi, output_digi_links);
0068 }