Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:09

0001 #include "PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h"
0002 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0003 #include "CommonTools/Egamma/interface/EffectiveAreas.h"
0004 #include "CommonTools/Utils/interface/StringObjectFunction.h"
0005 #include "CommonTools/Utils/interface/ThreadSafeFunctor.h"
0006 #include "RecoEgamma/EgammaTools/interface/EBEECutValues.h"
0007 
0008 class PhoGenericQuadraticRhoPtScaledCut : public CutApplicatorWithEventContentBase {
0009 public:
0010   PhoGenericQuadraticRhoPtScaledCut(const edm::ParameterSet& c);
0011 
0012   result_type operator()(const reco::PhotonPtr&) const final;
0013 
0014   void setConsumes(edm::ConsumesCollector&) final;
0015   void getEventContent(const edm::EventBase&) final;
0016 
0017   double value(const reco::CandidatePtr& cand) const final;
0018 
0019   CandidateType candidateType() const final { return PHOTON; }
0020 
0021 private:
0022   ThreadSafeFunctor<StringObjectFunction<reco::Photon>> varFunc_;
0023   bool lessThan_;
0024   //cut value is constTerm + linearRhoTerm_*rho + linearPtTerm*pt + quadraticPtTerm*pt*pt
0025   //note EBEECutValues & Effective areas are conceptually the same thing, both are eta
0026   //binned constants, just Effective areas have arbitary rather than barrel/endcap binnng
0027   EBEECutValues constTerm_;
0028   EffectiveAreas rhoEA_;
0029   EBEECutValues linearPtTerm_;
0030   EBEECutValues quadraticPtTerm_;
0031 
0032   edm::Handle<double> rhoHandle_;
0033 };
0034 
0035 DEFINE_EDM_PLUGIN(CutApplicatorFactory, PhoGenericQuadraticRhoPtScaledCut, "PhoGenericQuadraticRhoPtScaledCut");
0036 
0037 PhoGenericQuadraticRhoPtScaledCut::PhoGenericQuadraticRhoPtScaledCut(const edm::ParameterSet& params)
0038     : CutApplicatorWithEventContentBase(params),
0039       varFunc_(params.getParameter<std::string>("cutVariable")),
0040       lessThan_(params.getParameter<bool>("lessThan")),
0041       constTerm_(params, "constTerm"),
0042       rhoEA_(params.getParameter<edm::FileInPath>("effAreasConfigFile").fullPath(),
0043              params.getParameter<bool>("quadEAflag")),
0044       linearPtTerm_(params, "linearPtTerm"),
0045       quadraticPtTerm_(params, "quadraticPtTerm") {
0046   edm::InputTag rhoTag = params.getParameter<edm::InputTag>("rho");
0047   contentTags_.emplace("rho", rhoTag);
0048 }
0049 
0050 void PhoGenericQuadraticRhoPtScaledCut::setConsumes(edm::ConsumesCollector& cc) {
0051   auto rho = cc.consumes<double>(contentTags_["rho"]);
0052   contentTokens_.emplace("rho", rho);
0053 }
0054 
0055 void PhoGenericQuadraticRhoPtScaledCut::getEventContent(const edm::EventBase& ev) {
0056   ev.getByLabel(contentTags_["rho"], rhoHandle_);
0057 }
0058 
0059 CutApplicatorBase::result_type PhoGenericQuadraticRhoPtScaledCut::operator()(const reco::PhotonPtr& pho) const {
0060   const double rho = (*rhoHandle_);
0061 
0062   const float var = varFunc_(*pho);
0063 
0064   const float et = pho->et();
0065   const float absEta = std::abs(pho->superCluster()->eta());
0066   const float cutValue = constTerm_(pho) + rhoEA_.getLinearEA(absEta) * rho +
0067                          rhoEA_.getQuadraticEA(absEta) * rho * rho + linearPtTerm_(pho) * et +
0068                          quadraticPtTerm_(pho) * et * et;
0069   if (lessThan_)
0070     return var < cutValue;
0071   else
0072     return var >= cutValue;
0073 }
0074 
0075 double PhoGenericQuadraticRhoPtScaledCut::value(const reco::CandidatePtr& cand) const {
0076   reco::PhotonPtr pho(cand);
0077   return varFunc_(*pho);
0078 }