Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Groups
Pages
SimG4CMS
Calo
interface
HFShowerParam.h
Go to the documentation of this file.
1
#ifndef SimG4CMS_HFShowerParam_h
2
#define SimG4CMS_HFShowerParam_h
3
// File: HFShowerParam.h
5
// Description: Generates hits for HF with some parametrized information
7
8
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
9
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
10
#include "
Geometry/HcalCommonData/interface/HcalDDDSimConstants.h
"
11
#include "
CondFormats/GeometryObjects/interface/HcalSimulationParameters.h
"
12
#include "
SimG4CMS/Calo/interface/HFShowerLibrary.h
"
13
#include "
SimG4CMS/Calo/interface/HFFibre.h
"
14
#include "
SimG4CMS/Calo/interface/HFGflash.h
"
15
16
#include "G4ThreeVector.hh"
17
18
class
G4Step;
19
20
#include <TH1F.h>
21
#include <TH2F.h>
22
#include <string>
23
#include <vector>
24
25
class
HFShowerParam
{
26
public
:
27
HFShowerParam
(
const
std::string
&
name
,
28
const
HcalDDDSimConstants
* hcons,
29
const
HcalSimulationParameters
* hps,
30
edm::ParameterSet
const
&
p
);
31
virtual
~HFShowerParam
();
32
33
public
:
34
struct
Hit
{
35
Hit
() {}
36
G4ThreeVector
position
;
37
int
depth
;
38
double
time
;
39
double
edep
;
40
};
41
std::vector<Hit>
getHits
(
const
G4Step* aStep,
double
weight
,
bool
& isKilled);
42
43
private
:
44
const
HcalDDDSimConstants
*
hcalConstants_
;
45
std::unique_ptr<HFShowerLibrary>
showerLibrary_
;
46
std::unique_ptr<HFFibre>
fibre_
;
47
std::unique_ptr<HFGflash>
gflash_
;
48
bool
fillHisto_
;
49
double
pePerGeV_
,
edMin_
,
ref_index_
,
aperture_
,
attLMeanInv_
;
50
bool
trackEM_
,
equalizeTimeShift_
,
onlyLong_
,
applyFidCut_
,
parametrizeLast_
;
51
G4int
emPDG_
,
epPDG_
,
gammaPDG_
;
52
std::vector<double>
gpar_
;
53
TH1F *
em_long_1_
, *
em_lateral_1_
, *
em_long_2_
, *
em_lateral_2_
;
54
TH1F *
hzvem_
, *
hzvhad_
, *
em_long_1_tuned_
, *
em_long_gflash_
;
55
TH1F*
em_long_sl_
;
56
TH2F *
em_2d_1_
, *
em_2d_2_
;
57
};
58
59
#endif // HFShowerParam_h
HFShowerParam::gflash_
std::unique_ptr< HFGflash > gflash_
Definition:
HFShowerParam.h:47
HFShowerParam::edMin_
double edMin_
Definition:
HFShowerParam.h:49
HFShowerParam::HFShowerParam
HFShowerParam(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
Definition:
HFShowerParam.cc:29
HFShowerParam::aperture_
double aperture_
Definition:
HFShowerParam.h:49
HcalDDDSimConstants.h
HFShowerParam::em_long_sl_
TH1F * em_long_sl_
Definition:
HFShowerParam.h:55
HFShowerParam::gammaPDG_
G4int gammaPDG_
Definition:
HFShowerParam.h:51
MessageLogger.h
HFShowerParam::attLMeanInv_
double attLMeanInv_
Definition:
HFShowerParam.h:49
HFShowerParam::~HFShowerParam
virtual ~HFShowerParam()
Definition:
HFShowerParam.cc:106
HFShowerParam::pePerGeV_
double pePerGeV_
Definition:
HFShowerParam.h:49
HFShowerParam
Definition:
HFShowerParam.h:25
HFShowerParam::em_long_2_
TH1F * em_long_2_
Definition:
HFShowerParam.h:53
HFShowerParam::em_long_1_tuned_
TH1F * em_long_1_tuned_
Definition:
HFShowerParam.h:54
HcalDDDSimConstants
Definition:
HcalDDDSimConstants.h:24
HFShowerParam::epPDG_
G4int epPDG_
Definition:
HFShowerParam.h:51
HFShowerParam::em_long_1_
TH1F * em_long_1_
Definition:
HFShowerParam.h:53
HFShowerParam::fillHisto_
bool fillHisto_
Definition:
HFShowerParam.h:48
mergeVDriftHistosByStation.name
string name
Definition:
mergeVDriftHistosByStation.py:78
HFShowerParam::em_long_gflash_
TH1F * em_long_gflash_
Definition:
HFShowerParam.h:54
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
HFGflash.h
ParameterSet.h
HFShowerParam::onlyLong_
bool onlyLong_
Definition:
HFShowerParam.h:50
HFShowerParam::em_lateral_2_
TH1F * em_lateral_2_
Definition:
HFShowerParam.h:53
HcalSimulationParameters.h
HFShowerParam::Hit::edep
double edep
Definition:
HFShowerParam.h:39
HFShowerParam::gpar_
std::vector< double > gpar_
Definition:
HFShowerParam.h:52
HFShowerLibrary.h
HFShowerParam::fibre_
std::unique_ptr< HFFibre > fibre_
Definition:
HFShowerParam.h:46
HFShowerParam::getHits
std::vector< Hit > getHits(const G4Step *aStep, double weight, bool &isKilled)
Definition:
HFShowerParam.cc:108
HFShowerParam::em_lateral_1_
TH1F * em_lateral_1_
Definition:
HFShowerParam.h:53
HcalSimulationParameters
Definition:
HcalSimulationParameters.h:6
HFShowerParam::Hit::position
G4ThreeVector position
Definition:
HFShowerParam.h:36
HFFibre.h
HFShowerParam::emPDG_
G4int emPDG_
Definition:
HFShowerParam.h:51
HFShowerParam::equalizeTimeShift_
bool equalizeTimeShift_
Definition:
HFShowerParam.h:50
HFShowerParam::em_2d_2_
TH2F * em_2d_2_
Definition:
HFShowerParam.h:56
HFShowerParam::em_2d_1_
TH2F * em_2d_1_
Definition:
HFShowerParam.h:56
HFShowerParam::trackEM_
bool trackEM_
Definition:
HFShowerParam.h:50
HFShowerParam::parametrizeLast_
bool parametrizeLast_
Definition:
HFShowerParam.h:50
HFShowerParam::showerLibrary_
std::unique_ptr< HFShowerLibrary > showerLibrary_
Definition:
HFShowerParam.h:45
HFShowerParam::Hit::time
double time
Definition:
HFShowerParam.h:38
AlCaHLTBitMon_ParallelJobs.p
tuple p
Definition:
AlCaHLTBitMon_ParallelJobs.py:153
HFShowerParam::ref_index_
double ref_index_
Definition:
HFShowerParam.h:49
HFShowerParam::applyFidCut_
bool applyFidCut_
Definition:
HFShowerParam.h:50
HFShowerParam::Hit::depth
int depth
Definition:
HFShowerParam.h:37
HFShowerParam::hcalConstants_
const HcalDDDSimConstants * hcalConstants_
Definition:
HFShowerParam.h:44
edm::ParameterSet
Definition:
ParameterSet.h:47
HFShowerParam::Hit
Definition:
HFShowerParam.h:34
histoStyle.weight
int weight
Definition:
histoStyle.py:51
HFShowerParam::Hit::Hit
Hit()
Definition:
HFShowerParam.h:35
HFShowerParam::hzvhad_
TH1F * hzvhad_
Definition:
HFShowerParam.h:54
HFShowerParam::hzvem_
TH1F * hzvem_
Definition:
HFShowerParam.h:54
Generated for CMSSW Reference Manual by
1.8.5