FastSimulation
ShowerDevelopment
src
FastHFShowerLibrary.cc
Go to the documentation of this file.
1
// File: FastHFShowerLibrary.cc
3
// Description: Shower library for Very forward hadron calorimeter
5
6
#include "
FastSimulation/ShowerDevelopment/interface/FastHFShowerLibrary.h
"
7
#include "
FastSimulation/Event/interface/FSimEvent.h
"
8
#include "
FastSimulation/Event/interface/FSimTrack.h
"
9
#include "
FastSimulation/Utilities/interface/RandomEngineAndDistribution.h
"
10
#include "
SimG4CMS/Calo/interface/HFFibreFiducial.h
"
11
#include "
FWCore/Framework/interface/ESTransientHandle.h
"
12
#include "
FWCore/Framework/interface/ESHandle.h
"
13
#include "
Geometry/Records/interface/IdealGeometryRecord.h
"
14
#include "
Geometry/Records/interface/HcalParametersRcd.h
"
15
#include "
Geometry/Records/interface/HcalSimNumberingRecord.h
"
16
#include "
Geometry/HcalCommonData/interface/HcalDDDSimConstants.h
"
17
#include "
Geometry/HcalCommonData/interface/HcalSimulationConstants.h
"
18
#include "
FWCore/Utilities/interface/Exception.h
"
19
#include "
DataFormats/HcalDetId/interface/HcalDetId.h
"
20
#include "
DataFormats/HcalDetId/interface/HcalSubdetector.h
"
21
22
#include "Randomize.hh"
23
#include "CLHEP/Units/GlobalSystemOfUnits.h"
24
#include "CLHEP/Units/GlobalPhysicalConstants.h"
25
26
// Geant4 headers
27
#include "G4ParticleDefinition.hh"
28
#include "G4DynamicParticle.hh"
29
#include "G4DecayPhysics.hh"
30
#include "G4ParticleTable.hh"
31
#include "G4ParticleTypes.hh"
32
33
// STL headers
34
#include <memory>
35
36
#include <iostream>
37
#include <mutex>
38
#include <vector>
39
40
//#define DebugLog
41
42
static
std::once_flag
initializeOnce
;
43
44
FastHFShowerLibrary::FastHFShowerLibrary
(
edm::ParameterSet
const
&
p
) : fast(
p
) {
45
edm::ParameterSet
m_HS =
p
.getParameter<
edm::ParameterSet
>(
"HFShowerLibrary"
);
46
applyFidCut
= m_HS.
getParameter
<
bool
>(
"ApplyFiducialCut"
);
47
}
48
49
void
const
FastHFShowerLibrary::initHFShowerLibrary
(
const
edm::EventSetup
& iSetup) {
50
edm::LogInfo
(
"FastCalorimetry"
) <<
"initHFShowerLibrary::initialization"
;
51
52
edm::ESHandle<HcalDDDSimConstants>
hdc;
53
iSetup.
get
<
HcalSimNumberingRecord
>().
get
(hdc);
54
hcalConstants
= hdc.
product
();
55
56
edm::ESHandle<HcalSimulationConstants>
hdsc;
57
iSetup.
get
<
HcalSimNumberingRecord
>().
get
(hdsc);
58
const
HcalSimulationConstants
* hsps = hdsc.
product
();
59
60
std::string
name
=
"HcalHits"
;
61
numberingFromDDD
= std::make_unique<HcalNumberingFromDDD>(
hcalConstants
);
62
hfshower
= std::make_unique<HFShowerLibrary>(
name
,
hcalConstants
, hsps->
hcalsimpar
(),
fast
);
63
64
//only one thread can be allowed to setup the G4 physics table.
65
std::call_once(
initializeOnce
, []() {
66
// Geant4 particles
67
G4DecayPhysics
decays
;
68
decays
.ConstructParticle();
69
G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
70
partTable->SetReadiness();
71
});
72
//G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
73
//hfshower->initRun(partTable, hcalConstants); // init particle code
74
}
75
76
void
FastHFShowerLibrary::SetRandom
(
const
RandomEngineAndDistribution
* rnd) {
77
// define Geant4 engine per thread
78
G4Random::setTheEngine(&(rnd->
theEngine
()));
79
LogDebug
(
"FastHFShowerLibrary::recoHFShowerLibrary"
)
80
<<
"Begin of event "
<< G4UniformRand() <<
" "
<< rnd->
theEngine
().name() <<
" "
<< rnd->
theEngine
();
81
}
82
83
void
FastHFShowerLibrary::recoHFShowerLibrary
(
const
FSimTrack
& myTrack) {
84
#ifdef DebugLog
85
edm::LogInfo
(
"FastCalorimetry"
) <<
"FastHFShowerLibrary: recoHFShowerLibrary "
;
86
#endif
87
88
if
(!myTrack.
onVFcal
()) {
89
#ifdef DebugLog
90
edm::LogInfo
(
"FastCalorimetry"
) <<
"FastHFShowerLibrary: we should not be here "
;
91
#endif
92
}
93
94
hitMap
.clear();
95
double
eGen = 1000. * myTrack.
vfcalEntrance
().
e
();
// energy in [MeV]
96
double
delZv = (myTrack.
vfcalEntrance
().
vertex
().Z() > 0.0) ? 50.0 : -50.0;
97
G4ThreeVector
vertex
(10. * myTrack.
vfcalEntrance
().
vertex
().X(),
98
10. * myTrack.
vfcalEntrance
().
vertex
().Y(),
99
10. * myTrack.
vfcalEntrance
().
vertex
().Z() + delZv);
// in [mm]
100
101
G4ThreeVector direction(
102
myTrack.
vfcalEntrance
().
Vect
().X(), myTrack.
vfcalEntrance
().
Vect
().Y(), myTrack.
vfcalEntrance
().
Vect
().Z());
103
104
bool
ok
;
105
double
weight
= 1.0;
// rad. damage
106
int
parCode = myTrack.
type
();
107
double
tSlice = 0.1 *
vertex
.mag() / 29.98;
108
109
std::vector<HFShowerLibrary::Hit>
hits
=
110
hfshower
->fillHits(
vertex
, direction, parCode, eGen,
ok
,
weight
, tSlice,
false
);
111
112
for
(
unsigned
int
i
= 0;
i
<
hits
.size(); ++
i
) {
113
G4ThreeVector
pos
=
hits
[
i
].position;
114
int
depth
=
hits
[
i
].depth;
115
double
time
=
hits
[
i
].time;
116
if
(!
applyFidCut
|| (
HFFibreFiducial::PMTNumber
(
pos
) > 0)) {
117
// if (!applyFidCut || (applyFidCut && HFFibreFiducial::PMTNumber(pos)>0)) {
118
int
det = 5;
119
int
lay = 1;
120
uint32_t
id
= 0;
121
HcalNumberingFromDDD::HcalID
tmp
=
122
numberingFromDDD
->unitID(det,
math::XYZVectorD
(
pos
.x(),
pos
.y(),
pos
.z()),
depth
, lay);
123
modifyDepth
(
tmp
);
124
id
=
numberingScheme
.
getUnitID
(
tmp
);
125
126
CaloHitID
current_id(
id
,
time
, myTrack.
id
());
127
std::map<CaloHitID, float>::iterator cellitr;
128
cellitr =
hitMap
.find(current_id);
129
if
(cellitr ==
hitMap
.end()) {
130
hitMap
.insert(std::pair<CaloHitID, float>(current_id, 1.0));
131
}
else
{
132
cellitr->second += 1.0;
133
}
134
}
// end of isItinFidVolume check
135
}
// end loop over hits
136
}
137
138
void
FastHFShowerLibrary::modifyDepth
(
HcalNumberingFromDDD::HcalID
&
id
) {
139
if
(
id
.subdet ==
HcalForward
) {
140
int
ieta
= (
id
.zside == 0) ? -
id
.etaR :
id
.etaR;
141
if
(
hcalConstants
->
maxHFDepth
(
ieta
,
id
.phis) > 2) {
142
if
(
id
.
depth
<= 2) {
143
if
(G4UniformRand() > 0.5)
144
id
.depth += 2;
145
}
146
}
147
}
148
}
FastHFShowerLibrary::hcalConstants
const HcalDDDSimConstants * hcalConstants
Definition:
FastHFShowerLibrary.h:56
edm::ESHandle::product
T const * product() const
Definition:
ESHandle.h:86
FastHFShowerLibrary.h
FastHFShowerLibrary::fast
const edm::ParameterSet fast
Definition:
FastHFShowerLibrary.h:53
HcalNumberingScheme::getUnitID
virtual uint32_t getUnitID(const HcalNumberingFromDDD::HcalID &id)
Definition:
HcalNumberingScheme.cc:19
bJpsiMuMuTrigSettings_cff.decays
decays
Definition:
bJpsiMuMuTrigSettings_cff.py:15
mps_fire.i
i
Definition:
mps_fire.py:428
ESTransientHandle.h
HcalNumberingFromDDD::HcalID
Definition:
HcalNumberingFromDDD.h:21
hfClusterShapes_cfi.hits
hits
Definition:
hfClusterShapes_cfi.py:5
ESHandle.h
HcalSimNumberingRecord.h
RawParticle::vertex
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition:
RawParticle.h:320
FastHFShowerLibrary::recoHFShowerLibrary
void recoHFShowerLibrary(const FSimTrack &myTrack)
Definition:
FastHFShowerLibrary.cc:83
HFFibreFiducial.h
AlCaHLTBitMon_ParallelJobs.p
p
Definition:
AlCaHLTBitMon_ParallelJobs.py:153
pos
Definition:
PixelAliasList.h:18
FSimTrack.h
FastHFShowerLibrary::numberingScheme
HcalNumberingScheme numberingScheme
Definition:
FastHFShowerLibrary.h:57
edm::LogInfo
Log< level::Info, false > LogInfo
Definition:
MessageLogger.h:125
createJobs.tmp
tmp
align.sh
Definition:
createJobs.py:716
convertSQLiteXML.ok
bool ok
Definition:
convertSQLiteXML.py:98
RandomEngineAndDistribution.h
HFFibreFiducial::PMTNumber
int PMTNumber(const G4ThreeVector &pe_effect)
Definition:
HFFibreFiducial.cc:9
FastHFShowerLibrary::name
std::string name
Definition:
FastHFShowerLibrary.h:62
HcalSimulationConstants::hcalsimpar
const HcalSimulationParameters * hcalsimpar() const
Definition:
HcalSimulationConstants.h:20
FSimEvent.h
FastHFShowerLibrary::FastHFShowerLibrary
FastHFShowerLibrary(edm::ParameterSet const &p)
Definition:
FastHFShowerLibrary.cc:44
HcalParametersRcd.h
edm::EventSetup::get
T get() const
Definition:
EventSetup.h:80
FSimTrack::vfcalEntrance
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition:
FSimTrack.h:149
HcalSimulationConstants
Definition:
HcalSimulationConstants.h:15
edm::ESHandle
Definition:
DTSurvey.h:22
FSimTrack::onVFcal
int onVFcal() const
Definition:
FSimTrack.h:121
HcalDDDSimConstants::maxHFDepth
int maxHFDepth(const int &ieta, const int &iphi) const
Definition:
HcalDDDSimConstants.cc:655
LEDCalibrationChannels.depth
depth
Definition:
LEDCalibrationChannels.py:65
RawParticle::Vect
XYZVector Vect() const
the momentum threevector
Definition:
RawParticle.h:323
FastHFShowerLibrary::hitMap
std::map< CaloHitID, float > hitMap
Definition:
FastHFShowerLibrary.h:59
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
math::XYZVectorD
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition:
Vector3D.h:8
FSimTrack::id
int id() const
the index in FBaseSimEvent and other vectors
Definition:
FSimTrack.h:96
LEDCalibrationChannels.ieta
ieta
Definition:
LEDCalibrationChannels.py:63
bphysicsOniaDQM_cfi.vertex
vertex
Definition:
bphysicsOniaDQM_cfi.py:7
LogDebug
#define LogDebug(id)
Definition:
MessageLogger.h:223
edm::ParameterSet
Definition:
ParameterSet.h:47
HcalDetId.h
IdealGeometryRecord.h
CoreSimTrack::type
int type() const
particle type (HEP PDT convension)
Definition:
CoreSimTrack.h:22
edm::EventSetup
Definition:
EventSetup.h:57
HcalSubdetector.h
get
#define get
RandomEngineAndDistribution::theEngine
CLHEP::HepRandomEngine & theEngine() const
Definition:
RandomEngineAndDistribution.h:25
HcalForward
Definition:
HcalAssistant.h:36
FastHFShowerLibrary::modifyDepth
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition:
FastHFShowerLibrary.cc:138
RawParticle::e
double e() const
energy of the momentum
Definition:
RawParticle.h:305
initializeOnce
static std::once_flag initializeOnce
Definition:
FastHFShowerLibrary.cc:42
CaloHitID
Definition:
CaloHitID.h:11
HcalSimulationConstants.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
Exception.h
FastHFShowerLibrary::numberingFromDDD
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition:
FastHFShowerLibrary.h:55
FastHFShowerLibrary::hfshower
std::unique_ptr< HFShowerLibrary > hfshower
Definition:
FastHFShowerLibrary.h:54
FSimTrack
Definition:
FSimTrack.h:30
FastHFShowerLibrary::SetRandom
void SetRandom(const RandomEngineAndDistribution *)
Definition:
FastHFShowerLibrary.cc:76
ntuplemaker.time
time
Definition:
ntuplemaker.py:310
HcalSimNumberingRecord
Definition:
HcalSimNumberingRecord.h:25
FastHFShowerLibrary::applyFidCut
bool applyFidCut
Definition:
FastHFShowerLibrary.h:61
weight
Definition:
weight.py:1
FastHFShowerLibrary::initHFShowerLibrary
const void initHFShowerLibrary(const edm::EventSetup &)
Definition:
FastHFShowerLibrary.cc:49
HcalDDDSimConstants.h
RandomEngineAndDistribution
Definition:
RandomEngineAndDistribution.h:18
Generated for CMSSW Reference Manual by
1.8.16