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.
13 
14 #include <iostream>
15 
16 
18  reco_(conf.getParameter<bool>("correctForTimeslew"),
19  conf.getParameter<bool>("correctForPhaseContainment"),conf.getParameter<double>("correctionPhaseNS"),
20  conf.getParameter<int>("recoMethod"),
21  conf.getParameter<int>("lowGainOffset"),
22  conf.getParameter<double>("lowGainFrac")),
23  det_(DetId::Hcal),
24  inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
25  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed"))
26 {
27  std::string subd=conf.getParameter<std::string>("Subdetector");
28  if (!strcasecmp(subd.c_str(),"ZDC")) {
31  produces<ZDCRecHitCollection>();
32  } else if (!strcasecmp(subd.c_str(),"CALIB")) {
35  produces<HcalCalibRecHitCollection>();
36  } else {
37  std::cout << "ZdcSimpleReconstructor is not associated with a specific subdetector!" << std::endl;
38  }
39 
40 }
41 
43 }
45 
47  es.get<HcalLongRecoParamsRcd>().get(p);
49 }
50 
52  if (myobject) delete myobject;
53 }
55 {
56  // get conditions
58  eventSetup.get<HcalDbRecord>().get(conditions);
59  const HcalQIEShape* shape = conditions->getHcalShape (); // this one is generic
60  // define vectors to pass noiseTS and signalTS
61  std::vector<unsigned int> mySignalTS;
62  std::vector<unsigned int> myNoiseTS;
63 
66  e.getByLabel(inputLabel_,digi);
67 
68  // create empty output
69  std::auto_ptr<ZDCRecHitCollection> rec(new ZDCRecHitCollection);
70  rec->reserve(digi->size());
71  // run the algorithm
72  unsigned int toaddMem = 0;
73 
75  for (i=digi->begin(); i!=digi->end(); i++) {
76  HcalZDCDetId cell = i->id();
77  DetId detcell=(DetId)cell;
78  // rof 27.03.09: drop ZS marked and passed digis:
80  if (i->zsMarkAndPass()) continue;
81 
82 // get db values for signalTSs and noiseTSs
83  const HcalLongRecoParam* myParams = myobject->getValues(detcell);
84  mySignalTS.clear();
85  myNoiseTS.clear();
86  mySignalTS = myParams->signalTS();
87  myNoiseTS = myParams->noiseTS();
88 // warning: the PulseCorrection is not used by ZDC. If it gets a non-contingious set of
89 // signal TS, it may not work properly. Assume contiguous here....
90  unsigned int toadd = mySignalTS.size();
91  if(toaddMem != toadd) {
92  reco_.initPulseCorr(toadd);
93  toaddMem = toadd;
94  }
95  const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
96  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
97  HcalCoderDb coder (*channelCoder, *shape);
98  rec->push_back(reco_.reconstruct(*i,myNoiseTS,mySignalTS,coder,calibrations));
99  }
100  // return result
101  e.put(rec);
102  }
103 }
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
virtual void endRun(edm::Run &r, edm::EventSetup const &es)
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)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
virtual void beginRun(edm::Run &r, edm::EventSetup const &es)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
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
const Item * getValues(DetId fId) const
Definition: Run.h:33