CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalTBDigiProducer.cc
Go to the documentation of this file.
5 
15 
17  theParameterMap(new HcalTBSimParameterMap(ps)),
18  theHcalShape(new HcalShape()),
19  theHcalIntegratedShape(new CaloShapeIntegrator(theHcalShape)),
20  theHBHEResponse(new CaloHitResponse(theParameterMap, theHcalIntegratedShape)),
21  theHOResponse(new CaloHitResponse(theParameterMap, theHcalIntegratedShape)),
22  theAmplifier(0), theCoderFactory(0), theElectronicsSim(0),
23  theHitCorrection(0), theHBHEDigitizer(0), theHODigitizer(0), theHBHEHits(),
24  theHOHits(), thisPhaseShift(0) {
25  std::string const instance("simHcalDigis");
26  mixMod.produces<HBHEDigiCollection>(instance);
27  mixMod.produces<HODigiCollection>(instance);
28 
29  DetId detId(DetId::Hcal, 1);
30  bool syncPhase = (theParameterMap->simParameters(detId)).syncPhase();
31  doPhaseShift = !syncPhase;
32 
35 
36  bool doTimeSlew = ps.getParameter<bool>("doTimeSlew");
37  if(doTimeSlew) {
38  // no time slewing for HF
42  }
43 
44  bool doNoise = ps.getParameter<bool>("doNoise");
48 
51 
52  tunePhaseShift = ps.getUntrackedParameter<double>("tunePhaseShiftTB", 1.);
53  ecalTBInfoLabel = ps.getUntrackedParameter<std::string>("EcalTBInfoLabel","SimEcalTBG4Object");
54  edm::LogInfo("HcalSim") << "HcalTBDigiProducer initialized with doNoise = "
55  << doNoise << ", doTimeSlew = " << doTimeSlew
56  << " and doPhaseShift = " << doPhaseShift
57  << " tunePhasShift = " << tunePhaseShift;
58 
59 }
60 
62 
64  if (theHODigitizer) delete theHODigitizer;
65  if (theParameterMap) delete theParameterMap;
66  if (theHcalShape) delete theHcalShape;
68  if (theHBHEResponse) delete theHBHEResponse;
69  if (theHOResponse) delete theHOResponse;
71  if (theAmplifier) delete theAmplifier;
72  if (theCoderFactory) delete theCoderFactory;
74 }
75 
76 
78  // get the appropriate gains, noises, & widths for this event
80  eventSetup.get<HcalDbRecord>().get(conditions);
81  theAmplifier->setDbService(conditions.product());
82  theCoderFactory->setDbService(conditions.product());
83 
84  // get the correct geometry
85  checkGeometry(eventSetup);
86 
87  theHBHEHits.clear();
88  theHOHits.clear();
89  if (doPhaseShift) {
90 
91  edm::Handle<PEcalTBInfo> theEcalTBInfo;
92  e.getByLabel(ecalTBInfoLabel,theEcalTBInfo);
93  thisPhaseShift = theEcalTBInfo->phaseShift();
94 
95  DetId detIdHB(DetId::Hcal, 1);
96  setPhaseShift(detIdHB);
97  DetId detIdHO(DetId::Hcal, 3);
98  setPhaseShift(detIdHO);
99  }
100 
103 }
104 
105 void HcalTBDigiProducer::accumulateCaloHits(edm::Handle<std::vector<PCaloHit> > const& hcalHandle, int bunchCrossing) {
106 
107  LogDebug("HcalSim") << "HcalTBDigiProducer::accumulate trying to get SimHit";
108 
109  if(hcalHandle.isValid()) {
110  std::vector<PCaloHit> hits = *hcalHandle.product();
111  if(theHitCorrection != 0) {
113  }
114  LogDebug("HcalSim") << "HcalTBDigiProducer::accumulate Hits corrected";
115  theHBHEDigitizer->add(hits, bunchCrossing);
116  theHODigitizer->add(hits, bunchCrossing);
117  }
118 }
119 
121  // Step A: Get Inputs, and accumulate digis
122 
123  edm::InputTag hcalTag("g4SimHits", "HcalHits");
125  e.getByLabel(hcalTag, hcalHandle);
126 
127  accumulateCaloHits(hcalHandle, 0);
128 }
129 
131  // Step A: Get Inputs, and accumulate digis
132 
133  edm::InputTag hcalTag("g4SimHits", "HcalHits");
135  e.getByLabel(hcalTag, hcalHandle);
136 
137  accumulateCaloHits(hcalHandle, e.bunchCrossing());
138 }
139 
141  // Step B: Create empty output
142  std::auto_ptr<HBHEDigiCollection> hbheResult(new HBHEDigiCollection());
143  std::auto_ptr<HODigiCollection> hoResult(new HODigiCollection());
144  LogDebug("HcalSim") << "HcalTBDigiProducer::produce Empty collection created";
145  // Step C: Invoke the algorithm, getting back outputs.
146  theHBHEDigitizer->run(*hbheResult);
147  edm::LogInfo("HcalSim") << "HcalTBDigiProducer: HBHE digis : " << hbheResult->size();
148  theHODigitizer->run(*hoResult);
149  edm::LogInfo("HcalSim") << "HcalTBDigiProducer: HO digis : " << hoResult->size();
150 
151  // Step D: Put outputs into event
152  std::string const instance("simHcalDigis");
153  e.put(hbheResult, instance);
154  e.put(hoResult, instance);
155 
156 }
157 
159 
160  for (edm::PCaloHitContainer::const_iterator hitItr = hits.begin();
161  hitItr != hits.end(); ++hitItr) {
162  HcalSubdetector subdet = HcalDetId(hitItr->id()).subdet();
163  if(subdet == HcalBarrel || subdet == HcalEndcap) {
164  theHBHEHits.push_back(*hitItr);
165  } else if(subdet == HcalOuter) {
166  theHOHits.push_back(*hitItr);
167  } else {
168  edm::LogError("HcalSim") << "Bad HcalHit subdetector " << subdet;
169  }
170  }
171 }
172 
174 
175  // TODO find a way to avoid doing this every event
177  eventSetup.get<CaloGeometryRecord>().get(geometry);
178 
179  const CaloGeometry * pGeometry = &*geometry;
180 
181  // see if we need to update
182  if(pGeometry != theGeometry) {
183  theGeometry = pGeometry;
184  updateGeometry();
185  }
186 }
187 
189 
192 
193  // Get cells for HB and HE
194  hbheCells.clear();
196  std::vector<DetId> heCells = theGeometry->getValidDetIds(DetId::Hcal, HcalEndcap);
197  // combine HB & HE
198  hbheCells.insert(hbheCells.end(), heCells.begin(), heCells.end());
199 
200  // Get cells for HO
201  hoCells.clear();
203 
204  edm::LogInfo("HcalSim") << "HcalTBDigiProducer update Geometry with "
205  << hbheCells.size() << " cells in HB/HE and "
206  << hoCells.size() << " cells in HO";
207 
209  LogDebug("HcalSim") << "HcalTBDigiProducer: Set DetID's for HB/HE";
211  LogDebug("HcalSim") << "HcalTBDigiProducer: Set DetID's for HO";
212 }
213 
215 
217  if ( !parameters.syncPhase() ) {
218  int myDet = detId.subdetId();
219  double passPhaseShift = thisPhaseShift + tunePhaseShift;
220  if (myDet <= 2) {
221  theHBHEResponse->setPhaseShift(passPhaseShift);
222  } else {
223  theHOResponse->setPhaseShift(passPhaseShift);
224  }
225  }
226 }
#define LogDebug(id)
void setPhaseShift(const DetId &detId)
T getParameter(std::string const &) const
void setGeometry(const CaloGeometry *geometry)
geometry needed for time-of-flight
void setDbService(const HcalDbService *service)
the Producer will probably update this every event
T getUntrackedParameter(std::string const &, T const &) const
virtual void finalizeEvent(edm::Event &e, edm::EventSetup const &c) override
dictionary parameters
Definition: Parameters.py:2
void sortHits(const edm::PCaloHitContainer &hits)
fills the vectors for each subdetector
std::vector< PCaloHit > PCaloHitContainer
std::vector< DetId > hoCells
virtual const CaloSimParameters & simParameters(const DetId &id) const
CaloHitResponse * theHBHEResponse
void initializeHits()
static PFTauRenderPlugin instance
HBHEHitFilter theHBHEHitFilter
HODigitizer * theHODigitizer
void add(const std::vector< PCaloHit > &hits, int bunchCrossing)
HcalElectronicsSim * theElectronicsSim
HOHitFilter theHOHitFilter
HcalTBSimParameterMap * theParameterMap
Main class for Parameters in different subdetectors.
const CaloGeometry * theGeometry
void setPhaseShift(const double &thePhaseShift)
setting the phase shift for asynchronous trigger (e.g. test beams)
CaloTDigitizer< HBHEDigitizerTraits > HBHEDigitizer
std::vector< PCaloHit > theHBHEHits
CaloTDigitizer< HODigitizerTraits > HODigitizer
edm::SortedCollection< HODataFrame > HODigiCollection
CaloVShape * theHcalIntegratedShape
void setHitFilter(const CaloVHitFilter *filter)
if you want to reject hits, for example, from a certain subdetector, set this
Creates electronics signals from hits.
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
HBHEDigitizer * theHBHEDigitizer
void fillChargeSums(MixCollection< PCaloHit > &hits)
HcalSubdetector
Definition: HcalAssistant.h:32
void accumulateCaloHits(edm::Handle< std::vector< PCaloHit > > const &hits, int bunchCrossing)
void checkGeometry(const edm::EventSetup &eventSetup)
shaper for Hcal (not for HF)
Definition: HcalShape.h:15
std::vector< DetId > hbheCells
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
virtual void accumulate(edm::Event const &e, edm::EventSetup const &c) override
std::vector< PCaloHit > theHOHits
Definition: DetId.h:20
CaloHitResponse * theHOResponse
void run(MixCollection< PCaloHit > &, DigiCollection &)
turns hits into digis
HcalHitCorrection * theHitCorrection
void setDbService(const HcalDbService *service)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
HcalAmplifier * theAmplifier
std::string ecalTBInfoLabel
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:90
HcalCoderFactory * theCoderFactory
ESHandle< TrackerGeometry > geometry
CaloVShape * theHcalShape
virtual void initializeEvent(edm::Event const &e, edm::EventSetup const &c) override
void setHitCorrection(const CaloVHitCorrection *hitCorrection)
If you want to correct hits, for attenuation or delay, set this.
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
void setDetIds(const std::vector< DetId > &detIds)
edm::SortedCollection< HBHEDataFrame > HBHEDigiCollection
HcalTBDigiProducer(const edm::ParameterSet &ps, edm::EDProducer &mixMod)
bool syncPhase() const
choice of the ADC time alignment (synchronous for LHC, asynchronous for test beams) ...