CMS 3D CMS Logo

CastorSimpleReconstructor.cc
Go to the documentation of this file.
19 
21 #include <iostream>
22 
23 using namespace std;
24 
26  : reco_(conf.getParameter<int>("firstSample"),
27  conf.getParameter<int>("samplesToAdd"),
28  conf.getParameter<bool>("correctForTimeslew"),
29  conf.getParameter<bool>("correctForPhaseContainment"),
30  conf.getParameter<double>("correctionPhaseNS")),
31  det_(DetId::Hcal),
32  firstSample_(conf.getParameter<int>("firstSample")),
33  samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
34  maxADCvalue_(conf.getParameter<int>("maxADCvalue")),
35  tsFromDB_(conf.getParameter<bool>("tsFromDB")),
36  setSaturationFlag_(conf.getParameter<bool>("setSaturationFlag")),
37  doSaturationCorr_(conf.getParameter<bool>("doSaturationCorr")) {
38  tok_input_ = consumes<CastorDigiCollection>(conf.getParameter<edm::InputTag>("digiLabel"));
39 
40  std::string subd = conf.getParameter<std::string>("Subdetector");
41  if (!strcasecmp(subd.c_str(), "CASTOR")) {
42  det_ = DetId::Calo;
44  produces<CastorRecHitCollection>();
45  } else {
46  edm::LogWarning("CastorSimpleReconstructor")
47  << "CastorSimpleReconstructor is not associated with CASTOR subdetector!" << std::endl;
48  }
49 }
50 
52 
54  // get conditions
56  eventSetup.get<CastorDbRecord>().get(conditions);
57  const CastorQIEShape* shape = conditions->getCastorShape(); // this one is generic
58 
60 
61  // try to get the TS windows from the db
63  if (tsFromDB_) {
64  eventSetup.get<CastorRecoParamsRcd>().get(recoparams);
65  if (!recoparams.isValid()) {
66  tsFromDB_ = false;
67  edm::LogWarning("CastorSimpleReconstructor")
68  << "Could not handle the CastorRecoParamsRcd correctly, using parameters from cfg file from this event "
69  "onwards... These parameters could be wrong for this run... please check"
70  << std::endl;
71  }
72  }
73 
74  // try to get the saturation correction constants from the db
76  if (doSaturationCorr_) {
77  eventSetup.get<CastorSaturationCorrsRcd>().get(satcorr);
78  if (!satcorr.isValid()) {
79  doSaturationCorr_ = false;
80  edm::LogWarning("CastorSimpleReconstructor") << "Could not handle the CastorSaturationCorrsRcd correctly. We'll "
81  "not try the saturation correction from this event onwards..."
82  << std::endl;
83  }
84  }
85 
88  e.getByToken(tok_input_, digi);
89 
90  // create empty output
91  auto rec = std::make_unique<CastorRecHitCollection>();
92  // run the algorithm
94  for (i = digi->begin(); i != digi->end(); i++) {
95  HcalCastorDetId cell = i->id();
96  DetId detcell = (DetId)cell;
97  const CastorCalibrations& calibrations = conditions->getCastorCalibrations(cell);
98 
99  if (tsFromDB_) {
100  const CastorRecoParam* param_ts = recoparams->getValues(detcell.rawId());
101  reco_.resetTimeSamples(param_ts->firstSample(), param_ts->samplesToAdd());
102  }
103  const CastorQIECoder* channelCoder = conditions->getCastorCoder(cell);
104  CastorCoderDb coder(*channelCoder, *shape);
105 
106  // reconstruct the rechit
107  rec->push_back(reco_.reconstruct(*i, coder, calibrations));
108 
109  // set the saturation flag if needed
110  if (setSaturationFlag_) {
111  reco_.checkADCSaturation(rec->back(), *i, maxADCvalue_);
112 
113  //++++ Saturation Correction +++++
114  if (doSaturationCorr_ && rec->back().flagField(HcalCaloFlagLabels::ADCSaturationBit)) {
115  // get saturation correction value
116  const CastorSaturationCorr* saturationCorr = satcorr->getValues(detcell.rawId());
117  double satCorrConst = 1.;
118  satCorrConst = saturationCorr->getValue();
119  reco_.recoverADCSaturation(rec->back(), coder, calibrations, *i, maxADCvalue_, satCorrConst);
120  }
121  }
122  }
123  // return result
124  e.put(std::move(rec));
125  }
126 }
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
void resetTimeSamples(int f, int t)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< T >::const_iterator const_iterator
const CastorCalibrations & getCastorCalibrations(const HcalGenericDetId &fId) const
CastorSimpleReconstructor(const edm::ParameterSet &ps)
const Item * getValues(DetId fId, bool throwOnFail=true) const
CastorRecHit reconstruct(const CastorDataFrame &digi, const CastorCoder &coder, const CastorCalibrations &calibs) const
void produce(edm::Event &e, const edm::EventSetup &c) override
unsigned int samplesToAdd() const
const CastorQIEShape * getCastorShape() const
edm::EDGetTokenT< CastorDigiCollection > tok_input_
static const int SubdetectorId
const_iterator end() const
Definition: DetId.h:17
unsigned int firstSample() const
void recoverADCSaturation(CastorRecHit &rechit, const CastorCoder &coder, const CastorCalibrations &calibs, const CastorDataFrame &digi, const int &maxADCvalue, const double &satCorrConst) const
T get() const
Definition: EventSetup.h:73
void checkADCSaturation(CastorRecHit &rechit, const CastorDataFrame &digi, const int &maxADCvalue) const
bool isValid() const
Definition: ESHandle.h:44
def move(src, dest)
Definition: eostools.py:511
const_iterator begin() const
const CastorQIECoder * getCastorCoder(const HcalGenericDetId &fId) const