CMS 3D CMS Logo

Public Member Functions | Private Attributes

CastorSimpleReconstructor Class Reference

#include <CastorSimpleReconstructor.h>

Inheritance diagram for CastorSimpleReconstructor:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual void beginRun (edm::Run &r, edm::EventSetup const &es)
 CastorSimpleReconstructor (const edm::ParameterSet &ps)
virtual void produce (edm::Event &e, const edm::EventSetup &c)
virtual ~CastorSimpleReconstructor ()

Private Attributes

DetId::Detector det_
bool doSaturationCorr_
int firstSample_
edm::InputTag inputLabel_
int maxADCvalue_
CastorSimpleRecAlgo reco_
int samplesToAdd_
bool setSaturationFlag_
int subdet_
bool tsFromDB_

Detailed Description

Definition at line 13 of file CastorSimpleReconstructor.h.


Constructor & Destructor Documentation

CastorSimpleReconstructor::CastorSimpleReconstructor ( const edm::ParameterSet ps) [explicit]

Definition at line 26 of file CastorSimpleReconstructor.cc.

References DetId::Calo, det_, edm::ParameterSet::getParameter(), subdet_, and HcalCastorDetId::SubdetectorId.

                                                                               :
reco_(conf.getParameter<int>("firstSample"),conf.getParameter<int>("samplesToAdd"),conf.getParameter<bool>("correctForTimeslew"),
      conf.getParameter<bool>("correctForPhaseContainment"),conf.getParameter<double>("correctionPhaseNS")),
det_(DetId::Hcal),
inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
firstSample_(conf.getParameter<int>("firstSample")),
samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
tsFromDB_(conf.getParameter<bool>("tsFromDB")),
setSaturationFlag_(conf.getParameter<bool>("setSaturationFlag")),
maxADCvalue_(conf.getParameter<int>("maxADCvalue")),
doSaturationCorr_(conf.getParameter<bool>("doSaturationCorr"))
{
    std::string subd=conf.getParameter<std::string>("Subdetector");
    if (!strcasecmp(subd.c_str(),"CASTOR")) {
        det_=DetId::Calo;
        subdet_=HcalCastorDetId::SubdetectorId;
        produces<CastorRecHitCollection>();
    } else {
        edm::LogWarning("CastorSimpleReconstructor") << "CastorSimpleReconstructor is not associated with CASTOR subdetector!" << std::endl;
    }       
    
}
CastorSimpleReconstructor::~CastorSimpleReconstructor ( ) [virtual]

Definition at line 49 of file CastorSimpleReconstructor.cc.

                                                      {
}

Member Function Documentation

void CastorSimpleReconstructor::beginRun ( edm::Run r,
edm::EventSetup const &  es 
) [virtual]

Reimplemented from edm::EDProducer.

Definition at line 52 of file CastorSimpleReconstructor.cc.

                                                                          {
    
    
}
void CastorSimpleReconstructor::produce ( edm::Event e,
const edm::EventSetup c 
) [virtual]

Implements edm::EDProducer.

Definition at line 56 of file CastorSimpleReconstructor.cc.

References HcalCaloFlagLabels::ADCSaturationBit, DetId::Calo, CastorSimpleRecAlgo::checkADCSaturation(), det_, doSaturationCorr_, CastorRecoParam::firstSample(), edm::EventSetup::get(), edm::Event::getByLabel(), CastorSaturationCorr::getValue(), i, inputLabel_, edm::ESHandleBase::isValid(), maxADCvalue_, edm::Event::put(), DetId::rawId(), reco_, CastorSimpleRecAlgo::reconstruct(), CastorSimpleRecAlgo::recoverADCSaturation(), CastorSimpleRecAlgo::resetTimeSamples(), CastorRecoParam::samplesToAdd(), setSaturationFlag_, subdet_, HcalCastorDetId::SubdetectorId, and tsFromDB_.

{
    // get conditions
    edm::ESHandle<CastorDbService> conditions;
    eventSetup.get<CastorDbRecord>().get(conditions);
    const CastorQIEShape* shape = conditions->getCastorShape (); // this one is generic
    
    CastorCalibrations calibrations;
    
    // try to get the TS windows from the db
    edm::ESHandle<CastorRecoParams> recoparams;
    if (tsFromDB_) {
        eventSetup.get<CastorRecoParamsRcd>().get(recoparams);
        if (!recoparams.isValid()) { 
                tsFromDB_ = false;
                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;
        }
    }
    
    // try to get the saturation correction constants from the db
    edm::ESHandle<CastorSaturationCorrs> satcorr;
    if (doSaturationCorr_) {
        eventSetup.get<CastorSaturationCorrsRcd>().get(satcorr);
        if (!satcorr.isValid()) {
            doSaturationCorr_ = false;
            edm::LogWarning("CastorSimpleReconstructor") << "Could not handle the CastorSaturationCorrsRcd correctly. We'll not try the saturation correction from this event onwards..." << std::endl;
        }
    }
    
    
    if (det_==DetId::Calo && subdet_==HcalCastorDetId::SubdetectorId) {
        edm::Handle<CastorDigiCollection> digi;
        e.getByLabel(inputLabel_,digi);
        
        // create empty output
        std::auto_ptr<CastorRecHitCollection> rec(new CastorRecHitCollection);
        // run the algorithm
        CastorDigiCollection::const_iterator i;
        for (i=digi->begin(); i!=digi->end(); i++) {
            HcalCastorDetId cell = i->id();
            DetId detcell=(DetId)cell;    
            const CastorCalibrations& calibrations=conditions->getCastorCalibrations(cell);
            
            if (tsFromDB_) {
                const CastorRecoParam* param_ts = recoparams->getValues(detcell.rawId());
                reco_.resetTimeSamples(param_ts->firstSample(),param_ts->samplesToAdd());
            }          
            const CastorQIECoder* channelCoder = conditions->getCastorCoder (cell);
            CastorCoderDb coder (*channelCoder, *shape);
            
            // reconstruct the rechit
            rec->push_back(reco_.reconstruct(*i,coder,calibrations));
            
            // set the saturation flag if needed
            if (setSaturationFlag_) {
                reco_.checkADCSaturation(rec->back(),*i,maxADCvalue_);
                
                //++++ Saturation Correction +++++
                if (doSaturationCorr_ && rec->back().flagField(HcalCaloFlagLabels::ADCSaturationBit)) {
                    // get saturation correction value
                    const CastorSaturationCorr* saturationCorr = satcorr->getValues(detcell.rawId());
                    double satCorrConst = 1.;
                    satCorrConst = saturationCorr->getValue();
                    reco_.recoverADCSaturation(rec->back(),coder,calibrations,*i,maxADCvalue_,satCorrConst);
                }
            }
        }
        // return result
        e.put(rec);     
    }
}

Member Data Documentation

Definition at line 21 of file CastorSimpleReconstructor.h.

Referenced by CastorSimpleReconstructor(), and produce().

Definition at line 31 of file CastorSimpleReconstructor.h.

Referenced by produce().

Definition at line 26 of file CastorSimpleReconstructor.h.

Definition at line 24 of file CastorSimpleReconstructor.h.

Referenced by produce().

Definition at line 30 of file CastorSimpleReconstructor.h.

Referenced by produce().

Definition at line 20 of file CastorSimpleReconstructor.h.

Referenced by produce().

Definition at line 27 of file CastorSimpleReconstructor.h.

Definition at line 29 of file CastorSimpleReconstructor.h.

Referenced by produce().

Definition at line 22 of file CastorSimpleReconstructor.h.

Referenced by CastorSimpleReconstructor(), and produce().

Definition at line 28 of file CastorSimpleReconstructor.h.

Referenced by produce().