CMS 3D CMS Logo

CastorSimpleReconstructor.cc
Go to the documentation of this file.
25 
26 #include <iostream>
27 
29 public:
30  explicit CastorSimpleReconstructor(const edm::ParameterSet& ps);
31  ~CastorSimpleReconstructor() override;
32  void produce(edm::Event& e, const edm::EventSetup& c) override;
33 
34 private:
37  int subdet_;
38  // HcalOtherSubdetector subdetOther_;
43 
47  bool tsFromDB_;
50 };
51 
54 
55 using namespace std;
56 
58  : reco_(conf.getParameter<int>("firstSample"),
59  conf.getParameter<int>("samplesToAdd"),
60  conf.getParameter<bool>("correctForTimeslew"),
61  conf.getParameter<bool>("correctForPhaseContainment"),
62  conf.getParameter<double>("correctionPhaseNS")),
63  det_(DetId::Hcal),
64  firstSample_(conf.getParameter<int>("firstSample")),
65  samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
66  maxADCvalue_(conf.getParameter<int>("maxADCvalue")),
67  tsFromDB_(conf.getParameter<bool>("tsFromDB")),
68  setSaturationFlag_(conf.getParameter<bool>("setSaturationFlag")),
69  doSaturationCorr_(conf.getParameter<bool>("doSaturationCorr")) {
70  tok_input_ = consumes<CastorDigiCollection>(conf.getParameter<edm::InputTag>("digiLabel"));
71  tok_conditions_ = esConsumes<CastorDbService, CastorDbRecord>();
72 
73  std::string subd = conf.getParameter<std::string>("Subdetector");
74  if (!strcasecmp(subd.c_str(), "CASTOR")) {
75  det_ = DetId::Calo;
77  produces<CastorRecHitCollection>();
78  } else {
79  edm::LogWarning("CastorSimpleReconstructor")
80  << "CastorSimpleReconstructor is not associated with CASTOR subdetector!" << std::endl;
81  }
82  if (tsFromDB_) {
83  tok_recoParams_ = esConsumes<CastorRecoParams, CastorRecoParamsRcd>();
84  }
85  if (doSaturationCorr_) {
86  tok_satCorr_ = esConsumes<CastorSaturationCorrs, CastorSaturationCorrsRcd>();
87  }
88 }
89 
91 
93  // get conditions
95  const CastorQIEShape* shape = conditions->getCastorShape(); // this one is generic
96 
98 
99  // try to get the TS windows from the db
101  if (tsFromDB_) {
102  recoparams = eventSetup.getHandle(tok_recoParams_);
103  if (!recoparams.isValid()) {
104  tsFromDB_ = false;
105  edm::LogWarning("CastorSimpleReconstructor")
106  << "Could not handle the CastorRecoParamsRcd correctly, using parameters from cfg file from this event "
107  "onwards... These parameters could be wrong for this run... please check"
108  << std::endl;
109  }
110  }
111 
112  // try to get the saturation correction constants from the db
114  if (doSaturationCorr_) {
115  satcorr = eventSetup.getHandle(tok_satCorr_);
116  if (!satcorr.isValid()) {
117  doSaturationCorr_ = false;
118  edm::LogWarning("CastorSimpleReconstructor") << "Could not handle the CastorSaturationCorrsRcd correctly. We'll "
119  "not try the saturation correction from this event onwards..."
120  << std::endl;
121  }
122  }
123 
126  e.getByToken(tok_input_, digi);
127 
128  // create empty output
129  auto rec = std::make_unique<CastorRecHitCollection>();
130  // run the algorithm
132  for (i = digi->begin(); i != digi->end(); i++) {
133  HcalCastorDetId cell = i->id();
134  DetId detcell = (DetId)cell;
135  const CastorCalibrations& calibrations = conditions->getCastorCalibrations(cell);
136 
137  if (tsFromDB_) {
138  const CastorRecoParam* param_ts = recoparams->getValues(detcell.rawId());
139  reco_.resetTimeSamples(param_ts->firstSample(), param_ts->samplesToAdd());
140  }
141  const CastorQIECoder* channelCoder = conditions->getCastorCoder(cell);
142  CastorCoderDb coder(*channelCoder, *shape);
143 
144  // reconstruct the rechit
145  rec->push_back(reco_.reconstruct(*i, coder, calibrations));
146 
147  // set the saturation flag if needed
148  if (setSaturationFlag_) {
149  reco_.checkADCSaturation(rec->back(), *i, maxADCvalue_);
150 
151  //++++ Saturation Correction +++++
152  if (doSaturationCorr_ && rec->back().flagField(HcalCaloFlagLabels::ADCSaturationBit)) {
153  // get saturation correction value
154  const CastorSaturationCorr* saturationCorr = satcorr->getValues(detcell.rawId());
155  double satCorrConst = 1.;
156  satCorrConst = saturationCorr->getValue();
157  reco_.recoverADCSaturation(rec->back(), coder, calibrations, *i, maxADCvalue_, satCorrConst);
158  }
159  }
160  }
161  // return result
162  e.put(std::move(rec));
163  }
164 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void resetTimeSamples(int f, int t)
const Item * getValues(DetId fId, bool throwOnFail=true) const
void recoverADCSaturation(CastorRecHit &rechit, const CastorCoder &coder, const CastorCalibrations &calibs, const CastorDataFrame &digi, const int &maxADCvalue, const double &satCorrConst) const
std::vector< T >::const_iterator const_iterator
CastorSimpleReconstructor(const edm::ParameterSet &ps)
edm::ESGetToken< CastorSaturationCorrs, CastorSaturationCorrsRcd > tok_satCorr_
unsigned int samplesToAdd() const
unsigned int firstSample() const
void produce(edm::Event &e, const edm::EventSetup &c) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< CastorDigiCollection > tok_input_
const_iterator begin() const
static const int SubdetectorId
bool isValid() const
Definition: ESHandle.h:44
const_iterator end() const
edm::ESGetToken< CastorRecoParams, CastorRecoParamsRcd > tok_recoParams_
Definition: DetId.h:17
edm::ESGetToken< CastorDbService, CastorDbRecord > tok_conditions_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
Detector
Definition: DetId.h:24
void checkADCSaturation(CastorRecHit &rechit, const CastorDataFrame &digi, const int &maxADCvalue) const
CastorRecHit reconstruct(const CastorDataFrame &digi, const CastorCoder &coder, const CastorCalibrations &calibs) const
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511