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.
6 
16 
18  theParameterMap(new HcalTBSimParameterMap(ps)),
19  theHcalShape(new HcalShape()),
20  theHcalIntegratedShape(new CaloShapeIntegrator(theHcalShape)),
21  theHBHEResponse(new CaloHitResponse(theParameterMap, theHcalIntegratedShape)),
22  theHOResponse(new CaloHitResponse(theParameterMap, theHcalIntegratedShape)),
23  theAmplifier(0), theCoderFactory(0), theElectronicsSim(0),
24  theHitCorrection(0), theHBHEDigitizer(0), theHODigitizer(0), theHBHEHits(),
25  theHOHits(), thisPhaseShift(0) {
26  std::string const instance("simHcalDigis");
27  mixMod.produces<HBHEDigiCollection>(instance);
28  mixMod.produces<HODigiCollection>(instance);
29  iC.consumes<std::vector<PCaloHit> >(edm::InputTag("g4SimHits", "HcalHits"));
30 
31  DetId detId(DetId::Hcal, 1);
32  bool syncPhase = (theParameterMap->simParameters(detId)).syncPhase();
33  doPhaseShift = !syncPhase;
34 
37 
38  bool doTimeSlew = ps.getParameter<bool>("doTimeSlew");
39  if(doTimeSlew) {
40  // no time slewing for HF
44  }
45 
46  bool doNoise = ps.getParameter<bool>("doNoise");
47  bool dummy1 = false;
48  bool dummy2 = false; // extra arguments for premixing
49  theAmplifier = new HcalAmplifier(theParameterMap, doNoise, dummy1, dummy2);
52 
55 
56  tunePhaseShift = ps.getUntrackedParameter<double>("tunePhaseShiftTB", 1.);
57  ecalTBInfoLabel = ps.getUntrackedParameter<std::string>("EcalTBInfoLabel","SimEcalTBG4Object");
58  edm::LogInfo("HcalSim") << "HcalTBDigiProducer initialized with doNoise = "
59  << doNoise << ", doTimeSlew = " << doTimeSlew
60  << " and doPhaseShift = " << doPhaseShift
61  << " tunePhasShift = " << tunePhaseShift;
62 
63  if (doPhaseShift) {
65  }
66 }
67 
69 
71  if (theHODigitizer) delete theHODigitizer;
72  if (theParameterMap) delete theParameterMap;
73  if (theHcalShape) delete theHcalShape;
75  if (theHBHEResponse) delete theHBHEResponse;
76  if (theHOResponse) delete theHOResponse;
78  if (theAmplifier) delete theAmplifier;
79  if (theCoderFactory) delete theCoderFactory;
81 }
82 
83 
85  // get the appropriate gains, noises, & widths for this event
87  eventSetup.get<HcalDbRecord>().get(conditions);
88  theAmplifier->setDbService(conditions.product());
89  theCoderFactory->setDbService(conditions.product());
90 
91  // get the correct geometry
92  checkGeometry(eventSetup);
93 
94  theHBHEHits.clear();
95  theHOHits.clear();
96  if (doPhaseShift) {
97 
98  edm::Handle<PEcalTBInfo> theEcalTBInfo;
99  e.getByLabel(ecalTBInfoLabel,theEcalTBInfo);
100  thisPhaseShift = theEcalTBInfo->phaseShift();
101 
102  DetId detIdHB(DetId::Hcal, 1);
103  setPhaseShift(detIdHB);
104  DetId detIdHO(DetId::Hcal, 3);
105  setPhaseShift(detIdHO);
106  }
107 
110 }
111 
112 void HcalTBDigiProducer::accumulateCaloHits(edm::Handle<std::vector<PCaloHit> > const& hcalHandle, int bunchCrossing) {
113 
114  LogDebug("HcalSim") << "HcalTBDigiProducer::accumulate trying to get SimHit";
115 
116  if(hcalHandle.isValid()) {
117  std::vector<PCaloHit> hits = *hcalHandle.product();
118  if(theHitCorrection != 0) {
120  }
121  LogDebug("HcalSim") << "HcalTBDigiProducer::accumulate Hits corrected";
122  theHBHEDigitizer->add(hits, bunchCrossing);
123  theHODigitizer->add(hits, bunchCrossing);
124  }
125 }
126 
128  // Step A: Get Inputs, and accumulate digis
129 
130  edm::InputTag hcalTag("g4SimHits", "HcalHits");
132  e.getByLabel(hcalTag, hcalHandle);
133 
134  accumulateCaloHits(hcalHandle, 0);
135 }
136 
138  // Step A: Get Inputs, and accumulate digis
139 
140  edm::InputTag hcalTag("g4SimHits", "HcalHits");
142  e.getByLabel(hcalTag, hcalHandle);
143 
144  accumulateCaloHits(hcalHandle, e.bunchCrossing());
145 }
146 
148  // Step B: Create empty output
149  std::auto_ptr<HBHEDigiCollection> hbheResult(new HBHEDigiCollection());
150  std::auto_ptr<HODigiCollection> hoResult(new HODigiCollection());
151  LogDebug("HcalSim") << "HcalTBDigiProducer::produce Empty collection created";
152  // Step C: Invoke the algorithm, getting back outputs.
153  theHBHEDigitizer->run(*hbheResult);
154  edm::LogInfo("HcalSim") << "HcalTBDigiProducer: HBHE digis : " << hbheResult->size();
155  theHODigitizer->run(*hoResult);
156  edm::LogInfo("HcalSim") << "HcalTBDigiProducer: HO digis : " << hoResult->size();
157 
158  // Step D: Put outputs into event
159  std::string const instance("simHcalDigis");
160  e.put(hbheResult, instance);
161  e.put(hoResult, instance);
162 
163 }
164 
166 
167  for (edm::PCaloHitContainer::const_iterator hitItr = hits.begin();
168  hitItr != hits.end(); ++hitItr) {
169  HcalSubdetector subdet = HcalDetId(hitItr->id()).subdet();
170  if(subdet == HcalBarrel || subdet == HcalEndcap) {
171  theHBHEHits.push_back(*hitItr);
172  } else if(subdet == HcalOuter) {
173  theHOHits.push_back(*hitItr);
174  } else {
175  edm::LogError("HcalSim") << "Bad HcalHit subdetector " << subdet;
176  }
177  }
178 }
179 
181 
182  // TODO find a way to avoid doing this every event
184  eventSetup.get<CaloGeometryRecord>().get(geometry);
185 
186  const CaloGeometry * pGeometry = &*geometry;
187 
188  // see if we need to update
189  if(pGeometry != theGeometry) {
190  theGeometry = pGeometry;
191  updateGeometry();
192  }
193 }
194 
196 
199 
200  // Get cells for HB and HE
201  hbheCells.clear();
203  std::vector<DetId> heCells = theGeometry->getValidDetIds(DetId::Hcal, HcalEndcap);
204  // combine HB & HE
205  hbheCells.insert(hbheCells.end(), heCells.begin(), heCells.end());
206 
207  // Get cells for HO
208  hoCells.clear();
210 
211  edm::LogInfo("HcalSim") << "HcalTBDigiProducer update Geometry with "
212  << hbheCells.size() << " cells in HB/HE and "
213  << hoCells.size() << " cells in HO";
214 
216  LogDebug("HcalSim") << "HcalTBDigiProducer: Set DetID's for HB/HE";
218  LogDebug("HcalSim") << "HcalTBDigiProducer: Set DetID's for HO";
219 }
220 
222 
224  if ( !parameters.syncPhase() ) {
225  int myDet = detId.subdetId();
226  double passPhaseShift = thisPhaseShift + tunePhaseShift;
227  if (myDet <= 2) {
228  theHBHEResponse->setPhaseShift(passPhaseShift);
229  } else {
230  theHOResponse->setPhaseShift(passPhaseShift);
231  }
232  }
233 }
#define LogDebug(id)
void setPhaseShift(const DetId &detId)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
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
HcalTBDigiProducer(const edm::ParameterSet &ps, edm::one::EDProducerBase &mixMod, edm::ConsumesCollector &iC)
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:116
HBHEDigitizer * theHBHEDigitizer
void fillChargeSums(MixCollection< PCaloHit > &hits)
HcalSubdetector
Definition: HcalAssistant.h:31
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:37
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
virtual void accumulate(edm::Event const &e, edm::EventSetup const &c) override
std::vector< PCaloHit > theHOHits
Definition: DetId.h:18
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
bool syncPhase() const
choice of the ADC time alignment (synchronous for LHC, asynchronous for test beams) ...