CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
FastHFShowerLibrary Class Reference

#include <FastHFShowerLibrary.h>

Public Member Functions

 FastHFShowerLibrary (edm::ParameterSet const &p)
 
const std::map< CaloHitID, float > & getHitsMap ()
 
void const initHFShowerLibrary (const edm::EventSetup &)
 
void modifyDepth (HcalNumberingFromDDD::HcalID &id)
 
void recoHFShowerLibrary (const FSimTrack &myTrack)
 
void SetRandom (const RandomEngineAndDistribution *)
 
 ~FastHFShowerLibrary ()
 

Private Attributes

bool applyFidCut
 
const edm::ParameterSet fast
 
const HcalDDDSimConstantshcalConstants
 
std::unique_ptr< HFShowerLibraryhfshower
 
std::map< CaloHitID, float > hitMap
 
std::string name
 
std::unique_ptr< HcalNumberingFromDDDnumberingFromDDD
 
HcalNumberingScheme numberingScheme
 

Detailed Description

Definition at line 40 of file FastHFShowerLibrary.h.

Constructor & Destructor Documentation

FastHFShowerLibrary::FastHFShowerLibrary ( edm::ParameterSet const &  p)

Definition at line 43 of file FastHFShowerLibrary.cc.

References applyFidCut, and edm::ParameterSet::getParameter().

44  : fast(p) {
45  edm::ParameterSet m_HS = p.getParameter<edm::ParameterSet>("HFShowerLibrary");
46  applyFidCut = m_HS.getParameter<bool>("ApplyFiducialCut");
47 }
T getParameter(std::string const &) const
const edm::ParameterSet fast
FastHFShowerLibrary::~FastHFShowerLibrary ( )
inline

Definition at line 46 of file FastHFShowerLibrary.h.

References initHFShowerLibrary(), modifyDepth(), and recoHFShowerLibrary().

46 {;}

Member Function Documentation

const std::map<CaloHitID,float>& FastHFShowerLibrary::getHitsMap ( )
inline

Definition at line 53 of file FastHFShowerLibrary.h.

References hitMap, and SetRandom().

Referenced by CalorimetryManager::HDShowerSimulation(), and CalorimetryManager::reconstructTrack().

53 { return hitMap; };
std::map< CaloHitID, float > hitMap
void const FastHFShowerLibrary::initHFShowerLibrary ( const edm::EventSetup iSetup)

Definition at line 49 of file FastHFShowerLibrary.cc.

References fast, edm::EventSetup::get(), hcalConstants, hfshower, g4SimHits_cfi::HFShowerLibrary, initializeOnce, name, numberingFromDDD, edm::ESHandle< T >::product(), and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by FamosManager::setupGeometryAndField(), and ~FastHFShowerLibrary().

49  {
50 
51  edm::LogInfo("FastCalorimetry") << "initHFShowerLibrary::initialization";
52 
54  iSetup.get<IdealGeometryRecord>().get(cpv);
55 
57  iSetup.get<HcalSimNumberingRecord>().get(hdc);
58  hcalConstants = hdc.product();
59 
60  std::string name = "HcalHits";
62  hfshower.reset(new HFShowerLibrary(name,*cpv,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 }
const edm::ParameterSet fast
static std::once_flag initializeOnce
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
std::unique_ptr< HFShowerLibrary > hfshower
T get() const
Definition: EventSetup.h:71
const HcalDDDSimConstants * hcalConstants
T const * product() const
Definition: ESHandle.h:86
void FastHFShowerLibrary::modifyDepth ( HcalNumberingFromDDD::HcalID id)

Definition at line 144 of file FastHFShowerLibrary.cc.

References egammaForCoreTracking_cff::depth, hcalConstants, HcalForward, and HcalDDDSimConstants::maxHFDepth().

Referenced by recoHFShowerLibrary(), and ~FastHFShowerLibrary().

144  {
145  if (id.subdet == HcalForward) {
146  int ieta = (id.zside == 0) ? -id.etaR : id.etaR;
147  if (hcalConstants->maxHFDepth(ieta,id.phis) > 2) {
148  if (id.depth <= 2) {
149  if (G4UniformRand() > 0.5) id.depth += 2;
150  }
151  }
152  }
153 }
const HcalDDDSimConstants * hcalConstants
int maxHFDepth(const int &ieta, const int &iphi) const
void FastHFShowerLibrary::recoHFShowerLibrary ( const FSimTrack myTrack)

Definition at line 85 of file FastHFShowerLibrary.cc.

References applyFidCut, egammaForCoreTracking_cff::depth, RawParticle::e(), HcalNumberingScheme::getUnitID(), hfshower, hitMap, hfClusterShapes_cfi::hits, mps_fire::i, FSimTrack::id(), modifyDepth(), numberingFromDDD, numberingScheme, convertSQLiteXML::ok, FSimTrack::onVFcal(), HFFibreFiducial::PMTNumber(), ntuplemaker::time, tmp, CoreSimTrack::type(), RawParticle::Vect(), RawParticle::vertex(), and FSimTrack::vfcalEntrance().

Referenced by CalorimetryManager::HDShowerSimulation(), CalorimetryManager::reconstructTrack(), and ~FastHFShowerLibrary().

85  {
86 
87 #ifdef DebugLog
88  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: recoHFShowerLibrary ";
89 #endif
90 
91  if(!myTrack.onVFcal()) {
92 #ifdef DebugLog
93  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: we should not be here ";
94 #endif
95  }
96 
97  hitMap.clear();
98  double eGen = 1000.*myTrack.vfcalEntrance().e(); // energy in [MeV]
99  double delZv = (myTrack.vfcalEntrance().vertex().Z()>0.0) ? 50.0 : -50.0;
100  G4ThreeVector vertex( 10.*myTrack.vfcalEntrance().vertex().X(),
101  10.*myTrack.vfcalEntrance().vertex().Y(),
102  10.*myTrack.vfcalEntrance().vertex().Z()+delZv); // in [mm]
103 
104  G4ThreeVector direction(myTrack.vfcalEntrance().Vect().X(),
105  myTrack.vfcalEntrance().Vect().Y(),
106  myTrack.vfcalEntrance().Vect().Z());
107 
108  bool ok;
109  double weight = 1.0; // rad. damage
110  int parCode = myTrack.type();
111  double tSlice = 0.1*vertex.mag()/29.98;
112 
113  std::vector<HFShowerLibrary::Hit> hits =
114  hfshower->fillHits(vertex,direction,parCode,eGen,ok,weight,tSlice,false);
115 
116  for (unsigned int i=0; i<hits.size(); ++i) {
117  G4ThreeVector pos = hits[i].position;
118  int depth = hits[i].depth;
119  double time = hits[i].time;
120  if (!applyFidCut || (HFFibreFiducial::PMTNumber(pos)>0) ) {
121 // if (!applyFidCut || (applyFidCut && HFFibreFiducial::PMTNumber(pos)>0)) {
122  int det = 5;
123  int lay = 1;
124  uint32_t id = 0;
126  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(),pos.y(),
127  pos.z()), depth, lay);
128  modifyDepth(tmp);
129  id = numberingScheme.getUnitID(tmp);
130 
131  CaloHitID current_id(id,time,myTrack.id());
132  std::map<CaloHitID,float>::iterator cellitr;
133  cellitr = hitMap.find(current_id);
134  if(cellitr==hitMap.end()) {
135  hitMap.insert(std::pair<CaloHitID,float>(current_id,1.0));
136  } else {
137  cellitr->second += 1.0;
138  }
139  } // end of isItinFidVolume check
140  } // end loop over hits
141 
142 }
XYZVector Vect() const
the momentum threevector
Definition: RawParticle.h:342
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:144
virtual uint32_t getUnitID(const HcalNumberingFromDDD::HcalID &id)
Definition: weight.py:1
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
HcalNumberingScheme numberingScheme
double e() const
energy of the momentum
Definition: RawParticle.h:324
int onVFcal() const
Definition: FSimTrack.h:116
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:339
int PMTNumber(const G4ThreeVector &pe_effect)
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
std::unique_ptr< HFShowerLibrary > hfshower
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:91
std::map< CaloHitID, float > hitMap
void FastHFShowerLibrary::SetRandom ( const RandomEngineAndDistribution rnd)

Definition at line 76 of file FastHFShowerLibrary.cc.

References LogDebug, and RandomEngineAndDistribution::theEngine().

Referenced by getHitsMap(), CalorimetryManager::HDShowerSimulation(), and CalorimetryManager::initialize().

77 {
78  // define Geant4 engine per thread
79  G4Random::setTheEngine(&(rnd->theEngine()));
80  LogDebug("FastHFShowerLibrary::recoHFShowerLibrary")
81  << "Begin of event " << G4UniformRand() << " "
82  << rnd->theEngine().name() << " " << rnd->theEngine();
83 }
#define LogDebug(id)
CLHEP::HepRandomEngine & theEngine() const

Member Data Documentation

bool FastHFShowerLibrary::applyFidCut
private

Definition at line 67 of file FastHFShowerLibrary.h.

Referenced by FastHFShowerLibrary(), and recoHFShowerLibrary().

const edm::ParameterSet FastHFShowerLibrary::fast
private

Definition at line 59 of file FastHFShowerLibrary.h.

Referenced by initHFShowerLibrary().

const HcalDDDSimConstants* FastHFShowerLibrary::hcalConstants
private

Definition at line 62 of file FastHFShowerLibrary.h.

Referenced by initHFShowerLibrary(), and modifyDepth().

std::unique_ptr<HFShowerLibrary> FastHFShowerLibrary::hfshower
private

Definition at line 60 of file FastHFShowerLibrary.h.

Referenced by initHFShowerLibrary(), and recoHFShowerLibrary().

std::map<CaloHitID,float> FastHFShowerLibrary::hitMap
private

Definition at line 65 of file FastHFShowerLibrary.h.

Referenced by getHitsMap(), and recoHFShowerLibrary().

std::string FastHFShowerLibrary::name
private

Definition at line 68 of file FastHFShowerLibrary.h.

Referenced by ElectronMVAID.ElectronMVAID::__call__(), FWLite.ElectronMVAID::__call__(), dirstructure.Directory::__create_pie_image(), DisplayManager.DisplayManager::__del__(), dqm_interfaces.DirID::__eq__(), dirstructure.Directory::__get_full_path(), dirstructure.Comparison::__get_img_name(), dataset.Dataset::__getDataType(), dataset.Dataset::__getFileInfoList(), dirstructure.Comparison::__make_image(), core.autovars.NTupleVariable::__repr__(), core.autovars.NTupleObjectType::__repr__(), core.autovars.NTupleObject::__repr__(), core.autovars.NTupleCollection::__repr__(), dirstructure.Directory::__repr__(), dqm_interfaces.DirID::__repr__(), dirstructure.Comparison::__repr__(), config.Service::__setattr__(), config.CFG::__str__(), counter.Counter::__str__(), average.Average::__str__(), FWLite.WorkingPoints::_reformat_cut_definitions(), core.autovars.NTupleObjectType::addSubObjects(), core.autovars.NTupleObjectType::addVariables(), core.autovars.NTupleObjectType::allVars(), dirstructure.Directory::calcStats(), validation.Sample::digest(), python.rootplot.utilities.Hist::divide(), python.rootplot.utilities.Hist::divide_wilson(), DisplayManager.DisplayManager::Draw(), TreeCrawler.Package::dump(), core.autovars.NTupleVariable::fillBranch(), core.autovars.NTupleObject::fillBranches(), core.autovars.NTupleCollection::fillBranchesScalar(), core.autovars.NTupleCollection::fillBranchesVector(), core.autovars.NTupleCollection::get_cpp_declaration(), core.autovars.NTupleCollection::get_cpp_wrapper_class(), core.autovars.NTupleCollection::get_py_wrapper_class(), utils.StatisticalTest::get_status(), production_tasks.Task::getname(), dataset.CMSDataset::getPrimaryDatasetEntries(), dataset.PrivateDataset::getPrimaryDatasetEntries(), initHFShowerLibrary(), VIDSelectorBase.VIDSelectorBase::initialize(), core.autovars.NTupleVariable::makeBranch(), core.autovars.NTupleObject::makeBranches(), core.autovars.NTupleCollection::makeBranchesScalar(), core.autovars.NTupleCollection::makeBranchesVector(), dirstructure.Directory::print_report(), dataset.BaseDataset::printInfo(), dataset.Dataset::printInfo(), production_tasks.MonitorJobs::run(), python.rootplot.utilities.Hist::TGraph(), python.rootplot.utilities.Hist::TH1F(), Vispa.Views.PropertyView.Property::valueChanged(), counter.Counter::write(), and average.Average::write().

std::unique_ptr<HcalNumberingFromDDD> FastHFShowerLibrary::numberingFromDDD
private

Definition at line 61 of file FastHFShowerLibrary.h.

Referenced by initHFShowerLibrary(), and recoHFShowerLibrary().

HcalNumberingScheme FastHFShowerLibrary::numberingScheme
private

Definition at line 63 of file FastHFShowerLibrary.h.

Referenced by recoHFShowerLibrary().