CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZdcHitReconstructor.cc
Go to the documentation of this file.
1 #include "ZdcHitReconstructor.h"
15 
16 #include <iostream>
17 
18 /* Zdc Hit reconstructor allows for CaloRecHits with status words */
19 
21  reco_(conf.getParameter<bool>("correctForTimeslew"),
22  conf.getParameter<bool>("correctForPhaseContainment"),
23  conf.getParameter<double>("correctionPhaseNS"),
24  conf.getParameter<int>("recoMethod"),
25  conf.getParameter<int>("lowGainOffset"),
26  conf.getParameter<double>("lowGainFrac")),
27  det_(DetId::Hcal),
28  inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
29  correctTiming_(conf.getParameter<bool>("correctTiming")),
30  setNoiseFlags_(conf.getParameter<bool>("setNoiseFlags")),
31  setHSCPFlags_(conf.getParameter<bool>("setHSCPFlags")),
32  setSaturationFlags_(conf.getParameter<bool>("setSaturationFlags")),
33  setTimingTrustFlags_(conf.getParameter<bool>("setTimingTrustFlags")),
34  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
35  AuxTSvec_(conf.getParameter<std::vector<int> >("AuxTSvec"))
36 
37 {
38  std::sort(AuxTSvec_.begin(),AuxTSvec_.end()); // sort vector in ascending TS order
39  std::string subd=conf.getParameter<std::string>("Subdetector");
40 
42  {
43  const edm::ParameterSet& pssat = conf.getParameter<edm::ParameterSet>("saturationParameters");
44  saturationFlagSetter_ = new HcalADCSaturationFlag(pssat.getParameter<int>("maxADCvalue"));
45  }
46  if (!strcasecmp(subd.c_str(),"ZDC")) {
49  produces<ZDCRecHitCollection>();
50  } else if (!strcasecmp(subd.c_str(),"CALIB")) {
53  produces<HcalCalibRecHitCollection>();
54  } else {
55  std::cout << "ZdcHitReconstructor is not associated with a specific subdetector!" << std::endl;
56  }
57 
58 }
59 
61 }
63 
65  es.get<HcalLongRecoParamsRcd>().get(p);
67 }
68 
70  if (myobject) delete myobject;
71 }
73 {
74  // get conditions
76  eventSetup.get<HcalDbRecord>().get(conditions);
77  const HcalQIEShape* shape = conditions->getHcalShape (); // this one is generic
78 
80  eventSetup.get<HcalChannelQualityRcd>().get(p);
81  HcalChannelQuality* myqual = new HcalChannelQuality(*p.product());
82 
84  eventSetup.get<HcalSeverityLevelComputerRcd>().get(mycomputer);
85  const HcalSeverityLevelComputer* mySeverity = mycomputer.product();
86 
87  // define vectors to pass noiseTS and signalTS
88  std::vector<unsigned int> mySignalTS;
89  std::vector<unsigned int> myNoiseTS;
90 
93  e.getByLabel(inputLabel_,digi);
94 
95  // create empty output
96  std::auto_ptr<ZDCRecHitCollection> rec(new ZDCRecHitCollection);
97  rec->reserve(digi->size());
98  // run the algorithm
100  for (i=digi->begin(); i!=digi->end(); i++) {
101  HcalZDCDetId cell = i->id();
102  DetId detcell=(DetId)cell;
103  // check on cells to be ignored and dropped: (rof,20.Feb.09)
104  const HcalChannelStatus* mydigistatus=myqual->getValues(detcell.rawId());
105  if (mySeverity->dropChannel(mydigistatus->getValue() ) ) continue;
107  if (i->zsMarkAndPass()) continue;
108  const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
109  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
110  HcalCoderDb coder (*channelCoder, *shape);
111 
112 // get db values for signalTSs and noiseTSs
113  const HcalLongRecoParam* myParams = myobject->getValues(detcell);
114  mySignalTS.clear();
115  myNoiseTS.clear();
116  mySignalTS = myParams->signalTS();
117  myNoiseTS = myParams->noiseTS();
118 
119  rec->push_back(reco_.reconstruct(*i,myNoiseTS,mySignalTS,coder,calibrations));
120  (rec->back()).setFlags(0);
122  saturationFlagSetter_->setSaturationFlag(rec->back(),*i);
123 
124  // Set auxiliary flag with subset of digi information
125  // ZDC aux flag can store non-contiguous set of values
126  int auxflag=0;
127  for (unsigned int xx=0; xx<AuxTSvec_.size() && xx<4;++xx)
128  {
129  if (AuxTSvec_[xx]<0 || AuxTSvec_[xx]>9) continue; // don't allow
130  auxflag+=(i->sample(AuxTSvec_[xx]).adc())<<(7*xx); // store the time slices in the first 28 bits of aux, a set of 4 7-bit a dc values
131  }
132  // bits 28 and 29 are reserved for capid of the first time slice saved in aux
133  if (AuxTSvec_.size()>0)
134  auxflag+=((i->sample(AuxTSvec_[0]).capid())<<28);
135  (rec->back()).setAux(auxflag);
136  }
137  // return result
138  e.put(rec);
139  } // else if (det_==DetId::Calo...)
140 
141  delete myqual;
142 } // void HcalHitReconstructor::produce(...)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< unsigned int > signalTS() const
std::vector< T >::const_iterator const_iterator
virtual void beginRun(edm::Run &r, edm::EventSetup const &es)
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
ZDCRecHit reconstruct(const ZDCDataFrame &digi, const std::vector< unsigned int > &myNoiseTS, const std::vector< unsigned int > &mySignalTS, const HcalCoder &coder, const HcalCalibrations &calibs) const
ZdcSimpleRecAlgo reco_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
HcalOtherSubdetector subdetOther_
HcalLongRecoParams * myobject
bool dropChannel(const uint32_t &mystatus) const
void setSaturationFlag(HBHERecHit &rechit, const HBHEDataFrame &digi)
HcalADCSaturationFlag * saturationFlagSetter_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
tuple conf
Definition: dbtoconf.py:185
Definition: DetId.h:20
virtual void endRun(edm::Run &r, edm::EventSetup const &es)
static const int SubdetectorId
Definition: HcalZDCDetId.h:22
ZdcHitReconstructor(const edm::ParameterSet &ps)
virtual void produce(edm::Event &e, const edm::EventSetup &c)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
std::vector< int > AuxTSvec_
std::vector< unsigned int > noiseTS() const
tuple cout
Definition: gather_cfg.py:121
uint32_t getValue() const
const Item * getValues(DetId fId) const
Definition: Run.h:33