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