CMS 3D CMS Logo

SSDigitizerAlgorithm.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <cmath>
3 
5 
10 
11 // Geometry
15 
16 using namespace edm;
17 
18 namespace {
19  double nFactorial(int n);
20  double aScalingConstant(int N, int i);
21 } // namespace
22 namespace {
23  double nFactorial(int n) { return std::tgamma(n + 1); }
24  double aScalingConstant(int N, int i) {
25  return std::pow(-1, (double)i) * nFactorial(N) * nFactorial(N + 2) /
26  (nFactorial(N - i) * nFactorial(N + 2 - i) * nFactorial(i));
27  }
28 } // namespace
29 
31  if (use_LorentzAngle_DB_) { // Get Lorentz angle from DB record
32  es.get<SiPhase2OuterTrackerLorentzAngleSimRcd>().get(siPhase2OTLorentzAngle_);
33  }
34 
35  es.get<TrackerDigiGeometryRecord>().get(geom_);
36 }
37 
39  : Phase2TrackerDigitizerAlgorithm(conf.getParameter<ParameterSet>("AlgorithmCommon"),
40  conf.getParameter<ParameterSet>("SSDigitizerAlgorithm")),
41  hitDetectionMode_(conf.getParameter<ParameterSet>("SSDigitizerAlgorithm").getParameter<int>("HitDetectionMode")),
42  pulseShapeParameters_(conf.getParameter<ParameterSet>("SSDigitizerAlgorithm")
43  .getParameter<std::vector<double> >("PulseShapeParameters")),
44  deadTime_(conf.getParameter<ParameterSet>("SSDigitizerAlgorithm").getParameter<double>("CBCDeadTime")) {
45  pixelFlag_ = false;
46  LogDebug("SSDigitizerAlgorithm ") << "SSDigitizerAlgorithm constructed "
47  << "Configuration parameters:"
48  << "Threshold/Gain = "
49  << "threshold in electron Endcap = " << theThresholdInE_Endcap_
50  << "threshold in electron Barrel = " << theThresholdInE_Barrel_ << " "
51  << theElectronPerADC_ << " " << theAdcFullScale_ << " The delta cut-off is set to "
52  << tMax_ << " pix-inefficiency " << addPixelInefficiency_;
54 }
55 SSDigitizerAlgorithm::~SSDigitizerAlgorithm() { LogDebug("SSDigitizerAlgorithm") << "SSDigitizerAlgorithm deleted"; }
56 //
57 // -- Select the Hit for Digitization
58 //
59 bool SSDigitizerAlgorithm::select_hit(const PSimHit& hit, double tCorr, double& sigScale) const {
60  bool result = false;
62  result = select_hit_sampledMode(hit, tCorr, sigScale);
64  result = select_hit_latchedMode(hit, tCorr, sigScale);
65  else {
66  double toa = hit.tof() - tCorr;
67  result = (toa > theTofLowerCut_ && toa < theTofUpperCut_);
68  }
69  return result;
70 }
71 //
72 // -- Select Hits in Sampled Mode
73 //
74 bool SSDigitizerAlgorithm::select_hit_sampledMode(const PSimHit& hit, double tCorr, double& sigScale) const {
75  double toa = hit.tof() - tCorr;
76  double sampling_time = bx_time;
77 
78  DetId det_id = DetId(hit.detUnitId());
79  float theThresholdInE =
81 
82  sigScale = getSignalScale(sampling_time - toa);
83  return (sigScale * hit.energyLoss() / GeVperElectron_ > theThresholdInE);
84 }
85 //
86 // -- Select Hits in Hit Detection Mode
87 //
88 bool SSDigitizerAlgorithm::select_hit_latchedMode(const PSimHit& hit, double tCorr, double& sigScale) const {
89  float toa = hit.tof() - tCorr;
90  toa -= hit.eventId().bunchCrossing() * bx_time;
91 
92  float sampling_time = (-1) * (hit.eventId().bunchCrossing() + 1) * bx_time;
93 
94  DetId det_id = DetId(hit.detUnitId());
95  float theThresholdInE =
97 
98  bool lastPulse = true;
99  bool aboveThr = false;
100  for (float i = deadTime_; i <= bx_time; i++) {
101  sigScale = getSignalScale(sampling_time - toa + i);
102 
103  aboveThr = (sigScale * hit.energyLoss() / GeVperElectron_ > theThresholdInE);
104  if (!lastPulse && aboveThr)
105  return true;
106 
107  lastPulse = aboveThr;
108  }
109  return false;
110 }
112  constexpr size_t max_par = 6;
113  if (pulseShapeParameters_.size() < max_par)
114  return -1;
115  double xOffset = pulseShapeParameters_[0];
116  double tau = pulseShapeParameters_[1];
117  double r = pulseShapeParameters_[2];
118  double theta = pulseShapeParameters_[3];
119  int nTerms = static_cast<int>(pulseShapeParameters_[4]);
120 
121  double fN = 0;
122  double xx = x - xOffset;
123  if (xx < 0)
124  return 0;
125 
126  for (int i = 0; i < nTerms; i++) {
127  double angularTerm = 0;
128  double temporalTerm = 0;
129  double rTerm = std::pow(r, i) / (std::pow(tau, 2. * i) * ::nFactorial(i + 2));
130  for (int j = 0; j <= i; j++) {
131  angularTerm += std::pow(std::cos(theta), (double)(i - j)) * std::pow(std::sin(theta), (double)j);
132  temporalTerm += ::aScalingConstant(i, j) * std::pow(xx, (double)(i - j)) * std::pow(tau, (double)j);
133  }
134  double fi = rTerm * angularTerm * temporalTerm;
135 
136  fN += fi;
137  }
138  return fN;
139 }
140 double SSDigitizerAlgorithm::signalShape(double x) const {
141  double xOffset = pulseShapeParameters_[0];
142  double tau = pulseShapeParameters_[1];
143  double maxCharge = pulseShapeParameters_[5];
144 
145  double xx = x - xOffset;
146  return maxCharge * (std::exp(-xx / tau) * std::pow(xx / tau, 2.) * cbc3PulsePolarExpansion(x));
147 }
149  for (size_t i = 0; i < interpolationPoints; i++) {
150  float val = i / interpolationStep;
151 
152  pulseShapeVec_.push_back(signalShape(val));
153  }
154 }
155 double SSDigitizerAlgorithm::getSignalScale(double xval) const {
156  double res = 0.0;
157  int len = pulseShapeVec_.size();
158 
159  if (xval < 0.0 || xval * interpolationStep >= len)
160  return res;
161 
162  unsigned int lower = std::floor(xval) * interpolationStep;
163  unsigned int upper = std::ceil(xval) * interpolationStep;
164  for (size_t i = lower + 1; i < upper * interpolationStep; i++) {
165  float val = i * 0.1;
166  if (val > xval) {
167  res = pulseShapeVec_[i - 1];
168  break;
169  }
170  }
171  return res;
172 }
173 //
174 // -- Compare Signal with Threshold
175 //
177  float charge,
178  float thr) const {
179  return (charge >= thr);
180 }
SSDigitizerAlgorithm::storeSignalShape
void storeSignalShape()
Definition: SSDigitizerAlgorithm.cc:148
Phase2TrackerDigitizerAlgorithm::theThresholdInE_Barrel_
const float theThresholdInE_Barrel_
Definition: Phase2TrackerDigitizerAlgorithm.h:147
mps_fire.i
i
Definition: mps_fire.py:428
MessageLogger.h
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
SiPhase2OuterTrackerLorentzAngleSimRcd
Definition: SiPhase2OuterTrackerCondDataRecords.h:12
SSDigitizerAlgorithm::SSDigitizerAlgorithm
SSDigitizerAlgorithm(const edm::ParameterSet &conf)
Definition: SSDigitizerAlgorithm.cc:38
Phase2TrackerDigitizerAlgorithm::theTofUpperCut_
const float theTofUpperCut_
Definition: Phase2TrackerDigitizerAlgorithm.h:156
metsig::tau
Definition: SignAlgoResolutions.h:49
SSDigitizerAlgorithm::interpolationStep
static constexpr int interpolationStep
Definition: SSDigitizerAlgorithm.h:32
SSDigitizerAlgorithm::pulseShapeVec_
std::vector< double > pulseShapeVec_
Definition: SSDigitizerAlgorithm.h:27
edm
HLT enums.
Definition: AlignableModifier.h:19
SSDigitizerAlgorithm::isAboveThreshold
bool isAboveThreshold(const DigitizerUtility::SimHitInfo *hitInfo, float charge, float thr) const override
Definition: SSDigitizerAlgorithm.cc:176
SSDigitizerAlgorithm::interpolationPoints
static constexpr size_t interpolationPoints
Definition: SSDigitizerAlgorithm.h:31
SSDigitizerAlgorithm::bx_time
static constexpr float bx_time
Definition: SSDigitizerAlgorithm.h:30
SSDigitizerAlgorithm::~SSDigitizerAlgorithm
~SSDigitizerAlgorithm() override
Definition: SSDigitizerAlgorithm.cc:55
DDAxes::x
SiPhase2OuterTrackerLorentzAngleRcd.h
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Phase2TrackerDigitizerAlgorithm::theTofLowerCut_
const float theTofLowerCut_
Definition: Phase2TrackerDigitizerAlgorithm.h:155
DetId
Definition: DetId.h:17
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
edm::EventSetup::get
T get() const
Definition: EventSetup.h:87
reco::ceil
constexpr int32_t ceil(float num)
Definition: constexpr_cmath.h:7
Phase2TrackerDigitizerAlgorithm::theAdcFullScale_
const int theAdcFullScale_
Definition: Phase2TrackerDigitizerAlgorithm.h:141
SSDigitizerAlgorithm::signalShape
double signalShape(double x) const
Definition: SSDigitizerAlgorithm.cc:140
Service.h
TrackerDigiGeometryRecord
Definition: TrackerDigiGeometryRecord.h:15
SSDigitizerAlgorithm::pulseShapeParameters_
std::vector< double > pulseShapeParameters_
Definition: SSDigitizerAlgorithm.h:28
N
#define N
Definition: blowfish.cc:9
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
hcaldqm::quantity::fN
Definition: ValueQuantity.h:11
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
SSDigitizerAlgorithm::cbc3PulsePolarExpansion
double cbc3PulsePolarExpansion(double x) const
Definition: SSDigitizerAlgorithm.cc:111
SSDigitizerAlgorithm::select_hit_sampledMode
bool select_hit_sampledMode(const PSimHit &hit, double tCorr, double &sigScale) const
Definition: SSDigitizerAlgorithm.cc:74
TrackerDigiGeometryRecord.h
SSDigitizerAlgorithm.h
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
Phase2TrackerDigitizerAlgorithm
Definition: Phase2TrackerDigitizerAlgorithm.h:59
Phase2TrackerDigitizerAlgorithm::pixelFlag_
bool pixelFlag_
Definition: Phase2TrackerDigitizerAlgorithm.h:238
Phase2TrackerDigitizerAlgorithm::tMax_
const double tMax_
Definition: Phase2TrackerDigitizerAlgorithm.h:179
edm::get
T const & get(Event const &event, InputTag const &tag) noexcept(false)
Definition: Event.h:675
createfilelist.int
int
Definition: createfilelist.py:10
phase1PixelTopology::xOffset
constexpr int16_t xOffset
Definition: phase1PixelTopology.h:19
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
SSDigitizerAlgorithm::deadTime_
float deadTime_
Definition: SSDigitizerAlgorithm.h:29
SSDigitizerAlgorithm::select_hit
bool select_hit(const PSimHit &hit, double tCorr, double &sigScale) const override
Definition: SSDigitizerAlgorithm.cc:59
Phase2TrackerDigitizerAlgorithm::theElectronPerADC_
const float theElectronPerADC_
Definition: Phase2TrackerDigitizerAlgorithm.h:140
edm::EventSetup
Definition: EventSetup.h:58
Phase2TrackerDigitizerAlgorithm::GeVperElectron_
const float GeVperElectron_
Definition: Phase2TrackerDigitizerAlgorithm.h:126
res
Definition: Electron.h:6
Phase2TrackerDigitizerAlgorithm::theThresholdInE_Endcap_
const float theThresholdInE_Endcap_
Definition: Phase2TrackerDigitizerAlgorithm.h:146
alignCSCRings.r
r
Definition: alignCSCRings.py:93
SSDigitizerAlgorithm::select_hit_latchedMode
bool select_hit_latchedMode(const PSimHit &hit, double tCorr, double &sigScale) const
Definition: SSDigitizerAlgorithm.cc:88
DigitizerUtility::SimHitInfo
Definition: DigitizerUtility.h:14
heppy_batch.val
val
Definition: heppy_batch.py:351
std
Definition: JetResolutionObject.h:76
SSDigitizerAlgorithm::init
void init(const edm::EventSetup &es) override
Definition: SSDigitizerAlgorithm.cc:30
SSDigitizerAlgorithm::LatchedMode
Definition: SSDigitizerAlgorithm.h:18
SSDigitizerAlgorithm::SampledMode
Definition: SSDigitizerAlgorithm.h:18
SSDigitizerAlgorithm::hitDetectionMode_
int hitDetectionMode_
Definition: SSDigitizerAlgorithm.h:26
genVertex_cff.x
x
Definition: genVertex_cff.py:12
PixelGeomDetUnit.h
StripSubdetector::TOB
static constexpr auto TOB
Definition: StripSubdetector.h:18
pileupCalc.upper
upper
Definition: pileupCalc.py:214
SSDigitizerAlgorithm::getSignalScale
double getSignalScale(double xval) const
Definition: SSDigitizerAlgorithm.cc:155
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
mps_fire.result
result
Definition: mps_fire.py:311
ParameterSet.h
PSimHit
Definition: PSimHit.h:15
Phase2TrackerDigitizerAlgorithm::addPixelInefficiency_
const bool addPixelInefficiency_
Definition: Phase2TrackerDigitizerAlgorithm.h:166
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
StripSubdetector.h
geometryCSVtoXML.xx
xx
Definition: geometryCSVtoXML.py:19
hit
Definition: SiStripHitEffFromCalibTree.cc:88