![]() |
![]() |
#include <Work/EcalMIPRecHitFilter/src/EcalMIPRecHitFilter.cc>
Public Member Functions | |
EcalMIPRecHitFilter (const edm::ParameterSet &) | |
~EcalMIPRecHitFilter () | |
Private Member Functions | |
virtual void | beginJob () |
virtual void | endJob () |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
Private Attributes | |
edm::InputTag | EcalRecHitCollection_ |
std::vector< int > | maskedList_ |
double | minAmp1_ |
double | minAmp2_ |
double | minSingleAmp_ |
int | side_ |
Description: <one line="" class="" summary>="">
Implementation: <Notes on="" implementation>="">
Definition at line 55 of file EcalMIPRecHitFilter.cc.
EcalMIPRecHitFilter::EcalMIPRecHitFilter | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 79 of file EcalMIPRecHitFilter.cc.
References EcalRecHitCollection_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), maskedList_, minAmp1_, minAmp2_, minSingleAmp_, and side_.
{ //now do what ever initialization is needed minSingleAmp_ = iConfig.getUntrackedParameter<double>("SingleAmpMin", 0.108); minAmp1_ = iConfig.getUntrackedParameter<double>("AmpMinSeed", 0.063); minAmp2_ = iConfig.getUntrackedParameter<double>("AmpMin2", 0.045); maskedList_ = iConfig.getUntrackedParameter<std::vector<int> >("maskedChannels", maskedList_); //this is using the ashed index EcalRecHitCollection_ = iConfig.getParameter<edm::InputTag>("EcalRecHitCollection"); side_ = iConfig.getUntrackedParameter<int>("side", 3); }
EcalMIPRecHitFilter::~EcalMIPRecHitFilter | ( | ) |
Definition at line 91 of file EcalMIPRecHitFilter.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
void EcalMIPRecHitFilter::beginJob | ( | void | ) | [private, virtual] |
void EcalMIPRecHitFilter::endJob | ( | void | ) | [private, virtual] |
bool EcalMIPRecHitFilter::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
LASER and CALIB constants from the DB //PUT THRESHOLDS IN ADC... AGAIN.
LASER and CALIB constants from the DB
Implements HLTFilter.
Definition at line 106 of file EcalMIPRecHitFilter.cc.
References EcalRecHitCollection_, CaloRecHit::energy(), spr::find(), edm::EventSetup::get(), edm::Event::getByLabel(), EcalADCToGeVConstant::getEBValue(), EcalCondObjectContainer< T >::getMap(), ecalpyutils::hashedIndex(), EcalRecHit::id(), maskedList_, minAmp1_, minAmp2_, minSingleAmp_, query::result, side_, and edm::EventBase::time().
{ using namespace edm; // getting very basic uncalRH Handle<EcalRecHitCollection> recHits; //try { // iEvent.getByLabel(EcalRecHitCollection_, recHits); //} catch ( std::exception& ex) { // LogWarning("EcalMIPRecHitFilter") << EcalRecHitCollection_ << " not available"; //} iEvent.getByLabel(EcalRecHitCollection_, recHits); if (!recHits.isValid()){ LogWarning("EcalMIPRecHitFilter") << EcalRecHitCollection_ << " not available"; } ESHandle<CaloTopology> caloTopo; iSetup.get<CaloTopologyRecord>().get(caloTopo); // Intercalib constants edm::ESHandle<EcalIntercalibConstants> pIcal; iSetup.get<EcalIntercalibConstantsRcd>().get(pIcal); const EcalIntercalibConstants* ical = pIcal.product(); const EcalIntercalibConstantMap& icalMap=ical->getMap(); edm::ESHandle<EcalLaserDbService> pLaser; iSetup.get<EcalLaserDbRecord>().get( pLaser ); edm::ESHandle<EcalADCToGeVConstant> pAgc; iSetup.get<EcalADCToGeVConstantRcd>().get(pAgc); const EcalADCToGeVConstant* agc = pAgc.product(); //std::cout << "Global EB ADC->GeV scale: " << agc->getEBValue() << " GeV/ADC count" ; float adcconst = agc->getEBValue(); bool thereIsSignal = false; // loop on rechits for ( EcalRecHitCollection::const_iterator hitItr = recHits->begin(); hitItr != recHits->end(); ++hitItr ) { EcalRecHit hit = (*hitItr); // masking noisy channels //KEEP this for now, just in case a few show up std::vector<int>::iterator result; result = find( maskedList_.begin(), maskedList_.end(), EBDetId(hit.id()).hashedIndex() ); if (result != maskedList_.end()) // LogWarning("EcalFilter") << "skipping uncalRecHit for channel: " << ic << " with amplitude " << ampli_ ; continue; float ampli_ = hit.energy(); EBDetId ebDet = hit.id(); // find intercalib constant for this xtal EcalIntercalibConstantMap::const_iterator icalit=icalMap.find(ebDet); EcalIntercalibConstant icalconst = 1.; if( icalit!=icalMap.end() ){ icalconst = (*icalit); // LogDebug("EcalRecHitDebug") << "Found intercalib for xtal " << EBDetId(it->id()).ic() << " " << icalconst ; } else { //edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EBDetId(ebDet) << "! something wrong with EcalIntercalibConstants in your DB? " ; } float lasercalib = pLaser->getLaserCorrection( EBDetId(ebDet), iEvent.time() ); ampli_ /= (icalconst * lasercalib * adcconst); // seeking channels with signal and displaced jitter if (ampli_ >= minSingleAmp_ ) { //std::cout << " THIS AMPLITUDE WORKS " << ampli_ << std::endl; thereIsSignal = true; // LogWarning("EcalFilter") << "at evet: " << iEvent.id().event() // << " and run: " << iEvent.id().run() // << " there is OUT OF TIME signal at chanel: " << ic // << " with amplitude " << ampli_ << " and max at: " << jitter_; break; } //Check for more robust selection other than just single crystal cosmics if (ampli_ >= minAmp1_) { //std::cout << " THIS AMPLITUDE WORKS " << ampli_ << std::endl; std::vector<DetId> neighbors = caloTopo->getWindow(ebDet,side_,side_); float secondMin = 0.; for(std::vector<DetId>::const_iterator detitr = neighbors.begin(); detitr != neighbors.end(); ++detitr) { EcalRecHitCollection::const_iterator thishit = recHits->find((*detitr)); if (thishit == recHits->end()) { //LogWarning("EcalMIPRecHitFilter") << "No RecHit available, for "<< EBDetId(*detitr); continue; } if ((*thishit).id() != ebDet) { float thisamp = (*thishit).energy(); // find intercalib constant for this xtal EcalIntercalibConstantMap::const_iterator icalit2=icalMap.find((*thishit).id()); EcalIntercalibConstant icalconst2 = 1.; if( icalit2!=icalMap.end() ){ icalconst2 = (*icalit2); // LogDebug("EcalRecHitDebug") << "Found intercalib for xtal " << EBDetId(it->id()).ic() << " " << icalconst ; } else { //edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EBDetId(ebDet) << "! something wrong with EcalIntercalibConstants in your DB? " ; } float lasercalib2 = pLaser->getLaserCorrection( EBDetId((*thishit).id()), iEvent.time() ); thisamp /= (icalconst2 * lasercalib2 * adcconst); if (thisamp > secondMin) secondMin = thisamp; } } if (secondMin > minAmp2_ ) { thereIsSignal = true; break; } } } //std::cout << " Ok is There one of THEM " << thereIsSignal << std::endl; return thereIsSignal; }
Definition at line 67 of file EcalMIPRecHitFilter.cc.
Referenced by EcalMIPRecHitFilter(), and filter().
std::vector<int> EcalMIPRecHitFilter::maskedList_ [private] |
Definition at line 71 of file EcalMIPRecHitFilter.cc.
Referenced by EcalMIPRecHitFilter(), and filter().
double EcalMIPRecHitFilter::minAmp1_ [private] |
Definition at line 68 of file EcalMIPRecHitFilter.cc.
Referenced by EcalMIPRecHitFilter(), and filter().
double EcalMIPRecHitFilter::minAmp2_ [private] |
Definition at line 69 of file EcalMIPRecHitFilter.cc.
Referenced by EcalMIPRecHitFilter(), and filter().
double EcalMIPRecHitFilter::minSingleAmp_ [private] |
Definition at line 70 of file EcalMIPRecHitFilter.cc.
Referenced by EcalMIPRecHitFilter(), and filter().
int EcalMIPRecHitFilter::side_ [private] |
Definition at line 72 of file EcalMIPRecHitFilter.cc.
Referenced by EcalMIPRecHitFilter(), and filter().