CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZdcSimpleReconstructor.cc
Go to the documentation of this file.
12 
13 #include <iostream>
14 
15 
17  reco_(conf.getParameter<bool>("correctForTimeslew"),
18  conf.getParameter<bool>("correctForPhaseContainment"),conf.getParameter<double>("correctionPhaseNS"),
19  conf.getParameter<int>("recoMethod"),
20  conf.getParameter<int>("lowGainOffset"),
21  conf.getParameter<double>("lowGainFrac")),
22  det_(DetId::Hcal),
23  inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
24  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed"))
25 {
26  std::string subd=conf.getParameter<std::string>("Subdetector");
27  if (!strcasecmp(subd.c_str(),"ZDC")) {
30  produces<ZDCRecHitCollection>();
31  } else if (!strcasecmp(subd.c_str(),"CALIB")) {
34  produces<HcalCalibRecHitCollection>();
35  } else {
36  std::cout << "ZdcSimpleReconstructor is not associated with a specific subdetector!" << std::endl;
37  }
38 
39 }
40 
42 }
44 
46  es.get<HcalLongRecoParamsRcd>().get(p);
48 }
49 
51  delete myobject; myobject = 0;
52 }
54 {
55  // get conditions
57  eventSetup.get<HcalDbRecord>().get(conditions);
58  // define vectors to pass noiseTS and signalTS
59  std::vector<unsigned int> mySignalTS;
60  std::vector<unsigned int> myNoiseTS;
61 
64  e.getByLabel(inputLabel_,digi);
65 
66  // create empty output
67  std::auto_ptr<ZDCRecHitCollection> rec(new ZDCRecHitCollection);
68  rec->reserve(digi->size());
69  // run the algorithm
70  unsigned int toaddMem = 0;
71 
73  for (i=digi->begin(); i!=digi->end(); i++) {
74  HcalZDCDetId cell = i->id();
75  DetId detcell=(DetId)cell;
76  // rof 27.03.09: drop ZS marked and passed digis:
78  if (i->zsMarkAndPass()) continue;
79 
80 // get db values for signalTSs and noiseTSs
81  const HcalLongRecoParam* myParams = myobject->getValues(detcell);
82  mySignalTS.clear();
83  myNoiseTS.clear();
84  mySignalTS = myParams->signalTS();
85  myNoiseTS = myParams->noiseTS();
86 // warning: the PulseCorrection is not used by ZDC. If it gets a non-contingious set of
87 // signal TS, it may not work properly. Assume contiguous here....
88  unsigned int toadd = mySignalTS.size();
89  if(toaddMem != toadd) {
90  reco_.initPulseCorr(toadd);
91  toaddMem = toadd;
92  }
93  const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
94  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
95  const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
96  HcalCoderDb coder (*channelCoder, *shape);
97  rec->push_back(reco_.reconstruct(*i,myNoiseTS,mySignalTS,coder,calibrations));
98  }
99  // return result
100  e.put(rec);
101  }
102 }
virtual void produce(edm::Event &e, const edm::EventSetup &c)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
HcalOtherSubdetector subdetOther_
std::vector< unsigned int > signalTS() const
HcalLongRecoParams * myobject
std::vector< T >::const_iterator const_iterator
const Item * getValues(DetId fId, bool throwOnFail=true) const
virtual void endRun(edm::Run const &r, edm::EventSetup const &es) overridefinal
ZdcSimpleReconstructor(const edm::ParameterSet &ps)
ZDCRecHit reconstruct(const ZDCDataFrame &digi, const std::vector< unsigned int > &myNoiseTS, const std::vector< unsigned int > &mySignalTS, const HcalCoder &coder, const HcalCalibrations &calibs) const
void initPulseCorr(int toadd)
virtual void beginRun(edm::Run const &r, edm::EventSetup const &es) overridefinal
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
tuple conf
Definition: dbtoconf.py:185
Definition: DetId.h:20
static const int SubdetectorId
Definition: HcalZDCDetId.h:22
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
std::vector< unsigned int > noiseTS() const
tuple cout
Definition: gather_cfg.py:121
Definition: Run.h:36