CMS 3D CMS Logo

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

#include <GEMDigiModule.h>

Public Types

typedef edm::DetSet< GEMDigiSimLinkGEMDigiSimLinks
 
typedef edm::DetSet< StripDigiSimLinkStripDigiSimLinks
 

Public Member Functions

void fillDigis (int rollDetId, GEMDigiCollection &)
 
 GEMDigiModule (const edm::ParameterSet &)
 
const GEMDigiSimLinksgemDigiSimLinks () const
 
void setGeometry (const GEMGeometry *)
 
void simulate (const GEMEtaPartition *, const edm::PSimHitContainer &, CLHEP::HepRandomEngine *)
 
const StripDigiSimLinksstripDigiSimLinks () const
 
 ~GEMDigiModule ()
 

Private Member Functions

void addLinks (unsigned int strip, int bx)
 creates links from Digi to SimTrack More...
 
void addLinksWithPartId (unsigned int strip, int bx)
 

Private Attributes

DetectorHitMap detectorHitMap_
 
const GEMGeometrygeometry_
 
std::vector< std::unique_ptr< GEMDigiModel > > models
 
StripDigiSimLinks stripDigiSimLinks_
 
Strips strips_
 
GEMDigiSimLinks theGemDigiSimLinks_
 

Detailed Description

Base Class for the GEM strip response simulation

Author
Sven Dildick by Yechan Kang

Definition at line 35 of file GEMDigiModule.h.

Member Typedef Documentation

◆ GEMDigiSimLinks

Definition at line 42 of file GEMDigiModule.h.

◆ StripDigiSimLinks

Definition at line 41 of file GEMDigiModule.h.

Constructor & Destructor Documentation

◆ GEMDigiModule()

GEMDigiModule::GEMDigiModule ( const edm::ParameterSet config)

Definition at line 9 of file GEMDigiModule.cc.

9  {
10  bool simulateBkgNoise_(config.getParameter<bool>("simulateBkgNoise"));
11  bool simulateIntrinsicNoise_(config.getParameter<bool>("simulateIntrinsicNoise"));
12  if (simulateIntrinsicNoise_) {
13  models.push_back(std::make_unique<GEMNoiseModel>(config));
14  }
15  if (simulateBkgNoise_) {
16  models.push_back(std::make_unique<GEMBkgModel>(config));
17  }
18  models.push_back(std::make_unique<GEMSignalModel>(config));
19 }
Definition: config.py:1
Definition: models.py:1

◆ ~GEMDigiModule()

GEMDigiModule::~GEMDigiModule ( )
default

Member Function Documentation

◆ addLinks()

void GEMDigiModule::addLinks ( unsigned int  strip,
int  bx 
)
private

creates links from Digi to SimTrack

Definition at line 50 of file GEMDigiModule.cc.

References nano_mu_digi_cff::bx, ALCARECOTkAlJpsiMuMu_cff::charge, detectorHitMap_, f, edm::DetSet< T >::push_back(), nano_mu_digi_cff::strip, and stripDigiSimLinks_.

Referenced by fillDigis().

50  {
51  std::pair<unsigned int, int> digi(strip, bx);
52  auto channelHitItr = detectorHitMap_.equal_range(digi);
53 
54  // find the fraction contribution for each SimTrack
55  std::map<int, float> simTrackChargeMap;
56  std::map<int, EncodedEventId> eventIdMap;
57  float totalCharge(0.);
58  for (auto hitItr = channelHitItr.first; hitItr != channelHitItr.second; ++hitItr) {
59  const PSimHit* hit(hitItr->second);
60  // might be zero for unit tests and such
61  if (hit == nullptr)
62  continue;
63 
64  int simTrackId(hit->trackId());
65  //float charge = hit->getCharge();
66  const float charge(1.f);
67  auto chargeItr = simTrackChargeMap.find(simTrackId);
68  if (chargeItr == simTrackChargeMap.end()) {
69  simTrackChargeMap[simTrackId] = charge;
70  eventIdMap[simTrackId] = hit->eventId();
71  } else {
72  chargeItr->second += charge;
73  }
74  totalCharge += charge;
75  }
76 
77  for (const auto& charge : simTrackChargeMap) {
78  const int simTrackId(charge.first);
79  auto link(StripDigiSimLink(strip, simTrackId, eventIdMap[simTrackId], charge.second / totalCharge));
81  }
82 }
void push_back(const T &t)
Definition: DetSet.h:66
StripDigiSimLinks stripDigiSimLinks_
Definition: GEMDigiModule.h:64
double f[11][100]
DetectorHitMap detectorHitMap_
Definition: GEMDigiModule.h:63

◆ addLinksWithPartId()

void GEMDigiModule::addLinksWithPartId ( unsigned int  strip,
int  bx 
)
private

Definition at line 84 of file GEMDigiModule.cc.

References nano_mu_digi_cff::bx, detectorHitMap_, edm::DetSet< T >::push_back(), nano_mu_digi_cff::strip, and theGemDigiSimLinks_.

Referenced by fillDigis().

84  {
85  std::pair<unsigned int, int> digi(strip, bx);
86  std::pair<DetectorHitMap::iterator, DetectorHitMap::iterator> channelHitItr = detectorHitMap_.equal_range(digi);
87 
88  for (DetectorHitMap::iterator hitItr = channelHitItr.first; hitItr != channelHitItr.second; ++hitItr) {
89  const PSimHit* hit = (hitItr->second);
90  // might be zero for unit tests and such
91  if (hit == nullptr)
92  continue;
93 
94  theGemDigiSimLinks_.push_back(GEMDigiSimLink(strip, bx, hit->particleType(), hit->trackId(), hit->eventId()));
95  }
96 }
void push_back(const T &t)
Definition: DetSet.h:66
GEMDigiSimLinks theGemDigiSimLinks_
Definition: GEMDigiModule.h:65
DetectorHitMap detectorHitMap_
Definition: GEMDigiModule.h:63

◆ fillDigis()

void GEMDigiModule::fillDigis ( int  rollDetId,
GEMDigiCollection digis 
)

Definition at line 37 of file GEMDigiModule.cc.

References addLinks(), addLinksWithPartId(), ztail::d, and strips_.

37  {
38  for (const auto& d : strips_) {
39  if (d.second == -999)
40  continue;
41  // (strip, bx)
42  GEMDigi digi(d.first, d.second);
43  digis.insertDigi(GEMDetId(rollDetId), digi);
44  addLinks(d.first, d.second);
45  addLinksWithPartId(d.first, d.second);
46  }
47  strips_.clear();
48 }
void addLinks(unsigned int strip, int bx)
creates links from Digi to SimTrack
void addLinksWithPartId(unsigned int strip, int bx)
d
Definition: ztail.py:151

◆ gemDigiSimLinks()

const GEMDigiSimLinks& GEMDigiModule::gemDigiSimLinks ( ) const
inline

Definition at line 51 of file GEMDigiModule.h.

References theGemDigiSimLinks_.

51 { return theGemDigiSimLinks_; }
GEMDigiSimLinks theGemDigiSimLinks_
Definition: GEMDigiModule.h:65

◆ setGeometry()

void GEMDigiModule::setGeometry ( const GEMGeometry geom)

Definition at line 98 of file GEMDigiModule.cc.

References relativeConstraints::geom, and isotrackApplyRegressor::model.

98  {
99  for (auto&& model : models) {
100  model->setGeometry(geom);
101  }
102 }
Definition: models.py:1

◆ simulate()

void GEMDigiModule::simulate ( const GEMEtaPartition roll,
const edm::PSimHitContainer simHits,
CLHEP::HepRandomEngine *  engine 
)

Definition at line 23 of file GEMDigiModule.cc.

References edm::DetSet< T >::clear(), detectorHitMap_, isotrackApplyRegressor::model, nano_mu_digi_cff::roll, FastTrackerRecHitCombiner_cfi::simHits, stripDigiSimLinks_, strips_, and theGemDigiSimLinks_.

25  {
27  detectorHitMap_.clear();
28  stripDigiSimLinks_ = StripDigiSimLinks(roll->id().rawId());
30  theGemDigiSimLinks_ = GEMDigiSimLinks(roll->id().rawId());
31  for (auto&& model : models) {
32  model->simulate(roll, simHits, engine, strips_, detectorHitMap_);
33  }
34  return;
35 }
edm::DetSet< StripDigiSimLink > StripDigiSimLinks
Definition: GEMDigiModule.h:41
StripDigiSimLinks stripDigiSimLinks_
Definition: GEMDigiModule.h:64
Definition: models.py:1
GEMDigiSimLinks theGemDigiSimLinks_
Definition: GEMDigiModule.h:65
edm::DetSet< GEMDigiSimLink > GEMDigiSimLinks
Definition: GEMDigiModule.h:42
void clear()
Definition: DetSet.h:71
DetectorHitMap detectorHitMap_
Definition: GEMDigiModule.h:63

◆ stripDigiSimLinks()

const StripDigiSimLinks& GEMDigiModule::stripDigiSimLinks ( ) const
inline

Definition at line 50 of file GEMDigiModule.h.

References stripDigiSimLinks_.

50 { return stripDigiSimLinks_; }
StripDigiSimLinks stripDigiSimLinks_
Definition: GEMDigiModule.h:64

Member Data Documentation

◆ detectorHitMap_

DetectorHitMap GEMDigiModule::detectorHitMap_
private

Definition at line 63 of file GEMDigiModule.h.

Referenced by addLinks(), addLinksWithPartId(), and simulate().

◆ geometry_

const GEMGeometry* GEMDigiModule::geometry_
private

Definition at line 54 of file GEMDigiModule.h.

◆ models

std::vector<std::unique_ptr<GEMDigiModel> > GEMDigiModule::models
private

Definition at line 56 of file GEMDigiModule.h.

◆ stripDigiSimLinks_

StripDigiSimLinks GEMDigiModule::stripDigiSimLinks_
private

Definition at line 64 of file GEMDigiModule.h.

Referenced by addLinks(), simulate(), and stripDigiSimLinks().

◆ strips_

Strips GEMDigiModule::strips_
private

Definition at line 62 of file GEMDigiModule.h.

Referenced by fillDigis(), and simulate().

◆ theGemDigiSimLinks_

GEMDigiSimLinks GEMDigiModule::theGemDigiSimLinks_
private

Definition at line 65 of file GEMDigiModule.h.

Referenced by addLinksWithPartId(), gemDigiSimLinks(), and simulate().