Main Page
Namespaces
Classes
Package Documentation
SimCalorimetry
CastorSim
src
CastorSimParameters.cc
Go to the documentation of this file.
1
#include "CLHEP/Random/RandGaussQ.h"
2
#include "
CondFormats/CastorObjects/interface/CastorGain.h
"
3
#include "
CondFormats/CastorObjects/interface/CastorGainWidth.h
"
4
#include "
DataFormats/HcalDetId/interface/HcalDetId.h
"
5
#include "
DataFormats/HcalDetId/interface/HcalGenericDetId.h
"
6
#include "
DataFormats/HcalDetId/interface/HcalSubdetector.h
"
7
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
8
#include "
SimCalorimetry/CastorSim/src/CastorSimParameters.h
"
9
10
CastorSimParameters::CastorSimParameters
(
double
simHitToPhotoelectrons
,
11
double
photoelectronsToAnalog
,
12
double
samplingFactor
,
13
double
timePhase
,
14
bool
syncPhase
)
15
:
CaloSimParameters
(
16
simHitToPhotoelectrons, photoelectronsToAnalog, samplingFactor, timePhase, 6, 4,
false
, syncPhase),
17
theDbService(
nullptr
),
18
theSamplingFactor(samplingFactor),
19
nominalfCperPE(1) {}
20
21
CastorSimParameters::CastorSimParameters
(
const
edm::ParameterSet
&
p
)
22
:
CaloSimParameters
(p),
23
theDbService
(
nullptr
),
24
theSamplingFactor
(p.getParameter<double>(
"samplingFactor"
)),
25
nominalfCperPE
(p.getParameter<double>(
"photoelectronsToAnalog"
)) {}
26
27
double
CastorSimParameters::getNominalfCperPE
()
const
{
28
// return the nominal PMT gain value of CASTOR from the config file.
29
return
nominalfCperPE
;
30
}
31
32
double
CastorSimParameters::photoelectronsToAnalog
(
const
DetId
&detId)
const
{
33
// calculate factor (PMT gain) using sampling factor value & available
34
// electron gain
35
return
theSamplingFactor
/
fCtoGeV
(detId);
36
}
37
38
double
CastorSimParameters::fCtoGeV
(
const
DetId
&detId)
const
{
39
assert(
theDbService
!=
nullptr
);
40
HcalGenericDetId
hcalGenDetId(detId);
41
const
CastorGain
*gains =
theDbService
->
getGain
(hcalGenDetId);
42
const
CastorGainWidth
*gwidths =
theDbService
->
getGainWidth
(hcalGenDetId);
43
double
result
= 0.0;
44
if
(!gains || !gwidths) {
45
edm::LogError
(
"CastorAmplifier"
) <<
"Could not fetch HCAL conditions for channel "
<< hcalGenDetId;
46
}
else
{
47
// only one gain will be recorded per channel, so just use capID 0 for now
48
result = gains->
getValue
(0);
49
}
50
return
result
;
51
}
MessageLogger.h
CastorGain
Definition:
CastorGain.h:14
CastorSimParameters::getNominalfCperPE
double getNominalfCperPE() const
Definition:
CastorSimParameters.cc:27
AlCaHLTBitMon_ParallelJobs.p
p
Definition:
AlCaHLTBitMon_ParallelJobs.py:153
mps_fire.result
result
Definition:
mps_fire.py:283
CastorSimParameters::theDbService
const CastorDbService * theDbService
Definition:
CastorSimParameters.h:28
funct::false
false
Definition:
Factorize.h:36
HcalSubdetector.h
nullptr
#define nullptr
Definition:
GCC11Compatibility.h:37
CastorGainWidth.h
CaloSimParameters
Main class for Parameters in different subdetectors.
Definition:
CaloSimParameters.h:14
CastorSimParameters::CastorSimParameters
CastorSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog, double samplingFactor, double timePhase, bool syncPhase)
Definition:
CastorSimParameters.cc:10
HcalGenericDetId.h
edm::LogError
Definition:
MessageLogger.h:183
HcalDetId.h
CastorSimParameters::theSamplingFactor
double theSamplingFactor
Definition:
CastorSimParameters.h:29
hcalSimParameters_cfi.photoelectronsToAnalog
photoelectronsToAnalog
Definition:
hcalSimParameters_cfi.py:18
hcalSimParameters_cfi.simHitToPhotoelectrons
simHitToPhotoelectrons
Definition:
hcalSimParameters_cfi.py:19
CastorDbService::getGain
const CastorGain * getGain(const HcalGenericDetId &fId) const
Definition:
CastorDbService.cc:155
CastorGainWidth
Definition:
CastorGainWidth.h:13
CastorSimParameters.h
CastorSimParameters::fCtoGeV
double fCtoGeV(const DetId &detId) const
Definition:
CastorSimParameters.cc:38
HcalGenericDetId
Definition:
HcalGenericDetId.h:15
DetId
Definition:
DetId.h:18
CastorGain.h
CastorDbService::getGainWidth
const CastorGainWidth * getGainWidth(const HcalGenericDetId &fId) const
Definition:
CastorDbService.cc:162
ecalSimParameterMap_cff.samplingFactor
samplingFactor
Definition:
ecalSimParameterMap_cff.py:8
edm::ParameterSet
Definition:
ParameterSet.h:36
CastorSimParameters::nominalfCperPE
double nominalfCperPE
Definition:
CastorSimParameters.h:30
CastorGain::getValue
float getValue(int fCapId) const
get value for capId = 0..3
Definition:
CastorGain.h:19
ecalSimParameterMap_cff.syncPhase
syncPhase
Definition:
ecalSimParameterMap_cff.py:11
CaloSimParameters::photoelectronsToAnalog
double photoelectronsToAnalog() const
the factor which goes from photoelectrons to whatever gets read by ADCs
Definition:
CaloSimParameters.h:38
ecalSimParameterMap_cff.timePhase
timePhase
Definition:
ecalSimParameterMap_cff.py:9
Generated for CMSSW Reference Manual by
1.8.11