CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 (uint32_t &id)
 
void recoHFShowerLibrary (const FSimTrack &myTrack)
 
 ~FastHFShowerLibrary ()
 

Private Attributes

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

Detailed Description

Definition at line 40 of file FastHFShowerLibrary.h.

Constructor & Destructor Documentation

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

Definition at line 39 of file FastHFShowerLibrary.cc.

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

40  : fast(p) {
41  edm::ParameterSet m_HS = p.getParameter<edm::ParameterSet>("HFShowerLibrary");
42  applyFidCut = m_HS.getParameter<bool>("ApplyFiducialCut");
43 }
T getParameter(std::string const &) const
const edm::ParameterSet fast
FastHFShowerLibrary::~FastHFShowerLibrary ( )
inline

Definition at line 46 of file FastHFShowerLibrary.h.

46 {;}

Member Function Documentation

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

Definition at line 53 of file FastHFShowerLibrary.h.

References hitMap.

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

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

Definition at line 45 of file FastHFShowerLibrary.cc.

References fast, edm::EventSetup::get(), hcalConstants, hfshower, name, numberingFromDDD, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by FamosManager::setupGeometryAndField().

45  {
46 
47  edm::LogInfo("FastCalorimetry") << "initHFShowerLibrary::initialization";
48 
50  iSetup.get<IdealGeometryRecord>().get(cpv);
51 
53  iSetup.get<HcalSimNumberingRecord>().get(hdc);
54  hcalConstants = (HcalDDDSimConstants*)(&(*hdc));
55 
56  std::string name = "HcalHits";
58  hfshower.reset(new HFShowerLibrary(name,*cpv,fast));
59 
60  // Geant4 particles
61  G4DecayPhysics decays;
62  decays.ConstructParticle();
63  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
64  partTable->SetReadiness();
65 
66  hfshower->initRun(partTable, hcalConstants); // init particle code
67 }
const edm::ParameterSet fast
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
const T & get() const
Definition: EventSetup.h:56
std::unique_ptr< HFShowerLibrary > hfshower
HcalDDDSimConstants * hcalConstants
void FastHFShowerLibrary::modifyDepth ( uint32_t &  id)

Definition at line 125 of file FastHFShowerLibrary.cc.

References HcalDetId::depth(), eta, hcalConstants, HcalForward, HcalDetId::ieta(), HcalDetId::iphi(), HcalDDDSimConstants::maxHFDepth(), phi, DetId::rawId(), and HcalDetId::subdet().

Referenced by CalorimetryManager::loadFromHcal().

125  {
126 
127  HcalDetId hid(id);
128  if (hid.subdet() == HcalForward) {
129  int eta = hid.ieta();
130  int phi = hid.iphi();
131  if (hcalConstants->maxHFDepth(eta,phi) > 2) {
132  if (hid.depth() <= 2) {
133  int dep = (G4UniformRand() > 0.5) ? (2+hid.depth()) : hid.depth();
134  id = HcalDetId(HcalForward,eta,phi,dep).rawId();
135  }
136  }
137  }
138 }
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
HcalDDDSimConstants * hcalConstants
int maxHFDepth(int ieta, int iphi) const
void FastHFShowerLibrary::recoHFShowerLibrary ( const FSimTrack myTrack)

Definition at line 69 of file FastHFShowerLibrary.cc.

References applyFidCut, HLT_FULL_cff::depth, HcalNumberingScheme::getUnitID(), hfshower, hitMap, i, FSimTrack::id(), numberingFromDDD, numberingScheme, convertSQLiteXML::ok, FSimTrack::onVFcal(), HFFibreFiducial::PMTNumber(), tmp, CoreSimTrack::type(), RawParticle::vertex(), FSimTrack::vfcalEntrance(), puppiForMET_cff::weight, RawParticle::X(), RawParticle::Y(), and RawParticle::Z().

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

69  {
70 
71 #ifdef DebugLog
72  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: recoHFShowerLibrary ";
73 #endif
74 
75  if(!myTrack.onVFcal()) {
76 #ifdef DebugLog
77  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: we should not be here ";
78 #endif
79  }
80 
81  hitMap.clear();
82  double eGen = 1000.*myTrack.vfcalEntrance().e(); // energy in [MeV]
83  double delZv = (myTrack.vfcalEntrance().vertex().Z()>0.0) ? 50.0 : -50.0;
84  G4ThreeVector vertex( 10.*myTrack.vfcalEntrance().vertex().X(),
85  10.*myTrack.vfcalEntrance().vertex().Y(),
86  10.*myTrack.vfcalEntrance().vertex().Z()+delZv); // in [mm]
87 
88  G4ThreeVector direction(myTrack.vfcalEntrance().Vect().X(),
89  myTrack.vfcalEntrance().Vect().Y(),
90  myTrack.vfcalEntrance().Vect().Z());
91 
92  bool ok;
93  double weight = 1.0; // rad. damage
94  int parCode = myTrack.type();
95  double tSlice = 0.1*vertex.mag()/29.98;
96 
97  std::vector<HFShowerLibrary::Hit> hits =
98  hfshower->fillHits(vertex,direction,parCode,eGen,ok,weight,false,tSlice);
99 
100  for (unsigned int i=0; i<hits.size(); ++i) {
101  G4ThreeVector pos = hits[i].position;
102  int depth = hits[i].depth;
103  double time = hits[i].time;
104  if (!applyFidCut || (HFFibreFiducial::PMTNumber(pos)>0) ) {
105 // if (!applyFidCut || (applyFidCut && HFFibreFiducial::PMTNumber(pos)>0)) {
106  int det = 5;
107  int lay = 1;
108  uint32_t id = 0;
109  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, depth, lay);
110  id = numberingScheme.getUnitID(tmp);
111 
112  CaloHitID current_id(id,time,myTrack.id());
113  std::map<CaloHitID,float>::iterator cellitr;
114  cellitr = hitMap.find(current_id);
115  if(cellitr==hitMap.end()) {
116  hitMap.insert(std::pair<CaloHitID,float>(current_id,1.0));
117  } else {
118  cellitr->second += 1.0;
119  }
120  } // end of isItinFidVolume check
121  } // end loop over hits
122 
123 }
int i
Definition: DBlmapReader.cc:9
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:139
virtual uint32_t getUnitID(const HcalNumberingFromDDD::HcalID &id)
static int PMTNumber(const G4ThreeVector &pe_effect)
HcalNumberingScheme numberingScheme
double Y() const
y of vertex
Definition: RawParticle.h:275
double Z() const
z of vertex
Definition: RawParticle.h:276
int onVFcal() const
Definition: FSimTrack.h:111
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
double X() const
x of vertex
Definition: RawParticle.h:274
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
std::unique_ptr< HFShowerLibrary > hfshower
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:86
std::map< CaloHitID, float > hitMap

Member Data Documentation

bool FastHFShowerLibrary::applyFidCut
private

Definition at line 65 of file FastHFShowerLibrary.h.

Referenced by FastHFShowerLibrary(), and recoHFShowerLibrary().

const edm::ParameterSet FastHFShowerLibrary::fast
private

Definition at line 53 of file FastHFShowerLibrary.h.

Referenced by initHFShowerLibrary().

HcalDDDSimConstants* FastHFShowerLibrary::hcalConstants
private

Definition at line 60 of file FastHFShowerLibrary.h.

Referenced by initHFShowerLibrary(), and modifyDepth().

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

Definition at line 58 of file FastHFShowerLibrary.h.

Referenced by initHFShowerLibrary(), and recoHFShowerLibrary().

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

Definition at line 63 of file FastHFShowerLibrary.h.

Referenced by getHitsMap(), and recoHFShowerLibrary().

std::string FastHFShowerLibrary::name
private

Definition at line 66 of file FastHFShowerLibrary.h.

Referenced by ElectronMVAID.ElectronMVAID::__call__(), dirstructure.Directory::__create_pie_image(), dqm_interfaces.DirID::__eq__(), dirstructure.Directory::__get_full_path(), dirstructure.Comparison::__get_img_name(), dataset.Dataset::__getDataType(), dataset.Dataset::__getFileInfoList(), cuy.divideElement::__init__(), cuy.plotElement::__init__(), cuy.additionElement::__init__(), cuy.superimposeElement::__init__(), cuy.graphElement::__init__(), 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.CFG::__str__(), counter.Counter::__str__(), average.Average::__str__(), 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(), 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 59 of file FastHFShowerLibrary.h.

Referenced by initHFShowerLibrary(), and recoHFShowerLibrary().

HcalNumberingScheme FastHFShowerLibrary::numberingScheme
private

Definition at line 61 of file FastHFShowerLibrary.h.

Referenced by recoHFShowerLibrary().